File: nextorf.py

package info (click to toggle)
python-biopython 1.68%2Bdfsg-3
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 46,860 kB
  • ctags: 13,237
  • sloc: python: 160,306; xml: 93,216; ansic: 9,118; sql: 1,208; makefile: 155; sh: 63
file content (292 lines) | stat: -rwxr-xr-x 10,095 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
#!/usr/bin/env python
# Created: Thu Feb 15 14:22:12 2001
# Last changed: Time-stamp: <01/02/18 11:16:42 thomas>
# Copyright 2000 by Thomas Sicheritz-Ponten.  All rights reserved.
# This code is part of the Biopython distribution and governed by its
# license.  Please see the LICENSE file that should have been included
# as part of this package.
# Authors: Thomas Sicheritz-Ponten and Jan O. Andersson
# thomas@cbs.dtu.dk, http://www.cbs.dtu.dk/thomas
# Jan.O.Andersson@home.se
# File: nextorf.py

from __future__ import print_function

import re
import sys
import os
import getopt

from Bio import SeqIO
from Bio.Seq import Seq
from Bio import Alphabet
from Bio.Alphabet import IUPAC
from Bio.Data import IUPACData, CodonTable


class ProteinX(Alphabet.ProteinAlphabet):
    letters = IUPACData.extended_protein_letters + "X"

proteinX = ProteinX()


class MissingTable(object):
    def __init__(self, table):
        self._table = table

    def get(self, codon, stop_symbol):
        try:
            return self._table.get(codon, stop_symbol)
        except CodonTable.TranslationError:
            return 'X'


# Make the codon table given an existing table
def makeTableX(table):
    assert table.protein_alphabet == IUPAC.extended_protein
    return CodonTable.CodonTable(table.nucleotide_alphabet, proteinX,
                                 MissingTable(table.forward_table),
                                 table.back_table, table.start_codons,
                                 table.stop_codons)


class NextOrf(object):
    def __init__(self, file, options):
        self.options = options
        self.file = file
        self.genetic_code = int(self.options['table'])
        self.table = makeTableX(CodonTable.ambiguous_dna_by_id[self.genetic_code])
        self.counter = 0
        self.ReadFile()

    def ReadFile(self):
        handle = open(self.file)
        for record in SeqIO.parse(handle, "fasta"):
            self.header = record.id
            frame_coordinates = ''
            dir = self.options['strand']
            plus = dir in ['both', 'plus']
            minus = dir in ['both', 'minus']
            start, stop = int(self.options['start']), int(self.options['stop'])
            s = str(record.seq).upper()
            if stop > 0:
                s = s[start:stop]
            else:
                s = s[start:]
            self.seq = Seq(s, IUPAC.ambiguous_dna)
            self.length = len(self.seq)
            self.rseq = None
            CDS = []
            if plus:
                CDS.extend(self.GetCDS(self.seq))
            if minus:
                self.rseq = self.seq.reverse_complement()
                CDS.extend(self.GetCDS(self.rseq, strand=-1))
            self.Output(CDS)

    def ToFasta(self, header, seq):
        seq = re.sub('(............................................................)',
                     '\\1\n', seq)
        return '>%s\n%s' % (header, seq)

    def Gc(self, seq):
        d = {}
        for nt in 'ATGC':
            d[nt] = seq.count(nt)
        gc = d['G'] + d['C']
        if gc == 0:
            return 0
        return round(gc * 100.0 / (d['A'] + d['T'] + gc), 1)

    def Gc2(self, seq):
        l = len(seq)
        d = {}
        for nt in ['A', 'T', 'G', 'C']:
            d[nt] = [0, 0, 0]

        for i in range(0, l, 3):
            codon = seq[i:i + 3]
            if len(codon) < 3:
                codon += '  '
            for pos in range(0, 3):
                for nt in ['A', 'T', 'G', 'C']:
                    if codon[pos] == nt:
                        d[nt][pos] = d[nt][pos] + 1

        gc = {}
        gcall = 0
        nall = 0
        for i in range(0, 3):
            try:
                n = d['G'][i] + d['C'][i] + d['T'][i] + d['A'][i]
                gc[i] = (d['G'][i] + d['C'][i]) * 100.0 / n
            except:
                gc[i] = 0

            gcall = gcall + d['G'][i] + d['C'][i]
            nall += n

        gcall = 100.0 * gcall / nall
        res = '%.1f%%, %.1f%%, %.1f%%, %.1f%%' % (gcall, gc[0], gc[1], gc[2])
        return res

    def GetOrfCoordinates(self, seq):
        s = seq.data
        letters = []
        table = self.table
        get = self.table.forward_table.get
        n = len(seq)
        start_codons = self.table.start_codons
        stop_codons = self.table.stop_codons
