Package Bio :: Package SeqIO :: Module SffIO
[hide private]
[frames] | no frames]

Source Code for Module Bio.SeqIO.SffIO

   1  # Copyright 2009-2016 by Peter Cock.  All rights reserved. 
   2  # Based on code contributed and copyright 2009 by Jose Blanca (COMAV-UPV). 
   3  # 
   4  # This code is part of the Biopython distribution and governed by its 
   5  # license.  Please see the LICENSE file that should have been included 
   6  # as part of this package. 
   7  """Bio.SeqIO support for the binary Standard Flowgram Format (SFF) file format. 
   8   
   9  SFF was designed by 454 Life Sciences (Roche), the Whitehead Institute for 
  10  Biomedical Research and the Wellcome Trust Sanger Institute. SFF was also used 
  11  as the native output format from early versions of Ion Torrent's PGM platform 
  12  as well. You are expected to use this module via the Bio.SeqIO functions under 
  13  the format name "sff" (or "sff-trim" as described below). 
  14   
  15  For example, to iterate over the records in an SFF file, 
  16   
  17      >>> from Bio import SeqIO 
  18      >>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff"): 
  19      ...     print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
  20      ... 
  21      E3MFGYR02JWQ7T 265 tcagGGTCTACATGTTGGTT... 
  22      E3MFGYR02JA6IL 271 tcagTTTTTTTTGGAAAGGA... 
  23      E3MFGYR02JHD4H 310 tcagAAAGACAAGTGGTATC... 
  24      E3MFGYR02GFKUC 299 tcagCGGCCGGGCCTCTCAT... 
  25      E3MFGYR02FTGED 281 tcagTGGTAATGGGGGGAAA... 
  26      E3MFGYR02FR9G7 261 tcagCTCCGTAAGAAGGTGC... 
  27      E3MFGYR02GAZMS 278 tcagAAAGAAGTAAGGTAAA... 
  28      E3MFGYR02HHZ8O 221 tcagACTTTCTTCTTTACCG... 
  29      E3MFGYR02GPGB1 269 tcagAAGCAGTGGTATCAAC... 
  30      E3MFGYR02F7Z7G 219 tcagAATCATCCACTTTTTA... 
  31   
  32  Each SeqRecord object will contain all the annotation from the SFF file, 
  33  including the PHRED quality scores. 
  34   
  35      >>> print("%s %i" % (record.id, len(record))) 
  36      E3MFGYR02F7Z7G 219 
  37      >>> print("%s..." % record.seq[:10]) 
  38      tcagAATCAT... 
  39      >>> print("%r..." % (record.letter_annotations["phred_quality"][:10])) 
  40      [22, 21, 23, 28, 26, 15, 12, 21, 28, 21]... 
  41   
  42  Notice that the sequence is given in mixed case, the central upper case region 
  43  corresponds to the trimmed sequence. This matches the output of the Roche 
  44  tools (and the 3rd party tool sff_extract) for SFF to FASTA. 
  45   
  46      >>> print(record.annotations["clip_qual_left"]) 
  47      4 
  48      >>> print(record.annotations["clip_qual_right"]) 
  49      134 
  50      >>> print(record.seq[:4]) 
  51      tcag 
  52      >>> print("%s...%s" % (record.seq[4:20], record.seq[120:134])) 
  53      AATCATCCACTTTTTA...CAAAACACAAACAG 
  54      >>> print(record.seq[134:]) 
  55      atcttatcaacaaaactcaaagttcctaactgagacacgcaacaggggataagacaaggcacacaggggataggnnnnnnnnnnn 
  56   
  57  The annotations dictionary also contains any adapter clip positions 
  58  (usually zero), and information about the flows. e.g. 
  59   
  60      >>> len(record.annotations) 
  61      11 
  62      >>> print(record.annotations["flow_key"]) 
  63      TCAG 
  64      >>> print(record.annotations["flow_values"][:10]) 
  65      (83, 1, 128, 7, 4, 84, 6, 106, 3, 172) 
  66      >>> print(len(record.annotations["flow_values"])) 
  67      400 
  68      >>> print(record.annotations["flow_index"][:10]) 
  69      (1, 2, 3, 2, 2, 0, 3, 2, 3, 3) 
  70      >>> print(len(record.annotations["flow_index"])) 
  71      219 
  72   
  73  Note that to convert from a raw reading in flow_values to the corresponding 
  74  homopolymer stretch estimate, the value should be rounded to the nearest 100: 
  75   
  76      >>> print("%r..." % [int(round(value, -2)) // 100 
  77      ...                  for value in record.annotations["flow_values"][:10]]) 
  78      ... 
  79      [1, 0, 1, 0, 0, 1, 0, 1, 0, 2]... 
  80   
  81  If a read name is exactly 14 alphanumeric characters, the annotations 
  82  dictionary will also contain meta-data about the read extracted by 
  83  interpretting the name as a 454 Sequencing System "Universal" Accession 
  84  Number. Note that if a read name happens to be exactly 14 alphanumeric 
  85  characters but was not generated automatically, these annotation records 
  86  will contain nonsense information. 
  87   
  88      >>> print(record.annotations["region"]) 
  89      2 
  90      >>> print(record.annotations["time"]) 
  91      [2008, 1, 9, 16, 16, 0] 
  92      >>> print(record.annotations["coords"]) 
  93      (2434, 1658) 
  94   
  95  As a convenience method, you can read the file with SeqIO format name "sff-trim" 
  96  instead of "sff" to get just the trimmed sequences (without any annotation 
  97  except for the PHRED quality scores and anything encoded in the read names): 
  98   
  99      >>> from Bio import SeqIO 
 100      >>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff-trim"): 
 101      ...     print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
 102      ... 
 103      E3MFGYR02JWQ7T 260 GGTCTACATGTTGGTTAACC... 
 104      E3MFGYR02JA6IL 265 TTTTTTTTGGAAAGGAAAAC... 
 105      E3MFGYR02JHD4H 292 AAAGACAAGTGGTATCAACG... 
 106      E3MFGYR02GFKUC 295 CGGCCGGGCCTCTCATCGGT... 
 107      E3MFGYR02FTGED 277 TGGTAATGGGGGGAAATTTA... 
 108      E3MFGYR02FR9G7 256 CTCCGTAAGAAGGTGCTGCC... 
 109      E3MFGYR02GAZMS 271 AAAGAAGTAAGGTAAATAAC... 
 110      E3MFGYR02HHZ8O 150 ACTTTCTTCTTTACCGTAAC... 
 111      E3MFGYR02GPGB1 221 AAGCAGTGGTATCAACGCAG... 
 112      E3MFGYR02F7Z7G 130 AATCATCCACTTTTTAACGT... 
 113   
 114  Looking at the final record in more detail, note how this differs to the 
 115  example above: 
 116   
 117      >>> print("%s %i" % (record.id, len(record))) 
 118      E3MFGYR02F7Z7G 130 
 119      >>> print("%s..." % record.seq[:10]) 
 120      AATCATCCAC... 
 121      >>> print("%r..." % record.letter_annotations["phred_quality"][:10]) 
 122      [26, 15, 12, 21, 28, 21, 36, 28, 27, 27]... 
 123      >>> len(record.annotations) 
 124      3 
 125      >>> print(record.annotations["region"]) 
 126      2 
 127      >>> print(record.annotations["coords"]) 
 128      (2434, 1658) 
 129      >>> print(record.annotations["time"]) 
 130      [2008, 1, 9, 16, 16, 0] 
 131   
 132  You might use the Bio.SeqIO.convert() function to convert the (trimmed) SFF 
 133  reads into a FASTQ file (or a FASTA file and a QUAL file), e.g. 
 134   
 135      >>> from Bio import SeqIO 
 136      >>> try: 
 137      ...     from StringIO import StringIO # Python 2 
 138      ... except ImportError: 
 139      ...     from io import StringIO # Python 3 
 140      ... 
 141      >>> out_handle = StringIO() 
 142      >>> count = SeqIO.convert("Roche/E3MFGYR02_random_10_reads.sff", "sff", 
 143      ...                       out_handle, "fastq") 
 144      ... 
 145      >>> print("Converted %i records" % count) 
 146      Converted 10 records 
 147   
 148  The output FASTQ file would start like this: 
 149   
 150      >>> print("%s..." % out_handle.getvalue()[:50]) 
 151      @E3MFGYR02JWQ7T 
 152      tcagGGTCTACATGTTGGTTAACCCGTACTGATT... 
 153   
 154  Bio.SeqIO.index() provides memory efficient random access to the reads in an 
 155  SFF file by name. SFF files can include an index within the file, which can 
 156  be read in making this very fast. If the index is missing (or in a format not 
 157  yet supported in Biopython) the file is indexed by scanning all the reads - 
 158  which is a little slower. For example, 
 159   
 160      >>> from Bio import SeqIO 
 161      >>> reads = SeqIO.index("Roche/E3MFGYR02_random_10_reads.sff", "sff") 
 162      >>> record = reads["E3MFGYR02JHD4H"] 
 163      >>> print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
 164      E3MFGYR02JHD4H 310 tcagAAAGACAAGTGGTATC... 
 165      >>> reads.close() 
 166   
 167  Or, using the trimmed reads: 
 168   
 169      >>> from Bio import SeqIO 
 170      >>> reads = SeqIO.index("Roche/E3MFGYR02_random_10_reads.sff", "sff-trim") 
 171      >>> record = reads["E3MFGYR02JHD4H"] 
 172      >>> print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
 173      E3MFGYR02JHD4H 292 AAAGACAAGTGGTATCAACG... 
 174      >>> reads.close() 
 175   
 176  You can also use the Bio.SeqIO.write() function with the "sff" format. Note 
 177  that this requires all the flow information etc, and thus is probably only 
 178  useful for SeqRecord objects originally from reading another SFF file (and 
 179  not the trimmed SeqRecord objects from parsing an SFF file as "sff-trim"). 
 180   
 181  As an example, let's pretend this example SFF file represents some DNA which 
 182  was pre-amplified with a PCR primers AAAGANNNNN. The following script would 
 183  produce a sub-file containing all those reads whose post-quality clipping 
 184  region (i.e. the sequence after trimming) starts with AAAGA exactly (the non- 
 185  degenerate bit of this pretend primer): 
 186   
 187      >>> from Bio import SeqIO 
 188      >>> records = (record for record in 
 189      ...            SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff") 
 190      ...            if record.seq[record.annotations["clip_qual_left"]:].startswith("AAAGA")) 
 191      ... 
 192      >>> count = SeqIO.write(records, "temp_filtered.sff", "sff") 
 193      >>> print("Selected %i records" % count) 
 194      Selected 2 records 
 195   
 196  Of course, for an assembly you would probably want to remove these primers. 
 197  If you want FASTA or FASTQ output, you could just slice the SeqRecord. However, 
 198  if you want SFF output we have to preserve all the flow information - the trick 
 199  is just to adjust the left clip position! 
 200   
 201      >>> from Bio import SeqIO 
 202      >>> def filter_and_trim(records, primer): 
 203      ...     for record in records: 
 204      ...         if record.seq[record.annotations["clip_qual_left"]:].startswith(primer): 
 205      ...             record.annotations["clip_qual_left"] += len(primer) 
 206      ...             yield record 
 207      ... 
 208      >>> records = SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff") 
 209      >>> count = SeqIO.write(filter_and_trim(records, "AAAGA"), 
 210      ...                     "temp_filtered.sff", "sff") 
 211      ... 
 212      >>> print("Selected %i records" % count) 
 213      Selected 2 records 
 214   
 215  We can check the results, note the lower case clipped region now includes the "AAAGA" 
 216  sequence: 
 217   
 218      >>> for record in SeqIO.parse("temp_filtered.sff", "sff"): 
 219      ...     print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
 220      ... 
 221      E3MFGYR02JHD4H 310 tcagaaagaCAAGTGGTATC... 
 222      E3MFGYR02GAZMS 278 tcagaaagaAGTAAGGTAAA... 
 223      >>> for record in SeqIO.parse("temp_filtered.sff", "sff-trim"): 
 224      ...     print("%s %i %s..." % (record.id, len(record), record.seq[:20])) 
 225      ... 
 226      E3MFGYR02JHD4H 287 CAAGTGGTATCAACGCAGAG... 
 227      E3MFGYR02GAZMS 266 AGTAAGGTAAATAACAAACG... 
 228      >>> import os 
 229      >>> os.remove("temp_filtered.sff") 
 230   
 231  For a description of the file format, please see the Roche manuals and: 
 232  http://www.ncbi.nlm.nih.gov/Traces/trace.cgi?cmd=show&f=formats&m=doc&s=formats 
 233   
 234  """ 
 235   
 236  from __future__ import print_function 
 237   
 238  from Bio.SeqIO.Interfaces import SequenceWriter 
 239  from Bio import Alphabet 
 240  from Bio.Seq import Seq 
 241  from Bio.SeqRecord import SeqRecord 
 242  import struct 
 243  import sys 
 244  import re 
 245   
 246  from Bio._py3k import _bytes_to_string, _as_bytes 
 247   
 248  _null = b"\0" 
 249  _sff = b".sff" 
 250  _hsh = b".hsh" 
 251  _srt = b".srt" 
 252  _mft = b".mft" 
 253  _flag = b"\xff" 
 254   
 255   
