Mercurial > repos > peterjc > align_back_trans
annotate tools/align_back_trans/align_back_trans.py @ 4:c8469274d136 draft
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
author | peterjc |
---|---|
date | Wed, 10 May 2017 12:05:28 -0400 |
parents | de803005027f |
children | 2c32e8a8990f |
rev | line source |
---|---|
0 | 1 #!/usr/bin/env python |
2 """Back-translate a protein alignment to nucleotides | |
3 | |
4 This tool is a short Python script (using Biopython library functions) to | |
5 load a protein alignment, and matching nucleotide FASTA file of unaligned | |
6 sequences, in order to produce a codon aware nucleotide alignment - which | |
7 can be viewed as a back translation. | |
8 | |
9 The development repository for this tool is here: | |
10 | |
3 | 11 * https://github.com/peterjc/pico_galaxy/tree/master/tools/align_back_trans |
0 | 12 |
13 This tool is available with a Galaxy wrapper from the Galaxy Tool Shed at: | |
14 | |
15 * http://toolshed.g2.bx.psu.edu/view/peterjc/align_back_trans | |
16 | |
17 See accompanying text file for licence details (MIT licence). | |
18 """ | |
19 | |
20 import sys | |
4
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
21 |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
22 from Bio import AlignIO |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
23 from Bio import SeqIO |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
24 |
0 | 25 from Bio.Align import MultipleSeqAlignment |
4
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
26 from Bio.Alphabet import generic_protein |
0 | 27 from Bio.Data.CodonTable import ambiguous_generic_by_id |
4
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
28 from Bio.Seq import Seq |
0 | 29 |
30 if "-v" in sys.argv or "--version" in sys.argv: | |
3 | 31 print "v0.0.7" |
0 | 32 sys.exit(0) |
33 | |
34 | |
35 def check_trans(identifier, nuc, prot, table): | |
36 """Returns nucleotide sequence if works (can remove trailing stop)""" | |
37 if len(nuc) % 3: | |
3 | 38 sys.exit("Nucleotide sequence for %s is length %i (not a multiple of three)" |
1 | 39 % (identifier, len(nuc))) |
0 | 40 |
41 p = str(prot).upper().replace("*", "X") | |
42 t = str(nuc.translate(table)).upper().replace("*", "X") | |
43 if len(t) == len(p) + 1: | |
44 if str(nuc)[-3:].upper() in ambiguous_generic_by_id[table].stop_codons: | |
3 | 45 # Allow this... |
0 | 46 t = t[:-1] |
3 | 47 nuc = nuc[:-3] # edit return value |
0 | 48 if len(t) != len(p): |
2
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
49 err = ("Inconsistent lengths for %s, ungapped protein %i, " |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
50 "tripled %i vs ungapped nucleotide %i." % |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
51 (identifier, len(p), len(p) * 3, len(nuc))) |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
52 if t.endswith(p): |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
53 err += "\nThere are %i extra nucleotides at the start." % (len(t) - len(p)) |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
54 elif t.startswith(p): |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
55 err += "\nThere are %i extra nucleotides at the end." % (len(t) - len(p)) |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
56 elif p in t: |
3 | 57 # TODO - Calculate and report the number to trim at start and end? |
2
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
58 err += "\nHowever, protein sequence found within translated nucleotides." |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
59 elif p[1:] in t: |
9fbf29a8c12b
v0.0.6 use format_source; v0.0.5 more explicit error msg, citation info
peterjc
parents:
1
diff
changeset
|
60 err += "\nHowever, ignoring first amino acid, protein sequence found within translated nucleotides." |
3 | 61 sys.exit(err) |
0 | 62 |
63 if t == p: | |
64 return nuc | |
65 elif p.startswith("M") and "M" + t[1:] == p: | |
3 | 66 # Close, was there a start codon? |
0 | 67 if str(nuc[0:3]).upper() in ambiguous_generic_by_id[table].start_codons: |
68 return nuc | |
69 else: | |
3 | 70 sys.exit("Translation check failed for %s\n" |
0 | 71 "Would match if %s was a start codon (check correct table used)\n" |
72 % (identifier, nuc[0:3].upper())) | |
73 else: | |
3 | 74 # Allow * vs X here? e.g. internal stop codons |
75 m = "".join("." if x == y else "!" for (x, y) in zip(p, t)) | |
0 | 76 if len(prot) < 70: |
77 sys.stderr.write("Protein: %s\n" % p) | |
78 sys.stderr.write(" %s\n" % m) | |
79 sys.stderr.write("Translation: %s\n" % t) | |
80 else: | |
81 for offset in range(0, len(p), 60): | |
3 | 82 sys.stderr.write("Protein: %s\n" % p[offset:offset + 60]) |
83 sys.stderr.write(" %s\n" % m[offset:offset + 60]) | |
84 sys.stderr.write("Translation: %s\n\n" % t[offset:offset + 60]) | |
85 sys.exit("Translation check failed for %s\n" % identifier) | |
86 | |
0 | 87 |
88 def sequence_back_translate(aligned_protein_record, unaligned_nucleotide_record, gap, table=0): | |
3 | 89 # TODO - Separate arguments for protein gap and nucleotide gap? |
0 | 90 if not gap or len(gap) != 1: |
91 raise ValueError("Please supply a single gap character") | |
92 | |
93 alpha = unaligned_nucleotide_record.seq.alphabet | |
94 if hasattr(alpha, "gap_char"): | |
95 gap_codon = alpha.gap_char * 3 | |
96 assert len(gap_codon) == 3 | |
97 else: | |
98 from Bio.Alphabet import Gapped | |
99 alpha = Gapped(alpha, gap) | |
3 | 100 gap_codon = gap * 3 |
0 | 101 |
102 ungapped_protein = aligned_protein_record.seq.ungap(gap) | |
103 ungapped_nucleotide = unaligned_nucleotide_record.seq | |
104 if table: | |
105 ungapped_nucleotide = check_trans(aligned_protein_record.id, ungapped_nucleotide, ungapped_protein, table) | |
106 elif len(ungapped_protein) * 3 != len(ungapped_nucleotide): | |
3 | 107 sys.exit("Inconsistent lengths for %s, ungapped protein %i, " |
0 | 108 "tripled %i vs ungapped nucleotide %i" % |
109 (aligned_protein_record.id, | |
110 len(ungapped_protein), | |
111 len(ungapped_protein) * 3, | |
112 len(ungapped_nucleotide))) | |
113 | |
114 seq = [] | |
115 nuc = str(ungapped_nucleotide) | |
116 for amino_acid in aligned_protein_record.seq: | |
117 if amino_acid == gap: | |
118 seq.append(gap_codon) | |
119 else: | |
120 seq.append(nuc[:3]) | |
121 nuc = nuc[3:] | |
122 assert not nuc, "Nucleotide sequence for %r longer than protein %r" \ | |
123 % (unaligned_nucleotide_record.id, aligned_protein_record.id) | |
124 | |
3 | 125 aligned_nuc = unaligned_nucleotide_record[:] # copy for most annotation |
126 aligned_nuc.letter_annotation = {} # clear this | |
127 aligned_nuc.seq = Seq("".join(seq), alpha) # replace this | |
0 | 128 assert len(aligned_protein_record.seq) * 3 == len(aligned_nuc) |
129 return aligned_nuc | |
130 | |
3 | 131 |
0 | 132 def alignment_back_translate(protein_alignment, nucleotide_records, key_function=None, gap=None, table=0): |
133 """Thread nucleotide sequences onto a protein alignment.""" | |
3 | 134 # TODO - Separate arguments for protein and nucleotide gap characters? |
0 | 135 if gap is None: |
136 gap = "-" | |
137 | |
138 aligned = [] | |
4
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
139 try: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
140 if key_function is None: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
141 for protein in protein_alignment: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
142 nucleotide = nucleotide_records[protein.id] |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
143 aligned.append(sequence_back_translate(protein, nucleotide, gap, table)) |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
144 else: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
145 for protein in protein_alignment: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
146 nucleotide = nucleotide_records[key_function(protein.id)] |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
147 aligned.append(sequence_back_translate(protein, nucleotide, gap, table)) |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
148 except KeyError: |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
149 raise ValueError("Could not find nucleotide sequence for protein %r" |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
150 % protein.id) |
c8469274d136
v0.0.8 Using Biopython 1.67 from Tool Shed or Conda package
peterjc
parents:
3
diff
changeset
|
151 |
0 | 152 return MultipleSeqAlignment(aligned) |
153 | |
154 | |
155 if len(sys.argv) == 4: | |
156 align_format, prot_align_file, nuc_fasta_file = sys.argv[1:] | |
157 nuc_align_file = sys.stdout | |
158 table = 0 | |
159 elif len(sys.argv) == 5: | |
160 align_format, prot_align_file, nuc_fasta_file, nuc_align_file = sys.argv[1:] | |
161 table = 0 | |
162 elif len(sys.argv) == 6: | |
163 align_format, prot_align_file, nuc_fasta_file, nuc_align_file, table = sys.argv[1:] | |
164 else: | |
3 | 165 sys.exit("""This is a Python script for 'back-translating' a protein alignment, |
0 | 166 |
167 It requires three, four or five arguments: | |
168 - alignment format (e.g. fasta, clustal), | |
169 - aligned protein file (in specified format), | |
170 - unaligned nucleotide file (in fasta format). | |
171 - aligned nucleotiode output file (in same format), optional. | |
172 - NCBI translation table (0 for none), optional | |
173 | |
174 The nucleotide alignment is printed to stdout if no output filename is given. | |
175 | |
176 Example usage: | |
177 | |
178 $ python align_back_trans.py fasta demo_prot_align.fasta demo_nucs.fasta demo_nuc_align.fasta | |
179 | |
180 Warning: If the output file already exists, it will be overwritten. | |
181 | |
182 This script is available with sample data and a Galaxy wrapper here: | |
183 https://github.com/peterjc/pico_galaxy/tree/master/tools/align_back_trans | |
184 http://toolshed.g2.bx.psu.edu/view/peterjc/align_back_trans | |
185 """) | |
186 | |
187 try: | |
188 table = int(table) | |
3 | 189 except ValueError: |
190 sys.exit("Bad table argument %r" % table) | |
0 | 191 |
192 prot_align = AlignIO.read(prot_align_file, align_format, alphabet=generic_protein) | |
193 nuc_dict = SeqIO.index(nuc_fasta_file, "fasta") | |
194 nuc_align = alignment_back_translate(prot_align, nuc_dict, gap="-", table=table) | |
195 AlignIO.write(nuc_align, nuc_align_file, align_format) |