File: GfaIO.py

package info (click to toggle)
python-biopython 1.85%2Bdfsg-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 126,372 kB
  • sloc: xml: 1,047,995; python: 332,722; ansic: 16,944; sql: 1,208; makefile: 140; sh: 81
file content (213 lines) | stat: -rw-r--r-- 7,091 bytes parent folder | download
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
"""Bio.SeqIO support for the Graphical Fragment Assembly format.

This format is output by many assemblers and includes linkage information for
how the different sequences fit together, however, we just care about the
segment (sequence) information.

Documentation:
- Version 1.x: https://gfa-spec.github.io/GFA-spec/GFA1.html
- Version 2.0: https://gfa-spec.github.io/GFA-spec/GFA2.html
"""

import hashlib
import re
import warnings

from Bio import BiopythonWarning
from Bio.Seq import _UndefinedSequenceData
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord


from .Interfaces import _TextIOSource
from .Interfaces import SequenceIterator


def _check_tags(seq, tags):
    """Check a segment line's tags for inconsistencies (PRIVATE)."""
    for tag in tags:
        if tag[:2] == "LN":
            # Sequence length
            if len(seq) == 0:
                # No sequence data, set the sequence length
                seq._data = _UndefinedSequenceData(int(tag[5:]))
            elif int(tag[5:]) != len(seq):
                warnings.warn(
                    f"Segment line has incorrect length. Expected {tag[5:]} but got {len(seq)}.",
                    BiopythonWarning,
                )
        elif tag[:2] == "SH":
            # SHA256 checksum
            checksum = hashlib.sha256(str(seq).encode()).hexdigest()
            if checksum.upper() != tag[5:]:
                warnings.warn(
                    f"Segment line has incorrect checksum. Expected {tag[5:]} but got {checksum}.",
                    BiopythonWarning,
                )


def _tags_to_annotations(tags):
    """Build an annotations dictionary from a list of tags (PRIVATE)."""
    annotations = {}
    for tag in tags:
        parts = tag.split(":")
        if len(parts) < 3:
            raise ValueError(f"Segment line has invalid tag: {tag}.")
        if re.fullmatch(r"[A-Za-z][A-Za-z0-9]", parts[0]) is None:
            warnings.warn(
                f"Tag has invalid name: {parts[0]}. Are they tab delimited?",
                BiopythonWarning,
            )
        parts[2] = ":".join(parts[2:])  # tag value may contain : characters
        annotations[parts[0]] = (parts[1], parts[2])

        # Check type of the tag and raise warning on a mismatch. These RegExs
        # are part of the 1.0 standard.
        if parts[1] not in "AifZJHB":
            warnings.warn(f"Tag has invalid type: {parts[1]}", BiopythonWarning)
        elif parts[1] == "A" and re.fullmatch(r"[!-~]", parts[2]) is None:
            warnings.warn(
                f"Tag has incorrect type. Expected printable character, got {parts[2]}.",
                BiopythonWarning,
            )
        elif parts[1] == "i" and re.fullmatch(r"[-+]?[0-9]+", parts[2]) is None:
            warnings.warn(
                f"Tag has incorrect type. Expected signed integer, got {parts[2]}.",
                BiopythonWarning,
            )
        elif (
            parts[1] == "f"
            and re.fullmatch(r"[-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?", parts[2])
            is None
        ):
            warnings.warn(
                f"Tag has incorrect type. Expected float, got {parts[2]}.",
                BiopythonWarning,
            )
        elif parts[1] == "Z" and re.fullmatch(r"[ !-~]+", parts[2]) is None:
            warnings.warn(
                f"Tag has incorrect type. Expected printable string, got {parts[2]}.",
                BiopythonWarning,
            )
        elif parts[1] == "J" and re.fullmatch(r"[ !-~]+", parts[2]) is None:
            warnings.warn(
                f"Tag has incorrect type. Expected JSON excluding new-line and tab characters, got {parts[2]}.",
                BiopythonWarning,
            )
        elif parts[1] == "H" and re.fullmatch(r"[0-9A-F]+", parts[2]) is None:
            warnings.warn(
                f"Tag has incorrect type. Expected byte array in hex format, got {parts[2]}.",
                BiopythonWarning,
            )
        elif (
            parts[1] == "B"
            and re.fullmatch(
                r"[cCsSiIf](,[-+]?[0-9]*\.?[0-9]+([eE][-+]?[0-9]+)?)+", parts[2]
            )
            is None
        ):
            warnings.warn(
                f"Tag has incorrect type. Expected array of integers or floats, got {parts[2]}.",
                BiopythonWarning,
            )
    return annotations


class Gfa1Iterator(SequenceIterator):
    """Parser for GFA 1.x files.

    Documentation: https://gfa-spec.github.io/GFA-spec/GFA1.html
    """

    modes = "t"

    def __init__(
        self,
        source: _TextIOSource,
    ) -> None:
        """Iterate over a GFA file as SeqRecord objects.

        Arguments:
         - source - input stream opened in text mode, or a path to a file
        """
        super().__init__(source, fmt="GFA 1.0")

    def __next__(self):
        for line in self.stream:
            if line == "\n":
                warnings.warn("GFA data has a blank line.", BiopythonWarning)
                continue

            fields = line.strip("\n").split("\t")
            if fields[0] == "S":
                break
        else:
            raise StopIteration
        if len(fields) < 3:
            raise ValueError(
                f"Segment line must have name and sequence fields: {line}."
            )

        if fields[2] == "*":
            seq = Seq(None, length=0)
        else:
            seq = Seq(fields[2])

        tags = fields[3:]
        _check_tags(seq, tags)
        annotations = _tags_to_annotations(tags)

        return SeqRecord(seq, id=fields[1], name=fields[1], annotations=annotations)


class Gfa2Iterator(SequenceIterator):
    """Parser for GFA 2.0 files.

    Documentation for version 2: https://gfa-spec.github.io/GFA-spec/GFA2.html
    """

    modes = "t"

    def __init__(
        self,
        source: _TextIOSource,
    ) -> None:
        """Iterate over a GFA file as SeqRecord objects.

        Arguments:
         - source - input stream opened in text mode, or a path to a file
        """
        super().__init__(source, fmt="GFA 2.0")

    def __next__(self):
        for line in self.stream:
            if line == "\n":
                warnings.warn("GFA data has a blank line.", BiopythonWarning)
                continue

            fields = line.strip("\n").split("\t")
            if fields[0] == "S":
                break
        else:
            raise StopIteration
        if len(fields) < 4:
            raise ValueError(
                f"Segment line must have name, length, and sequence fields: {line}."
            )
        try:
            int(fields[2])
        except ValueError:
            raise ValueError(
                f"Segment line must have an integer length: {line}."
            ) from None

        if fields[3] == "*":
            seq = Seq(None, length=0)
        else:
            seq = Seq(fields[3])

        tags = fields[4:]
        _check_tags(seq, tags)
        annotations = _tags_to_annotations(tags)

        return SeqRecord(seq, id=fields[1], name=fields[1], annotations=annotations)