/
clustal.py
251 lines (220 loc) · 9.11 KB
/
clustal.py
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
# Copyright 2006-2016 by Peter Cock. All rights reserved.
#
# This file is part of the Biopython distribution and governed by your
# choice of the "Biopython License Agreement" or the "BSD 3-Clause License".
# Please see the LICENSE file that should have been included as part of this
# package.
"""Bio.Align support for "clustal" output from CLUSTAL W and other tools.
You are expected to use this module via the Bio.Align functions (or the
Bio.SeqIO functions if you are interested in the sequences only).
"""
import Bio
from Bio.Align import Alignment
from Bio.Align import interfaces
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord
class AlignmentWriter(interfaces.AlignmentWriter):
"""Clustalw alignment writer."""
fmt = "Clustal"
def write_header(self, stream, alignments):
"""Use this to write the file header."""
try:
metadata = alignments.metadata
program = metadata["Program"]
except (AttributeError, KeyError):
program = "Biopython"
version = Bio.__version__
else:
version = metadata.get("Version", "")
line = f"{program} {version} multiple sequence alignment\n"
stream.write(line)
stream.write("\n")
stream.write("\n")
def format_alignment(self, alignment):
"""Return a string with a single alignment in the Clustal format."""
nseqs, length = alignment.shape
if nseqs == 0:
raise ValueError("Must have at least one sequence")
if length == 0:
raise ValueError("Non-empty sequences are required")
try:
column_annotations = alignment.column_annotations
except AttributeError:
consensus = None
else:
consensus = column_annotations.get("clustal_consensus")
gapped_sequences = list(alignment)
names = []
for i, sequence in enumerate(alignment.sequences):
try:
name = sequence.id
except AttributeError:
name = "sequence_%d" % i # Clustal format doesn't allow an empty string
else:
# when we output, we do a nice 80 column output, although
# this may result in truncation of the ids. Also, make sure
# we don't get any spaces in the record identifier when output
# in the file by replacing them with underscores.
name = name[:30].replace(" ", "_")
name = name.ljust(36)
names.append(name)
lines = []
start = 0
while start != length:
# calculate the number of letters to show, which will
# be less if we are at the end of the alignment.
stop = start + 50
if stop > length:
stop = length
for name, gapped_sequence in zip(names, gapped_sequences):
line = f"{name}{gapped_sequence[start:stop]}\n"
lines.append(line)
# now we need to print out the star info, if we've got it
if consensus is not None:
line = " " * 36 + consensus[start:stop] + "\n"
lines.append(line)
lines.append("\n")
start = stop
lines.append("\n")
return "".join(lines)
class AlignmentIterator(interfaces.AlignmentIterator):
"""Clustalw alignment iterator."""
fmt = "Clustal"
def _read_header(self, stream):
try:
line = next(stream)
except StopIteration:
raise ValueError("Empty file.") from None
self.metadata = {}
# Whitelisted programs we know about
words = line.split()
known_programs = [
"CLUSTAL",
"PROBCONS",
"MUSCLE",
"MSAPROBS",
"Kalign",
"Biopython",
]
program = words[0]
if program not in known_programs:
raise ValueError(
"%s is not known to generate CLUSTAL files: %s"
% (program, ", ".join(known_programs))
)
self.metadata["Program"] = program
# find the clustal version in the header line
for word in words:
if word[0] == "(" and word[-1] == ")":
word = word[1:-1]
if word[0].isdigit():
self.metadata["Version"] = word
break
def _read_next_alignment(self, stream):
# If the alignment contains entries with the same sequence
# identifier (not a good idea - but seems possible), then this
# dictionary based parser will merge their sequences. Fix this?
ids = []
aligned_seqs = []
consensus = ""
index = None # Used to extract the consensus
# Use the first block to get the sequence identifiers
for line in stream:
if line.startswith(" "):
# Sequence consensus line...
assert len(ids) > 0
assert index is not None
length = len(aligned_seq) # noqa: F821
consensus = line[index : index + length]
break
elif line.strip():
# Sequences identifier...
fields = line.split()
# We expect there to be two fields, there can be an optional
# "sequence number" field containing the letter count.
if len(fields) < 2 or len(fields) > 3:
raise ValueError("Could not parse line:\n%s" % line)
seqid, aligned_seq = fields[:2]
ids.append(seqid)
aligned_seqs.append(aligned_seq)
# Record the sequence position to get the consensus
if index is None:
index = line.find(aligned_seq, len(seqid))
if len(fields) == 3:
# This MAY be an old style file with a letter count...
try:
count = int(fields[2])
except ValueError:
raise ValueError(
"Could not parse line, bad sequence count:\n%s" % line
) from None
if len(aligned_seq) - aligned_seq.count("-") != count:
raise ValueError(
"Could not parse line, incorrect sequence count:\n%s" % line
) from None
else:
# no consensus line
if index:
break
else:
return
assert index is not None
# Confirm all same length
length = len(aligned_seqs[0])
for aligned_seq in aligned_seqs:
assert len(aligned_seq) == length
if consensus:
assert len(consensus) == length
n = len(aligned_seqs)
i = 0
# Loop over any remaining blocks...
for line in stream:
if line.startswith(" "): # Sequence consensus line
assert index is not None
length = len(aligned_seq)
consensus += line[index : index + length]
elif not line.strip(): # Blank line
continue
else:
seqid = ids[i]
# Sequences identifier...
fields = line.split()
# We expect there to be two fields, there can be an optional
# "sequence number" field containing the letter count.
if len(fields) < 2 or len(fields) > 3:
raise ValueError("Could not parse line:\n%s" % line)
assert seqid == fields[0]
aligned_seq = fields[1]
aligned_seqs[i] += aligned_seq
if len(fields) == 3:
# This MAY be an old style file with a letter count...
try:
count = int(fields[2])
except ValueError:
raise ValueError(
"Could not parse line, bad sequence number:\n%s" % line
) from None
if len(aligned_seqs[i]) - aligned_seqs[i].count("-") != count:
raise ValueError(
"Could not parse line, incorrect sequence count:\n%s" % line
) from None
i += 1
if i == n:
i = 0
aligned_seqs = [s.encode() for s in aligned_seqs]
seqs, coordinates = Alignment.parse_printed_alignment(aligned_seqs)
records = [
SeqRecord(Seq(seq), id=seqid, description="")
for (seqid, seq) in zip(ids, seqs)
]
alignment = Alignment(records, coordinates)
if consensus:
columns = alignment.length
if len(consensus) != columns:
raise ValueError(
"Alignment has %i columns, consensus length is %i, '%s'"
% (columns, len(consensus), consensus)
)
alignment.column_annotations = {}
alignment.column_annotations["clustal_consensus"] = consensus
return alignment