256 -def _check_mode(handle):
257 """Ensure handle not opened in text mode. 258 259 Ensures mode is not set for Universal new line 260 and ensures mode is binary for Windows 261 """ 262 # TODO - Does this need to be stricter under Python 3? 263 mode = "" 264 if hasattr(handle, "mode"): 265 mode = handle.mode 266 if mode == 1: 267 # gzip.open(...) does this, fine 268 return 269 mode = str(mode) 270 271 if mode and "U" in mode.upper(): 272 raise ValueError("SFF files must NOT be opened in universal new " 273 "lines mode. Binary mode is recommended (although " 274 "on Unix the default mode is also fine).") 275 elif mode and "B" not in mode.upper() \ 276 and sys.platform == "win32": 277 raise ValueError("SFF files must be opened in binary mode on Windows")
278 279
280 -def _sff_file_header(handle):
281 """Read in an SFF file header (PRIVATE). 282 283 Assumes the handle is at the start of the file, will read forwards 284 though the header and leave the handle pointing at the first record. 285 Returns a tuple of values from the header (header_length, index_offset, 286 index_length, number_of_reads, flows_per_read, flow_chars, key_sequence) 287 288 >>> with open("Roche/greek.sff", "rb") as handle: 289 ... values = _sff_file_header(handle) 290 ... 291 >>> print(values[0]) 292 840 293 >>> print(values[1]) 294 65040 295 >>> print(values[2]) 296 256 297 >>> print(values[3]) 298 24 299 >>> print(values[4]) 300 800 301 >>> values[-1] 302 'TCAG' 303 304 """ 305 _check_mode(handle) 306 # file header (part one) 307 # use big endiean encdoing > 308 # magic_number I 309 # version 4B 310 # index_offset Q 311 # index_length I 312 # number_of_reads I 313 # header_length H 314 # key_length H 315 # number_of_flows_per_read H 316 # flowgram_format_code B 317 # [rest of file header depends on the number of flows and how many keys] 318 fmt = '>4s4BQIIHHHB' 319 assert 31 == struct.calcsize(fmt) 320 data = handle.read(31) 321 if not data: 322 raise ValueError("Empty file.") 323 elif len(data) < 31: 324 raise ValueError("File too small to hold a valid SFF header.") 325 magic_number, ver0, ver1, ver2, ver3, index_offset, index_length, \ 326 number_of_reads, header_length, key_length, number_of_flows_per_read, \ 327 flowgram_format = struct.unpack(fmt, data) 328 if magic_number in [_hsh, _srt, _mft]: 329 # Probably user error, calling Bio.SeqIO.parse() twice! 330 raise ValueError("Handle seems to be at SFF index block, not start") 331 if magic_number != _sff: # 779314790 332 raise ValueError("SFF file did not start '.sff', but %r" % magic_number) 333 if (ver0, ver1, ver2, ver3) != (0, 0, 0, 1): 334 raise ValueError("Unsupported SFF version in header, %i.%i.%i.%i" 335 % (ver0, ver1, ver2, ver3)) 336 if flowgram_format != 1: 337 raise ValueError("Flowgram format code %i not supported" 338 % flowgram_format) 339 if (index_offset != 0) ^ (index_length != 0): 340 raise ValueError("Index offset %i but index length %i" 341 % (index_offset, index_length)) 342 flow_chars = _bytes_to_string(handle.read(number_of_flows_per_read)) 343 key_sequence = _bytes_to_string(handle.read(key_length)) 344 # According to the spec, the header_length field should be the total number 345 # of bytes required by this set of header fields, and should be equal to 346 # "31 + number_of_flows_per_read + key_length" rounded up to the next value 347 # divisible by 8. 348 assert header_length % 8 == 0 349 padding = header_length - number_of_flows_per_read - key_length - 31 350 assert 0 <= padding < 8, padding 351 if handle.read(padding).count(_null) != padding: 352 import warnings 353 from Bio import BiopythonParserWarning 354 warnings.warn("Your SFF file is invalid, post header %i byte " 355 "null padding region contained data." % padding, 356 BiopythonParserWarning) 357 return header_length, index_offset, index_length, \ 358 number_of_reads, number_of_flows_per_read, \ 359 flow_chars, key_sequence
360 361
362 -def _sff_do_slow_index(handle):
363 """Generates an index by scanning though all the reads in an SFF file (PRIVATE). 364 365 This is a slow but generic approach if we can't parse the provided index 366 (if present). 367 368 Will use the handle seek/tell functions. 369 """ 370 handle.seek(0) 371 header_length, index_offset, index_length, number_of_reads, \ 372 number_of_flows_per_read, flow_chars, key_sequence \ 373 = _sff_file_header(handle) 374 # Now on to the reads... 375 read_header_fmt = '>2HI4H' 376 read_header_size = struct.calcsize(read_header_fmt) 377 # NOTE - assuming flowgram_format==1, which means struct type H 378 read_flow_fmt = ">%iH" % number_of_flows_per_read 379 read_flow_size = struct.calcsize(read_flow_fmt) 380 assert 1 == struct.calcsize(">B") 381 assert 1 == struct.calcsize(">s") 382 assert 1 == struct.calcsize(">c") 383 assert read_header_size % 8 == 0 # Important for padding calc later! 384 for read in range(number_of_reads): 385 record_offset = handle.tell() 386 if record_offset == index_offset: 387 # Found index block within reads, ignore it: 388 offset = index_offset + index_length 389 if offset % 8: 390 offset += 8 - (offset % 8) 391 assert offset % 8 == 0 392 handle.seek(offset) 393 record_offset = offset 394 # assert record_offset%8 == 0 # Worth checking, but slow 395 # First the fixed header 396 data = handle.read(read_header_size) 397 read_header_length, name_length, seq_len, clip_qual_left, \ 398 clip_qual_right, clip_adapter_left, clip_adapter_right \ 399 = struct.unpack(read_header_fmt, data) 400 if read_header_length < 10 or read_header_length % 8 != 0: 401 raise ValueError("Malformed read header, says length is %i:\n%r" 402 % (read_header_length, data)) 403 # now the name and any padding (remainder of header) 404 name = _bytes_to_string(handle.read(name_length)) 405 padding = read_header_length - read_header_size - name_length 406 if handle.read(padding).count(_null) != padding: 407 import warnings 408 from Bio import BiopythonParserWarning 409 warnings.warn("Your SFF file is invalid, post name %i byte " 410 "padding region contained data" % padding, 411 BiopythonParserWarning) 412 assert record_offset + read_header_length == handle.tell() 413 # now the flowgram values, flowgram index, bases and qualities 414 size = read_flow_size + 3 * seq_len 415 handle.seek(size, 1) 416 # now any padding... 417 padding = size % 8 418 if padding: 419 padding = 8 - padding 420 if handle.read(padding).count(_null) != padding: 421 import warnings 422 from Bio import BiopythonParserWarning 423 warnings.warn("Your SFF file is invalid, post quality %i " 424 "byte padding region contained data" % padding, 425 BiopythonParserWarning) 426 # print("%s %s %i" % (read, name, record_offset)) 427 yield name, record_offset 428 if handle.tell() % 8 != 0: 429 raise ValueError( 430 "After scanning reads, did not end on a multiple of 8")
431 432
433 -def _sff_find_roche_index(handle):
434 """Locate any existing Roche style XML meta data and read index (PRIVATE). 435 436 Makes a number of hard coded assumptions based on reverse engineered SFF 437 files from Roche 454 machines. 438 439 Returns a tuple of read count, SFF "index" offset and size, XML offset 440 and size, and the actual read index offset and size. 441 442 Raises a ValueError for unsupported or non-Roche index blocks. 443 """ 444 handle.seek(0) 445 header_length, index_offset, index_length, number_of_reads, \ 446 number_of_flows_per_read, flow_chars, key_sequence \ 447 = _sff_file_header(handle) 448 assert handle.tell() == header_length 449 if not index_offset or not index_offset: 450 raise ValueError("No index present in this SFF file") 451 # Now jump to the header... 452 handle.seek(index_offset) 453 fmt = ">4s4B" 454 fmt_size = struct.calcsize(fmt) 455 data = handle.read(fmt_size) 456 if not data: 457 raise ValueError("Premature end of file? Expected index of size %i at offest %i, found nothing" 458 % (index_length, index_offset)) 459 if len(data) < fmt_size: 460 raise ValueError("Premature end of file? Expected index of size %i at offest %i, found %r" 461 % (index_length, index_offset, data)) 462 magic_number, ver0, ver1, ver2, ver3 = struct.unpack(fmt, data) 463 if magic_number == _mft: # 778921588 464 # Roche 454 manifest index 465 # This is typical from raw Roche 454 SFF files (2009), and includes 466 # both an XML manifest and the sorted index. 467 if (ver0, ver1, ver2, ver3) != (49, 46, 48, 48): 468 # This is "1.00" as a string 469 raise ValueError("Unsupported version in .mft index header, %i.%i.%i.%i" 470 % (ver0, ver1, ver2, ver3)) 471 fmt2 = ">LL" 472 fmt2_size = struct.calcsize(fmt2) 473 xml_size, data_size = struct.unpack(fmt2, handle.read(fmt2_size)) 474 if index_length != fmt_size + fmt2_size + xml_size + data_size: 475 raise ValueError("Problem understanding .mft index header, %i != %i + %i + %i + %i" 476 % (index_length, fmt_size, fmt2_size, xml_size, data_size)) 477 return number_of_reads, header_length, \ 478 index_offset, index_length, \ 479 index_offset + fmt_size + fmt2_size, xml_size, \ 480 index_offset + fmt_size + fmt2_size + xml_size, data_size 481 elif magic_number == _srt: # 779317876 482 # Roche 454 sorted index 483 # I've had this from Roche tool sfffile when the read identifiers 484 # had nonstandard lengths and there was no XML manifest. 485 if (ver0, ver1, ver2, ver3) != (49, 46, 48, 48): 486 # This is "1.00" as a string 487 raise ValueError("Unsupported version in .srt index header, %i.%i.%i.%i" 488 % (ver0, ver1, ver2, ver3)) 489 data = handle.read(4) 490 if data != _null * 4: 491 raise ValueError( 492 "Did not find expected null four bytes in .srt index") 493 return number_of_reads, header_length, \ 494 index_offset, index_length, \ 495 0, 0, \ 496 index_offset + fmt_size + 4, index_length - fmt_size - 4 497 elif magic_number == _hsh: 498 raise ValueError("Hash table style indexes (.hsh) in SFF files are " 499 "not (yet) supported") 500 else: 501 raise ValueError("Unknown magic number %r in SFF index header:\n%r" 502 % (magic_number, data))
503 504
505 -def ReadRocheXmlManifest(handle):
506 """Reads any Roche style XML manifest data in the SFF "index". 507 508 The SFF file format allows for multiple different index blocks, and Roche 509 took advantage of this to define their own index block which also embeds 510 an XML manifest string. This is not a publicly documented extension to 511 the SFF file format, this was reverse engineered. 512 513 The handle should be to an SFF file opened in binary mode. This function 514 will use the handle seek/tell functions and leave the handle in an 515 arbitrary location. 516 517 Any XML manifest found is returned as a Python string, which you can then 518 parse as appropriate, or reuse when writing out SFF files with the 519 SffWriter class. 520 521 Returns a string, or raises a ValueError if an Roche manifest could not be 522 found. 523 """ 524 number_of_reads, header_length, index_offset, index_length, xml_offset, \ 525 xml_size, read_index_offset, read_index_size = _sff_find_roche_index( 526 handle) 527 if not xml_offset or not xml_size: 528 raise ValueError("No XML manifest found") 529 handle.seek(xml_offset) 530 return _bytes_to_string(handle.read(xml_size))
531 532 533 # This is a generator function!
534 -def _sff_read_roche_index(handle):
535 """Reads any existing Roche style read index provided in the SFF file (PRIVATE). 536 537 Will use the handle seek/tell functions. 538 539 This works on ".srt1.00" and ".mft1.00" style Roche SFF index blocks. 540 541 Roche SFF indices use base 255 not 256, meaning we see bytes in range the 542 range 0 to 254 only. This appears to be so that byte 0xFF (character 255) 543 can be used as a marker character to separate entries (required if the 544 read name lengths vary). 545 546 Note that since only four bytes are used for the read offset, this is 547 limited to 255^4 bytes (nearly 4GB). If you try to use the Roche sfffile 548 tool to combine SFF files beyound this limit, they issue a warning and 549 omit the index (and manifest). 550 """ 551 number_of_reads, header_length, index_offset, index_length, xml_offset, \ 552 xml_size, read_index_offset, read_index_size = _sff_find_roche_index( 553 handle) 554 # Now parse the read index... 555 handle.seek(read_index_offset) 556 fmt = ">5B" 557 for read in range(number_of_reads): 558 # TODO - Be more aware of when the index should end? 559 data = handle.read(6) 560 while True: 561 more = handle.read(1) 562 if not more: 563 raise ValueError("Premature end of file!") 564 data += more 565 if more == _flag: 566 break 567 assert data[-1:] == _flag, data[-1:] 568 name = _bytes_to_string(data[:-6]) 569 off4, off3, off2, off1, off0 = struct.unpack(fmt, data[-6:-1]) 570 offset = off0 + 255 * off1 + 65025 * off2 + 16581375 * off3 571 if off4: 572 # Could in theory be used as a fifth piece of offset information, 573 # i.e. offset =+ 4228250625L*off4, but testing the Roche tools this 574 # is not the case. They simple don't support such large indexes. 575 raise ValueError("Expected a null terminator to the read name.") 576 yield name, offset 577 if handle.tell() != read_index_offset + read_index_size: 578 raise ValueError("Problem with index length? %i vs %i" 579 % (handle.tell(), read_index_offset + read_index_size))
580 581 582 _valid_UAN_read_name = re.compile(r'^[a-zA-Z0-9]{14}$') 583 584
585 -def _sff_read_seq_record(handle, number_of_flows_per_read, flow_chars, 586 key_sequence, alphabet, trim=False):
587 """Parse the next read in the file, return data as a SeqRecord (PRIVATE).""" 588 # Now on to the reads... 589 # the read header format (fixed part): 590 # read_header_length H 591 # name_length H 592 # seq_len I 593 # clip_qual_left H 594 # clip_qual_right H 595 # clip_adapter_left H 596 # clip_adapter_right H 597 # [rest of read header depends on the name length etc] 598 read_header_fmt = '>2HI4H' 599 read_header_size = struct.calcsize(read_header_fmt) 600 read_flow_fmt = ">%iH" % number_of_flows_per_read 601 read_flow_size = struct.calcsize(read_flow_fmt) 602 603 read_header_length, name_length, seq_len, clip_qual_left, \ 604 clip_qual_right, clip_adapter_left, clip_adapter_right \ 605 = struct.unpack(read_header_fmt, handle.read(read_header_size)) 606 if clip_qual_left: 607 clip_qual_left -= 1 # python counting 608 if clip_adapter_left: 609 clip_adapter_left -= 1 # python counting 610 if read_header_length < 10 or read_header_length % 8 != 0: 611 raise ValueError("Malformed read header, says length is %i" 612 % read_header_length) 613 # now the name and any padding (remainder of header) 614 name = _bytes_to_string(handle.read(name_length)) 615 padding = read_header_length - read_header_size - name_length 616 if handle.read(padding).count(_null) != padding: 617 import warnings 618 from Bio import BiopythonParserWarning 619 warnings.warn("Your SFF file is invalid, post name %i " 620 "byte padding region contained data" % padding, 621 BiopythonParserWarning) 622 # now the flowgram values, flowgram index, bases and qualities 623 # NOTE - assuming flowgram_format==1, which means struct type H 624 flow_values = handle.read(read_flow_size) # unpack later if needed 625 temp_fmt = ">%iB" % seq_len # used for flow index and quals 626 flow_index = handle.read(seq_len) # unpack later if needed 627 seq = _bytes_to_string(handle.read(seq_len)) # TODO - Use bytes in Seq? 628 quals = list(struct.unpack(temp_fmt, handle.read(seq_len))) 629 # now any padding... 630 padding = (read_flow_size + seq_len * 3) % 8 631 if padding: 632 padding = 8 - padding 633 if handle.read(padding).count(_null) != padding: 634 import warnings 635 from Bio import BiopythonParserWarning 636 warnings.warn("Your SFF file is invalid, post quality %i " 637 "byte padding region contained data" % padding, 638 BiopythonParserWarning) 639 # Follow Roche and apply most aggressive of qual and adapter clipping. 640 # Note Roche seems to ignore adapter clip fields when writing SFF, 641 # and uses just the quality clipping values for any clipping. 642 clip_left = max(clip_qual_left, clip_adapter_left) 643 # Right clipping of zero means no clipping 644 if clip_qual_right: 645 if clip_adapter_right: 646 clip_right = min(clip_qual_right, clip_adapter_right) 647 else: 648 # Typical case with Roche SFF files 649 clip_right = clip_qual_right 650 elif clip_adapter_right: 651 clip_right = clip_adapter_right 652 else: 653 clip_right = seq_len 654 # Now build a SeqRecord 655 if trim: 656 if clip_left >= clip_right: 657 # Raise an error? 658 import warnings 659 from Bio import BiopythonParserWarning 660 warnings.warn("Overlapping clip values in SFF record, trimmed to nothing", 661 BiopythonParserWarning) 662 seq = "" 663 quals = [] 664 else: 665 seq = seq[clip_left:clip_right].upper() 666 quals = quals[clip_left:clip_right] 667 # Don't record the clipping values, flow etc, they make no sense now: 668 annotations = {} 669 else: 670 if clip_left >= clip_right: 671 import warnings 672 from Bio import BiopythonParserWarning 673 warnings.warn("Overlapping clip values in SFF record", BiopythonParserWarning) 674 seq = seq.lower() 675 else: 676 # This use of mixed case mimics the Roche SFF tool's FASTA output 677 seq = seq[:clip_left].lower() + \ 678 seq[clip_left:clip_right].upper() + \ 679 seq[clip_right:].lower() 680 annotations = {"flow_values": struct.unpack(read_flow_fmt, flow_values), 681 "flow_index": struct.unpack(temp_fmt, flow_index), 682 "flow_chars": flow_chars, 683 "flow_key": key_sequence, 684 "clip_qual_left": clip_qual_left, 685 "clip_qual_right": clip_qual_right, 686 "clip_adapter_left": clip_adapter_left, 687 "clip_adapter_right": clip_adapter_right} 688 if re.match(_valid_UAN_read_name, name): 689 annotations["time"] = _get_read_time(name) 690 annotations["region"] = _get_read_region(name) 691 annotations["coords"] = _get_read_xy(name) 692 record = SeqRecord(Seq(seq, alphabet), 693 id=name, 694 name=name, 695 description="", 696 annotations=annotations) 697 # Dirty trick to speed up this line: 698 # record.letter_annotations["phred_quality"] = quals 699 dict.__setitem__(record._per_letter_annotations, 700 "phred_quality", quals) 701 # Return the record and then continue... 702 return record
703 704 705 _powers_of_36 = [36 ** i for i in range(6)] 706 707
708 -def _string_as_base_36(string):
709 """Interpret a string as a base-36 number as per 454 manual.""" 710 total = 0 711 for c, power in zip(string[::-1], _powers_of_36): 712 # For reference: ord('0') = 48, ord('9') = 57 713 # For reference: ord('A') = 65, ord('Z') = 90 714 # For reference: ord('a') = 97, ord('z') = 122 715 if 48 <= ord(c) <= 57: 716 val = ord(c) - 22 # equivalent to: - ord('0') + 26 717 elif 65 <= ord(c) <= 90: 718 val = ord(c) - 65 719 elif 97 <= ord(c) <= 122: 720 val = ord(c) - 97 721 else: 722 # Invalid character 723 val = 0 724 total += val * power 725 return total
726 727
728 -def _get_read_xy(read_name):
729 """Extract coordinates from last 5 characters of read name.""" 730 number = _string_as_base_36(read_name[9:]) 731 return divmod(number, 4096)
732 733 734 _time_denominators = [13 * 32 * 24 * 60 * 60, 735 32 * 24 * 60 * 60, 736 24 * 60 * 60, 737 60 * 60, 738 60] 739 740
741 -def _get_read_time(read_name):
742 """Extract time from first 6 characters of read name.""" 743 time_list = [] 744 remainder = _string_as_base_36(read_name[:6]) 745 for denominator in _time_denominators: 746 this_term, remainder = divmod(remainder, denominator) 747 time_list.append(this_term) 748 time_list.append(remainder) 749 time_list[0] += 2000 750 return time_list
751 752
753 -def _get_read_region(read_name):
754 """Extract region from read name.""" 755 return int(read_name[8])
756 757
758 -def _sff_read_raw_record(handle, number_of_flows_per_read):
759 """Extract the next read in the file as a raw (bytes) string (PRIVATE).""" 760 read_header_fmt = '>2HI' 761 read_header_size = struct.calcsize(read_header_fmt) 762 read_flow_fmt = ">%iH" % number_of_flows_per_read 763 read_flow_size = struct.calcsize(read_flow_fmt) 764 765 raw = handle.read(read_header_size) 766 read_header_length, name_length, seq_len \ 767 = struct.unpack(read_header_fmt, raw) 768 if read_header_length < 10 or read_header_length % 8 != 0: 769 raise ValueError("Malformed read header, says length is %i" 770 % read_header_length) 771 # now the four clip values (4H = 8 bytes), and read name 772 raw += handle.read(8 + name_length) 773 # and any padding (remainder of header) 774 padding = read_header_length - read_header_size - 8 - name_length 775 pad = handle.read(padding) 776 if pad.count(_null) != padding: 777 import warnings 778 from Bio import BiopythonParserWarning 779 warnings.warn("Your SFF file is invalid, post name %i " 780 "byte padding region contained data" % padding, 781 BiopythonParserWarning) 782 raw += pad 783 # now the flowgram values, flowgram index, bases and qualities 784 raw += handle.read(read_flow_size + seq_len * 3) 785 padding = (read_flow_size + seq_len * 3) % 8 786 # now any padding... 787 if padding: 788 padding = 8 - padding 789 pad = handle.read(padding) 790 if pad.count(_null) != padding: 791 import warnings 792 from Bio import BiopythonParserWarning 793 warnings.warn("Your SFF file is invalid, post quality %i " 794 "byte padding region contained data" % padding, 795 BiopythonParserWarning) 796 raw += pad 797 # Return the raw bytes 798 return raw
799 800
801 -class _AddTellHandle(object):
802 """Wrapper for handles which do not support the tell method (PRIVATE). 803 804 Intended for use with things like network handles where tell (and reverse 805 seek) are not supported. The SFF file needs to track the current offset in 806 order to deal with the index block. 807 """ 808
809 - def __init__(self, handle):
810 self._handle = handle 811 self._offset = 0
812
813 - def read(self, length):
814 data = self._handle.read(length) 815 self._offset += len(data) 816 return data
817
818 - def tell(self):
819 return self._offset
820
821 - def seek(self, offset):
822 if offset < self._offset: 823 raise RuntimeError("Can't seek backwards") 824 self._handle.read(offset - self._offset)
825
826 - def close(self):
827 return self._handle.close()
828 829 830 # This is a generator function!
831 -def SffIterator(handle, alphabet=Alphabet.generic_dna, trim=False):
832 """Iterate over Standard Flowgram Format (SFF) reads (as SeqRecord objects). 833 834 - handle - input file, an SFF file, e.g. from Roche 454 sequencing. 835 This must NOT be opened in universal read lines mode! 836 - alphabet - optional alphabet, defaults to generic DNA. 837 - trim - should the sequences be trimmed? 838 839 The resulting SeqRecord objects should match those from a paired FASTA 840 and QUAL file converted from the SFF file using the Roche 454 tool 841 ssfinfo. i.e. The sequence will be mixed case, with the trim regions 842 shown in lower case. 843 844 This function is used internally via the Bio.SeqIO functions: 845 846 >>> from Bio import SeqIO 847 >>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff"): 848 ... print("%s %i" % (record.id, len(record))) 849 ... 850 E3MFGYR02JWQ7T 265 851 E3MFGYR02JA6IL 271 852 E3MFGYR02JHD4H 310 853 E3MFGYR02GFKUC 299 854 E3MFGYR02FTGED 281 855 E3MFGYR02FR9G7 261 856 E3MFGYR02GAZMS 278 857 E3MFGYR02HHZ8O 221 858 E3MFGYR02GPGB1 269 859 E3MFGYR02F7Z7G 219 860 861 You can also call it directly: 862 863 >>> with open("Roche/E3MFGYR02_random_10_reads.sff", "rb") as handle: 864 ... for record in SffIterator(handle): 865 ... print("%s %i" % (record.id, len(record))) 866 ... 867 E3MFGYR02JWQ7T 265 868 E3MFGYR02JA6IL 271 869 E3MFGYR02JHD4H 310 870 E3MFGYR02GFKUC 299 871 E3MFGYR02FTGED 281 872 E3MFGYR02FR9G7 261 873 E3MFGYR02GAZMS 278 874 E3MFGYR02HHZ8O 221 875 E3MFGYR02GPGB1 269 876 E3MFGYR02F7Z7G 219 877 878 Or, with the trim option: 879 880 >>> with open("Roche/E3MFGYR02_random_10_reads.sff", "rb") as handle: 881 ... for record in SffIterator(handle, trim=True): 882 ... print("%s %i" % (record.id, len(record))) 883 ... 884 E3MFGYR02JWQ7T 260 885 E3MFGYR02JA6IL 265 886 E3MFGYR02JHD4H 292 887 E3MFGYR02GFKUC 295 888 E3MFGYR02FTGED 277 889 E3MFGYR02FR9G7 256 890 E3MFGYR02GAZMS 271 891 E3MFGYR02HHZ8O 150 892 E3MFGYR02GPGB1 221 893 E3MFGYR02F7Z7G 130 894 895 """ 896 if isinstance(Alphabet._get_base_alphabet(alphabet), 897 Alphabet.ProteinAlphabet): 898 raise ValueError("Invalid alphabet, SFF files do not hold proteins.") 899 if isinstance(Alphabet._get_base_alphabet(alphabet), 900 Alphabet.RNAAlphabet): 901 raise ValueError("Invalid alphabet, SFF files do not hold RNA.") 902 try: 903 if 0 != handle.tell(): 904 raise ValueError("Not at start of file, offset %i" % handle.tell()) 905 except AttributeError: 906 # Probably a network handle or something like that 907 handle = _AddTellHandle(handle) 908 header_length, index_offset, index_length, number_of_reads, \ 909 number_of_flows_per_read, flow_chars, key_sequence \ 910 = _sff_file_header(handle) 911 # Now on to the reads... 912 # the read header format (fixed part): 913 # read_header_length H 914 # name_length H 915 # seq_len I 916 # clip_qual_left H 917 # clip_qual_right H 918 # clip_adapter_left H 919 # clip_adapter_right H 920 # [rest of read header depends on the name length etc] 921 read_header_fmt = '>2HI4H' 922 read_header_size = struct.calcsize(read_header_fmt) 923 read_flow_fmt = ">%iH" % number_of_flows_per_read 924 read_flow_size = struct.calcsize(read_flow_fmt) 925 assert 1 == struct.calcsize(">B") 926 assert 1 == struct.calcsize(">s") 927 assert 1 == struct.calcsize(">c") 928 assert read_header_size % 8 == 0 # Important for padding calc later! 929 # The spec allows for the index block to be before or even in the middle 930 # of the reads. We can check that if we keep track of our position 931 # in the file... 932 for read in range(number_of_reads): 933 if index_offset and handle.tell() == index_offset: 934 offset = index_offset + index_length 935 if offset % 8: 936 offset += 8 - (offset % 8) 937 assert offset % 8 == 0 938 handle.seek(offset) 939 # Now that we've done this, we don't need to do it again. Clear 940 # the index_offset so we can skip extra handle.tell() calls: 941 index_offset = 0 942 yield _sff_read_seq_record(handle, 943 number_of_flows_per_read, 944 flow_chars, 945 key_sequence, 946 alphabet, 947 trim) 948 _check_eof(handle, index_offset, index_length)
949 950
951 -def _check_eof(handle, index_offset, index_length):
952 """Check final padding is OK (8 byte alignment) and file ends (PRIVATE). 953 954 Will attempt to spot apparent SFF file concatenation and give an error. 955 956 Will not attempt to seek, only moves the handle forward. 957 """ 958 offset = handle.tell() 959 extra = b"" 960 padding = 0 961 962 if index_offset and offset <= index_offset: 963 # Index block then end of file... 964 if offset < index_offset: 965 raise ValueError("Gap of %i bytes after final record end %i, " 966 "before %i where index starts?" 967 % (index_offset - offset, offset, index_offset)) 968 # Doing read to jump the index rather than a seek 969 # in case this is a network handle or similar 970 handle.read(index_offset + index_length - offset) 971 offset = index_offset + index_length 972 assert offset == handle.tell(), \ 973 "Wanted %i, got %i, index is %i to %i" \ 974 % (offset, handle.tell(), index_offset, index_offset + index_length) 975 976 if offset % 8: 977 padding = 8 - (offset % 8) 978 extra = handle.read(padding) 979 980 if padding >= 4 and extra[-4:] == _sff: 981 # Seen this in one user supplied file, should have been 982 # four bytes of null padding but was actually .sff and 983 # the start of a new concatenated SFF file! 984 raise ValueError("Your SFF file is invalid, post index %i byte " 985 "null padding region ended '.sff' which could " 986 "be the start of a concatenated SFF file? " 987 "See offset %i" % (padding, offset)) 988 if padding and not extra: 989 # TODO - Is this error harmless enough to just ignore? 990 import warnings 991 from Bio import BiopythonParserWarning 992 warnings.warn("Your SFF file is technically invalid as it is missing " 993 "a terminal %i byte null padding region." % padding, 994 BiopythonParserWarning) 995 return 996 if extra.count(_null) != padding: 997 import warnings 998 from Bio import BiopythonParserWarning 999 warnings.warn("Your SFF file is invalid, post index %i byte " 1000 "null padding region contained data: %r" 1001 % (padding, extra), BiopythonParserWarning) 1002 1003 offset = handle.tell() 1004 assert offset % 8 == 0, \ 1005 "Wanted offset %i %% 8 = %i to be zero" % (offset, offset % 8) 1006 # Should now be at the end of the file... 1007 extra = handle.read(4) 1008 if extra == _sff: 1009 raise ValueError("Additional data at end of SFF file, " 1010 "perhaps multiple SFF files concatenated? " 1011 "See offset %i" % offset) 1012 elif extra: 1013 raise ValueError("Additional data at end of SFF file, " 1014 "see offset %i" % offset)
1015 1016 1017 # This is a generator function!
1018 -def _SffTrimIterator(handle, alphabet=Alphabet.generic_dna):
1019 """Iterate over SFF reads (as SeqRecord objects) with trimming (PRIVATE).""" 1020 return SffIterator(handle, alphabet, trim=True)
1021 1022
1023 -class SffWriter(SequenceWriter):
1024 """SFF file writer.""" 1025
1026 - def __init__(self, handle, index=True, xml=None):
1027 """Creates the writer object. 1028 1029 - handle - Output handle, ideally in binary write mode. 1030 - index - Boolean argument, should we try and write an index? 1031 - xml - Optional string argument, xml manifest to be recorded in the index 1032 block (see function ReadRocheXmlManifest for reading this data). 1033 """ 1034 _check_mode(handle) 1035 self.handle = handle 1036 self._xml = xml 1037 if index: 1038 self._index = [] 1039 else: 1040 self._index = None
1041
1042 - def write_file(self, records):
1043 """Use this to write an entire file containing the given records.""" 1044 try: 1045 self._number_of_reads = len(records) 1046 except TypeError: 1047 self._number_of_reads = 0 # dummy value 1048 if not hasattr(self.handle, "seek") \ 1049 or not hasattr(self.handle, "tell"): 1050 raise ValueError("A handle with a seek/tell methods is " 1051 "required in order to record the total " 1052 "record count in the file header (once it " 1053 "is known at the end).") 1054 if self._index is not None and \ 1055 not (hasattr(self.handle, "seek") and hasattr(self.handle, "tell")): 1056 import warnings 1057 warnings.warn("A handle with a seek/tell methods is required in " 1058 "order to record an SFF index.") 1059 self._index = None 1060 self._index_start = 0 1061 self._index_length = 0 1062 if not hasattr(records, "next"): 1063 records = iter(records) 1064 # Get the first record in order to find the flow information 1065 # we will need for the header. 1066 try: 1067 record = next(records) 1068 except StopIteration: 1069 record = None 1070 if record is None: 1071 # No records -> empty SFF file (or an error)? 1072 # We can't write a header without the flow information. 1073 # return 0 1074 raise ValueError("Must have at least one sequence") 1075 try: 1076 self._key_sequence = _as_bytes(record.annotations["flow_key"]) 1077 self._flow_chars = _as_bytes(record.annotations["flow_chars"]) 1078 self._number_of_flows_per_read = len(self._flow_chars) 1079 except KeyError: 1080 raise ValueError("Missing SFF flow information") 1081 self.write_header() 1082 self.write_record(record) 1083 count = 1 1084 for record in records: 1085 self.write_record(record) 1086 count += 1 1087 if self._number_of_reads == 0: 1088 # Must go back and record the record count... 1089 offset = self.handle.tell() 1090 self.handle.seek(0) 1091 self._number_of_reads = count 1092 self.write_header() 1093 self.handle.seek(offset) # not essential? 1094 else: 1095 assert count == self._number_of_reads 1096 if self._index is not None: 1097 self._write_index() 1098 return count
1099
1100 - def _write_index(self):
1101 assert len(self._index) == self._number_of_reads 1102 handle = self.handle 1103 self._index.sort() 1104 self._index_start = handle.tell() # need for header 1105 # XML... 1106 if self._xml is not None: 1107 xml = _as_bytes(self._xml) 1108 else: 1109 from Bio import __version__ 1110 xml = "<!-- This file was output with Biopython %s -->\n" % __version__ 1111 xml += "<!-- This XML and index block attempts to mimic Roche SFF files -->\n" 1112 xml += "<!-- This file may be a combination of multiple SFF files etc -->\n" 1113 xml = _as_bytes(xml) 1114 xml_len = len(xml) 1115 # Write to the file... 1116 fmt = ">I4BLL" 1117 fmt_size = struct.calcsize(fmt) 1118 handle.write(_null * fmt_size + xml) # fill this later 1119 fmt2 = ">6B" 1120 assert 6 == struct.calcsize(fmt2) 1121 self._index.sort() 1122 index_len = 0 # don't know yet! 1123 for name, offset in self._index: 1124 # Roche files record the offsets using base 255 not 256. 1125 # See comments for parsing the index block. There may be a faster 1126 # way to code this, but we can't easily use shifts due to odd base 1127 off3 = offset 1128 off0 = off3 % 255 1129 off3 -= off0 1130 off1 = off3 % 65025 1131 off3 -= off1 1132 off2 = off3 % 16581375 1133 off3 -= off2 1134 assert offset == off0 + off1 + off2 + off3, \ 1135 "%i -> %i %i %i %i" % (offset, off0, off1, off2, off3) 1136 off3, off2, off1, off0 = off3 // 16581375, off2 // 65025, \ 1137 off1 // 255, off0 1138 assert off0 < 255 and off1 < 255 and off2 < 255 and off3 < 255, \ 1139 "%i -> %i %i %i %i" % (offset, off0, off1, off2, off3) 1140 handle.write(name + struct.pack(fmt2, 0, 1141 off3, off2, off1, off0, 255)) 1142 index_len += len(name) + 6 1143 # Note any padding in not included: 1144 self._index_length = fmt_size + xml_len + index_len # need for header 1145 # Pad out to an 8 byte boundary (although I have noticed some 1146 # real Roche SFF files neglect to do this depsite their manual 1147 # suggesting this padding should be there): 1148 if self._index_length % 8: 1149 padding = 8 - (self._index_length % 8) 1150 handle.write(_null * padding) 1151 else: 1152 padding = 0 1153 offset = handle.tell() 1154 assert offset == self._index_start + self._index_length + padding, \ 1155 "%i vs %i + %i + %i" % (offset, self._index_start, 1156 self._index_length, padding) 1157 # Must now go back and update the index header with index size... 1158 handle.seek(self._index_start) 1159 handle.write(struct.pack(fmt, 778921588, # magic number 1160 49, 46, 48, 48, # Roche index version, "1.00" 1161 xml_len, index_len) + xml) 1162 # Must now go back and update the header... 1163 handle.seek(0) 1164 self.write_header() 1165 handle.seek(offset) # not essential?
1166
1167 - def write_header(self):
1168 # Do header... 1169 key_length = len(self._key_sequence) 1170 # file header (part one) 1171 # use big endiean encdoing > 1172 # magic_number I 1173 # version 4B 1174 # index_offset Q 1175 # index_length I 1176 # number_of_reads I 1177 # header_length H 1178 # key_length H 1179 # number_of_flows_per_read H 1180 # flowgram_format_code B 1181 # [rest of file header depends on the number of flows and how many keys] 1182 fmt = '>I4BQIIHHHB%is%is' % ( 1183 self._number_of_flows_per_read, key_length) 1184 # According to the spec, the header_length field should be the total 1185 # number of bytes required by this set of header fields, and should be 1186 # equal to "31 + number_of_flows_per_read + key_length" rounded up to 1187 # the next value divisible by 8. 1188 if struct.calcsize(fmt) % 8 == 0: 1189 padding = 0 1190 else: 1191 padding = 8 - (struct.calcsize(fmt) % 8) 1192 header_length = struct.calcsize(fmt) + padding 1193 assert header_length % 8 == 0 1194 header = struct.pack(fmt, 779314790, # magic number 0x2E736666 1195 0, 0, 0, 1, # version 1196 self._index_start, self._index_length, 1197 self._number_of_reads, 1198 header_length, key_length, 1199 self._number_of_flows_per_read, 1200 1, # the only flowgram format code we support 1201 self._flow_chars, self._key_sequence) 1202 self.handle.write(header + _null * padding)
1203
1204 - def write_record(self, record):
1205 """Write a single additional record to the output file. 1206 1207 This assumes the header has been done. 1208 """ 1209 # Basics 1210 name = _as_bytes(record.id) 1211 name_len = len(name) 1212 seq = _as_bytes(str(record.seq).upper()) 1213 seq_len = len(seq) 1214 # Qualities 1215 try: 1216 quals = record.letter_annotations["phred_quality"] 1217 except KeyError: 1218 raise ValueError("Missing PHRED qualities information for %s" % record.id) 1219 # Flow 1220 try: 1221 flow_values = record.annotations["flow_values"] 1222 flow_index = record.annotations["flow_index"] 1223 if self._key_sequence != _as_bytes(record.annotations["flow_key"]) \ 1224 or self._flow_chars != _as_bytes(record.annotations["flow_chars"]): 1225 raise ValueError("Records have inconsistent SFF flow data") 1226 except KeyError: 1227 raise ValueError("Missing SFF flow information for %s" % record.id) 1228 except AttributeError: 1229 raise ValueError("Header not written yet?") 1230 # Clipping 1231 try: 1232 clip_qual_left = record.annotations["clip_qual_left"] 1233 if clip_qual_left < 0: 1234 raise ValueError("Negative SFF clip_qual_left value for %s" % record.id) 1235 if clip_qual_left: 1236 clip_qual_left += 1 1237 clip_qual_right = record.annotations["clip_qual_right"] 1238 if clip_qual_right < 0: 1239 raise ValueError("Negative SFF clip_qual_right value for %s" % record.id) 1240 clip_adapter_left = record.annotations["clip_adapter_left"] 1241 if clip_adapter_left < 0: 1242 raise ValueError("Negative SFF clip_adapter_left value for %s" % record.id) 1243 if clip_adapter_left: 1244 clip_adapter_left += 1 1245 clip_adapter_right = record.annotations["clip_adapter_right"] 1246 if clip_adapter_right < 0: 1247 raise ValueError("Negative SFF clip_adapter_right value for %s" % record.id) 1248 except KeyError: 1249 raise ValueError("Missing SFF clipping information for %s" % record.id) 1250 1251 # Capture information for index 1252 if self._index is not None: 1253 offset = self.handle.tell() 1254 # Check the position of the final record (before sort by name) 1255 # Using a four-digit base 255 number, so the upper bound is 1256 # 254*(1)+254*(255)+254*(255**2)+254*(255**3) = 4228250624 1257 # or equivalently it overflows at 255**4 = 4228250625 1258 if offset > 4228250624: 1259 import warnings 1260 warnings.warn("Read %s has file offset %i, which is too large " 1261 "to store in the Roche SFF index structure. No " 1262 "index block will be recorded." % (name, offset)) 1263 # No point recoring the offsets now 1264 self._index = None 1265 else: 1266 self._index.append((name, self.handle.tell())) 1267 1268 # the read header format (fixed part): 1269 # read_header_length H 1270 # name_length H 1271 # seq_len I 1272 # clip_qual_left H 1273 # clip_qual_right H 1274 # clip_adapter_left H 1275 # clip_adapter_right H 1276 # [rest of read header depends on the name length etc] 1277 # name 1278 # flow values 1279 # flow index 1280 # sequence 1281 # padding 1282 read_header_fmt = '>2HI4H%is' % name_len 1283 if struct.calcsize(read_header_fmt) % 8 == 0: 1284 padding = 0 1285 else: 1286 padding = 8 - (struct.calcsize(read_header_fmt) % 8) 1287 read_header_length = struct.calcsize(read_header_fmt) + padding 1288 assert read_header_length % 8 == 0 1289 data = struct.pack(read_header_fmt, 1290 read_header_length, 1291 name_len, seq_len, 1292 clip_qual_left, clip_qual_right, 1293 clip_adapter_left, clip_adapter_right, 1294 name) + _null * padding 1295 assert len(data) == read_header_length 1296 # now the flowgram values, flowgram index, bases and qualities 1297 # NOTE - assuming flowgram_format==1, which means struct type H 1298 read_flow_fmt = ">%iH" % self._number_of_flows_per_read 1299 read_flow_size = struct.calcsize(read_flow_fmt) 1300 temp_fmt = ">%iB" % seq_len # used for flow index and quals 1301 data += struct.pack(read_flow_fmt, *flow_values) \ 1302 + struct.pack(temp_fmt, *flow_index) \ 1303 + seq \ 1304 + struct.pack(temp_fmt, *quals) 1305 # now any final padding... 1306 padding = (read_flow_size + seq_len * 3) % 8 1307 if padding: 1308 padding = 8 - padding 1309 self.handle.write(data + _null * padding)
1310 1311 1312 if __name__ == "__main__": 1313 from Bio._utils import run_doctest 1314 1315 run_doctest(verbose=0) 1316