#        print('Start codons %s' % start_codons)
#        print('Stop codons %s' % stop_codons)
        frame_coordinates = []
        for frame in range(0, 3):
            coordinates = []
            for i in range(0 + frame, n - n % 3, 3):
                codon = s[i:i + 3]
                if codon in start_codons:
                    coordinates.append((i + 1, 1, codon))
                elif codon in stop_codons:
                    coordinates.append((i + 1, 0, codon))
            frame_coordinates.append(coordinates)
        return frame_coordinates

    def GetCDS(self, seq, strand=1):
        frame_coordinates = self.GetOrfCoordinates(seq)
        START, STOP = 1, 0
        so = self.options
        nostart = so['nostart']
        minlength, maxlength = int(so['minlength']), int(so['maxlength'])
        CDS = []
        f = 0
        for frame in frame_coordinates:
            f += 1
            start_site = 0
            if nostart == '1':
                start_site = 1
            frame.append((self.length, 0, 'XXX'))
            for pos, codon_type, codon in frame:
                if codon_type == START:
                    if start_site == 0:
                        start_site = pos
                elif codon_type == STOP:
                    if start_site == 0:
                        continue
#                    if codon == 'XXX': print('do something')
                    stop = pos + 2
#                    print("stop")
                    length = stop - start_site + 1
                    if length >= minlength and length <= maxlength:
                        if nostart == '1' and start_site == 1:
                            start_site = start_site + f - 1
                        if codon == 'XXX':
                            stop = start_site \
                                + 3 * ((int((stop - 1) - start_site) // 3))
                        s = seq[start_site - 1:stop]
                        CDS.append((start_site, stop, length, s, strand * f))
                        start_site = 0
                        if nostart == '1':
                            start_site = stop + 1
                    elif length < minlength or length > maxlength:
                        start_site = 0
                        if nostart == '1':
                            start_site = stop + 1
                    del stop
        return CDS

    def Output(self, CDS):
        out = self.options['output']
        seqs = (self.seq, self.rseq)
        n = len(self.seq)
        for start, stop, length, subs, strand in CDS:
            self.counter += 1
            if strand > 0:
                head = 'orf_%s:%s:%d:%d:%d' % (self.counter, self.header,
                                               strand, start, stop)
            if strand < 0:
                head = 'orf_%s:%s:%d:%d:%d' % (self.counter, self.header,
                                               strand,
                                               n - stop + 1,
                                               n - start + 1)
            if self.options['gc']:
                head = '%s:%s' % (head, self.Gc2(subs.data))

            if out == 'aa':
                orf = subs.translate(table=self.genetic_code)
                print(self.ToFasta(head, orf.data))
            elif out == 'nt':
                print(self.ToFasta(head, subs.data))
            elif out == 'pos':
                print(head)


def help():
    global options
    print('Usage: %s (<options>) <FASTA file>' % sys.argv[0])
    print("")
    print('Options:                                                       default')
    print('--start       Start position in sequence                             0')
    print('--stop        Stop position in sequence            (end of seqence)')
    print('--minlength   Minimum length of orf in bp                          100')
    print('--maxlength   Maximum length of orf in bp, default           100000000')
    print('--strand      Strand to analyse [both, plus, minus]               both')
    print('--frame       Frame to analyse [1 2 3]                             all')
    print('--noframe     Ignore start codons [0 1]                              0')
    print('--output      Output to generate [aa nt pos]                        aa')
    print('--gc          Creates GC statistics of ORF [0 1]                     0')
    print('--table       Genetic code to use (see below)                        1')

#    for a,b in options.items():
#        print("\t%s %s" % (a, b)
#    print("")
    print("\nNCBI's Codon Tables:")
    for key, table in CodonTable.ambiguous_dna_by_id.items():
        print('\t%s %s' % (key, table._codon_table.names[0]))
    print('\ne.g.')
    print('./nextorf.py --minlength 5 --strand plus --output nt --gc 1 testjan.fas')
    sys.exit(0)


options = {
    'start': 0,
    'stop': 0,
    'minlength': 100,
    'maxlength': 100000000,
    'strand': 'both',
    'output': 'aa',
    'frames': [1, 2, 3],
    'gc': 0,
    'nostart': 0,
    'table': 1,
}

if __name__ == '__main__':
    args = sys.argv[1:]
    show_help = len(sys.argv) <= 1

    shorts = 'hv'
    longs = [x + '=' for x in options] + ['help']

    optlist, args = getopt.getopt(args, shorts, longs)
    if show_help:
        help()

    for arg in optlist:
        if arg[0] == '-h' or arg[0] == '--help':
            help()
            sys.exit(0)
        for key in options:
            if arg[1].lower() == 'no':
                arg[1] = 0
            elif arg[1].lower() == 'yes':
                arg[1] = 1

            if arg[0][2:] == key:
                options[key] = arg[1]

        if arg[0] == '-v':
            print('OPTIONS %s' % options)

    file = args[0]
    nextorf = NextOrf(file, options)