File size: 18,265 Bytes
b7731cd
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
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
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
# Copyright 2011, 2012 by Andrew Sczesnak.  All rights reserved.
# Revisions Copyright 2011, 2017 by Peter Cock.  All rights reserved.
# Revisions Copyright 2014, 2015 by Adam Novak.  All rights reserved.
# Revisions Copyright 2015, 2017 by Blaise Li.  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 the "maf" multiple alignment format.

The Multiple Alignment Format, described by UCSC, stores a series of
multiple alignments in a single file. It is suitable for whole-genome
to whole-genome alignments, metadata such as source chromosome, start
position, size, and strand can be stored.

See http://genome.ucsc.edu/FAQ/FAQformat.html#format5

You are expected to use this module via the Bio.Align functions.

Coordinates in the MAF format are defined in terms of zero-based start
positions (like Python) and aligning region sizes.

A minimal aligned region of length one and starting at first position in the
source sequence would have ``start == 0`` and ``size == 1``.

As we can see on this example, ``start + size`` will give one more than the
zero-based end position. We can therefore manipulate ``start`` and
``start + size`` as python list slice boundaries.
"""
import shlex
import itertools


from Bio.Align import Alignment
from Bio.Align import interfaces
from Bio.Seq import Seq, reverse_complement
from Bio.SeqRecord import SeqRecord


class AlignmentWriter(interfaces.AlignmentWriter):
    """Accepts Alignment objects, writes a MAF file."""

    fmt = "MAF"

    def _write_trackline(self, metadata):
        stream = self.stream
        stream.write("track")
        for key, value in metadata.items():
            if key in ("name", "description", "frames"):
                pass
            elif key == "mafDot":
                if value not in ("on", "off"):
                    raise ValueError(
                        "mafDot value must be 'on' or 'off' (received '%s')" % value
                    )
            elif key == "visibility":
                if value not in ("dense", "pack", "full"):
                    raise ValueError(
                        "visibility value must be 'dense', 'pack', or 'full' (received '%s')"
                        % value
                    )
            elif key == "speciesOrder":
                value = " ".join(value)
            else:
                continue
            if " " in value:
                value = '"%s"' % value
            stream.write(f" {key}={value}")
        stream.write("\n")

    def write_header(self, alignments):
        """Write the MAF header."""
        stream = self.stream
        try:
            metadata = alignments.metadata
        except AttributeError:
            metadata = {"MAF Version": "1"}
        track_keys = (
            "name",
            "description",
            "frames",
            "mafDot",
            "visibility",
            "speciesOrder",
        )
        for key in track_keys:
            if key in metadata:
                self._write_trackline(metadata)
                break
        stream.write("##maf")
        for key, value in metadata.items():
            if key in track_keys:
                continue
            if key == "Comments":
                continue
            if key == "MAF Version":
                if value != "1":
                    raise ValueError("MAF version must be 1")
                key = "version"
            elif key == "Scoring":
                key = "scoring"
            elif key == "Program":
                key = "program"
            else:
                raise ValueError("Unexpected key '%s' for header" % key)
            stream.write(f" {key}={value}")
        stream.write("\n")
        comments = metadata.get("Comments")
        if comments is not None:
            for comment in comments:
                stream.write(f"# {comment}\n")
            stream.write("\n")

    def _format_score_line(self, alignment, annotations):
        try:
            score = alignment.score
        except AttributeError:
            line = "a"
        else:
            line = f"a score={score:.6f}"
        value = annotations.get("pass")
        if value is not None:
            line += f" pass={value}"
        return line + "\n"

    def format_alignment(self, alignment):
        """Return a string with a single alignment formatted as a MAF block."""
        if not isinstance(alignment, Alignment):
            raise TypeError("Expected an Alignment object")
        try:
            alignment_annotations = alignment.annotations
        except AttributeError:
            alignment_annotations = {}
        lines = []
        line = self._format_score_line(alignment, alignment_annotations)
        lines.append(line)
        name_width = 0
        start_width = 0
        size_width = 0
        length_width = 0
        n = len(alignment.sequences)
        for i in range(n):
            record = alignment.sequences[i]
            coordinates = alignment.coordinates[i]
            try:
                name = record.id
            except AttributeError:
                name = "sequence_%d" % i
            start = coordinates[0]
            end = coordinates[-1]
            length = len(record)
            if start < end:
                size = end - start
            else:
                size = start - end
                start = length - start
            name_width = max(name_width, len(name))
            start_width = max(start_width, len(str(start)))
            size_width = max(size_width, len(str(size)))
            length_width = max(length_width, len(str(length)))
        for i, empty in enumerate(alignment_annotations.get("empty", [])):
            record, segment, status = empty
            try:
                name = record.id
            except AttributeError:
                name = "sequence_%d" % (i + n)
            name_width = max(name_width, len(name))
            start, end = segment
            length = len(record.seq)
            if start <= end:
                size = end - start
            else:
                size = start - end
                start = length - start
            start_width = max(start_width, len(str(start)))
            size_width = max(size_width, len(str(size)))
            length_width = max(length_width, len(str(length)))
        quality_width = name_width + start_width + size_width + length_width + 5
        for i in range(n):
            record = alignment.sequences[i]
            coordinates = alignment.coordinates[i]
            try:
                record_id = record.id
            except AttributeError:
                record_id = "sequence_%d" % i
            start = coordinates[0]
            end = coordinates[-1]
            length = len(record)
            if start < end:
                size = end - start
                strand = "+"
            else:
                size = start - end
                start = length - start
                strand = "-"
            text = alignment[i]
            name = record_id.ljust(name_width)
            start = str(start).rjust(start_width)
            size = str(size).rjust(size_width)
            length = str(length).rjust(length_width)
            line = f"s {name} {start} {size} {strand} {length} {text}\n"
            lines.append(line)
            try:
                annotations = record.annotations
            except AttributeError:
                annotations = None
            if annotations is not None:
                quality = annotations.get("quality")
                if quality is not None:
                    gapped_quality = ""
                    j = 0
                    for letter in text:
                        if letter == "-":
                            gapped_quality += "-"
                        else:
                            gapped_quality += quality[j]
                            j += 1
                    name = record_id.ljust(quality_width)
                    line = f"q {name} {gapped_quality}\n"
                    lines.append(line)
                try:
                    leftStatus = annotations["leftStatus"]
                    leftCount = annotations["leftCount"]
                    rightStatus = annotations["rightStatus"]
                    rightCount = annotations["rightCount"]
                except KeyError:
                    pass
                else:
                    name = record_id.ljust(name_width)
                    line = f"i {name} {leftStatus} {leftCount} {rightStatus} {rightCount}\n"
                    lines.append(line)
        for i, empty in enumerate(alignment_annotations.get("empty", [])):
            record, segment, status = empty
            try:
                name = record.id
            except AttributeError:
                name = "sequence_%d" % (i + n)
            name = name.ljust(name_width)
            start, end = segment
            length = len(record.seq)
            if start <= end:
                size = end - start
                strand = "+"
            else:
                size = start - end
                start = length - start
                strand = "-"
            start = str(start).rjust(start_width)
            size = str(size).rjust(size_width)
            length = str(length).rjust(length_width)
            line = f"e {name} {start} {size} {strand} {length} {status}\n"
            lines.append(line)
        lines.append("\n")
        return "".join(lines)


class AlignmentIterator(interfaces.AlignmentIterator):
    """Alignment iterator for Multiple Alignment Format files.

    The file may contain multiple concatenated alignments, which are loaded
    and returned incrementally.

    File meta-data are stored in the ``.metadata`` attribute of the returned
    iterator.  Alignment annotations are stored in the ``.annotations``
    attribute of the ``Alignment`` object, except for the alignment score,
    which is stored as an attribute.  Sequence information of empty parts in
    the alignment block (sequences that connect the previous alignment block to
    the next alignment block, but do not align to the current alignment block)
    is stored in the alignment annotations under the ``"empty"`` key.
    Annotations specific to each line in the alignment are stored in the
    ``.annotations`` attribute of the corresponding sequence record.
    """

    fmt = "MAF"

    status_characters = ("C", "I", "N", "n", "M", "T")
    empty_status_characters = ("C", "I", "M", "n")

    def _read_header(self, stream):
        metadata = {}
        line = next(stream)
        if line.startswith("track "):
            words = shlex.split(line)
            for word in words[1:]:
                key, value = word.split("=")
                if key in ("name", "description", "frames"):
                    pass
                elif key == "mafDot":
                    if value not in ("on", "off"):
                        raise ValueError(
                            "Variable mafDot in track line has unexpected value '%s'"
                            % value
                        )
                elif key == "visibility":
                    if value not in ("dense", "pack", "full"):
                        raise ValueError(
                            "Variable visibility in track line has unexpected value '%s'"
                            % value
                        )
                elif key == "speciesOrder":
                    value = value.split()
                else:
                    raise ValueError("Unexpected variable '%s' in track line" % key)
                metadata[key] = value
            line = next(stream)
        words = line.split()
        if words[0] != "##maf":
            raise ValueError("header line does not start with ##maf")
        for word in words[1:]:
            key, value = word.split("=")
            if key == "version":
                key = "MAF Version"
            elif key == "scoring":
                key = "Scoring"
            elif key == "program":
                key = "Porgram"
            else:
                raise ValueError("Unexpected variable '%s' in header line" % key)
            metadata[key] = value
        if metadata.get("MAF Version") != "1":
            raise ValueError("MAF version must be 1")
        comments = []
        for line in stream:
            if line.strip():
                if not line.startswith("#"):
                    self._line = line
                    break
                comment = line[1:].strip()
                comments.append(comment)
        else:
            self._close()
        if comments:
            metadata["Comments"] = comments
        self.metadata = metadata

    def _read_next_alignment(self, stream):
        line = self._line
        if line is None:
            return
        lines = itertools.chain([line], stream)
        alignment = self._create_alignment(lines)
        return alignment

    def _create_alignment(self, lines):
        records = []
        strands = []
        column_annotations = {}
        aligned_sequences = []
        annotations = {}
        line = next(lines)
        assert line.startswith("a")
        words = line[1:].split()
        for word in words:
            key, value = word.split("=")
            if key == "score":
                score = float(value)
            elif key == "pass":
                value = int(value)
                if value <= 0:
                    raise ValueError("pass value must be positive (found %d)" % value)
                annotations["pass"] = value
            else:
                raise ValueError("Unknown annotation variable '%s'" % key)

        for line in lines:
            if line.startswith("a"):
                self._line = line
                break
            elif line.startswith("s "):
                words = line.strip().split()
                if len(words) != 7:
                    raise ValueError(
                        "Error parsing alignment - 's' line must have 7 fields"
                    )
                src = words[1]
                start = int(words[2])
                size = int(words[3])
                strand = words[4]
                srcSize = int(words[5])
                text = words[6]
                for gap_char in ".=_":
                    text = text.replace(gap_char, "-")
                aligned_sequences.append(text)
                sequence = text.replace("-", "")
                if len(sequence) != size:
                    raise ValueError(
                        "sequence size is incorrect (found %d, expected %d)"
                        % (len(sequence), size)
                    )
                if strand == "-":
                    sequence = reverse_complement(sequence)
                    start = srcSize - start - size
                seq = Seq({start: sequence}, length=srcSize)
                record = SeqRecord(seq, id=src, name="", description="")
                records.append(record)
                strands.append(strand)
            elif line.startswith("i "):
                words = line.strip().split()
                assert len(words) == 6
                assert words[1] == src  # from the previous "s" line
                leftStatus = words[2]
                leftCount = int(words[3])
                rightStatus = words[4]
                rightCount = int(words[5])
                assert leftStatus in AlignmentIterator.status_characters
                assert rightStatus in AlignmentIterator.status_characters
                record.annotations["leftStatus"] = leftStatus
                record.annotations["leftCount"] = leftCount
                record.annotations["rightStatus"] = rightStatus
                record.annotations["rightCount"] = rightCount
            elif line.startswith("e"):
                words = line[1:].split()
                assert len(words) == 6
                src = words[0]
                start = int(words[1])
                size = int(words[2])
                strand = words[3]
                srcSize = int(words[4])
                status = words[5]
                assert status in AlignmentIterator.empty_status_characters
                sequence = Seq(None, length=srcSize)
                record = SeqRecord(sequence, id=src, name="", description="")
                end = start + size
                if strand == "+":
                    segment = (start, end)
                else:
                    segment = (srcSize - start, srcSize - end)
                empty = (record, segment, status)
                annotation = annotations.get("empty")
                if annotation is None:
                    annotation = []
                    annotations["empty"] = annotation
                annotation.append(empty)
            elif line.startswith("q "):
                words = line.strip().split()
                assert len(words) == 3
                assert words[1] == src  # from the previous "s" line
                value = words[2].replace("-", "")
                record.annotations["quality"] = value
            elif not line.strip():
                # reached the end of the alignment, but keep reading until we
                # find the next alignment
                continue
            else:
                raise ValueError(f"Error parsing alignment - unexpected line:\n{line}")
        else:
            self._line = None
        coordinates = Alignment.infer_coordinates(aligned_sequences)
        for record, strand, row in zip(records, strands, coordinates):
            if strand == "-":
                row[:] = row[-1] - row[0] - row
            start = record.seq.defined_ranges[0][0]
            row += start
        alignment = Alignment(records, coordinates)
        if annotations is not None:
            alignment.annotations = annotations
        if column_annotations is not None:
            alignment.column_annotations = column_annotations
        if score is not None:
            alignment.score = score
        return alignment