Bio.SeqIO.SffIO module

Bio.SeqIO support for the binary Standard Flowgram Format (SFF) file format.

SFF was designed by 454 Life Sciences (Roche), the Whitehead Institute for Biomedical Research and the Wellcome Trust Sanger Institute. SFF was also used as the native output format from early versions of Ion Torrent’s PGM platform as well. You are expected to use this module via the Bio.SeqIO functions under the format name “sff” (or “sff-trim” as described below).

For example, to iterate over the records in an SFF file,

>>> from Bio import SeqIO
>>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff"):
...     print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
...
E3MFGYR02JWQ7T 265 tcagGGTCTACATGTTGGTT...
E3MFGYR02JA6IL 271 tcagTTTTTTTTGGAAAGGA...
E3MFGYR02JHD4H 310 tcagAAAGACAAGTGGTATC...
E3MFGYR02GFKUC 299 tcagCGGCCGGGCCTCTCAT...
E3MFGYR02FTGED 281 tcagTGGTAATGGGGGGAAA...
E3MFGYR02FR9G7 261 tcagCTCCGTAAGAAGGTGC...
E3MFGYR02GAZMS 278 tcagAAAGAAGTAAGGTAAA...
E3MFGYR02HHZ8O 221 tcagACTTTCTTCTTTACCG...
E3MFGYR02GPGB1 269 tcagAAGCAGTGGTATCAAC...
E3MFGYR02F7Z7G 219 tcagAATCATCCACTTTTTA...

Each SeqRecord object will contain all the annotation from the SFF file, including the PHRED quality scores.

>>> print("%s %i" % (record.id, len(record)))
E3MFGYR02F7Z7G 219
>>> print("%s..." % record.seq[:10])
tcagAATCAT...
>>> print("%r..." % (record.letter_annotations["phred_quality"][:10]))
[22, 21, 23, 28, 26, 15, 12, 21, 28, 21]...

Notice that the sequence is given in mixed case, the central upper case region corresponds to the trimmed sequence. This matches the output of the Roche tools (and the 3rd party tool sff_extract) for SFF to FASTA.

>>> print(record.annotations["clip_qual_left"])
4
>>> print(record.annotations["clip_qual_right"])
134
>>> print(record.seq[:4])
tcag
>>> print("%s...%s" % (record.seq[4:20], record.seq[120:134]))
AATCATCCACTTTTTA...CAAAACACAAACAG
>>> print(record.seq[134:])
atcttatcaacaaaactcaaagttcctaactgagacacgcaacaggggataagacaaggcacacaggggataggnnnnnnnnnnn

The annotations dictionary also contains any adapter clip positions (usually zero), and information about the flows. e.g.

>>> len(record.annotations)
12
>>> print(record.annotations["flow_key"])
TCAG
>>> print(record.annotations["flow_values"][:10])
(83, 1, 128, 7, 4, 84, 6, 106, 3, 172)
>>> print(len(record.annotations["flow_values"]))
400
>>> print(record.annotations["flow_index"][:10])
(1, 2, 3, 2, 2, 0, 3, 2, 3, 3)
>>> print(len(record.annotations["flow_index"]))
219

Note that to convert from a raw reading in flow_values to the corresponding homopolymer stretch estimate, the value should be rounded to the nearest 100:

>>> print("%r..." % [int(round(value, -2)) // 100
...                  for value in record.annotations["flow_values"][:10]])
...
[1, 0, 1, 0, 0, 1, 0, 1, 0, 2]...

If a read name is exactly 14 alphanumeric characters, the annotations dictionary will also contain meta-data about the read extracted by interpreting the name as a 454 Sequencing System “Universal” Accession Number. Note that if a read name happens to be exactly 14 alphanumeric characters but was not generated automatically, these annotation records will contain nonsense information.

>>> print(record.annotations["region"])
2
>>> print(record.annotations["time"])
[2008, 1, 9, 16, 16, 0]
>>> print(record.annotations["coords"])
(2434, 1658)

As a convenience method, you can read the file with SeqIO format name “sff-trim” instead of “sff” to get just the trimmed sequences (without any annotation except for the PHRED quality scores and anything encoded in the read names):

>>> from Bio import SeqIO
>>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff-trim"):
...     print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
...
E3MFGYR02JWQ7T 260 GGTCTACATGTTGGTTAACC...
E3MFGYR02JA6IL 265 TTTTTTTTGGAAAGGAAAAC...
E3MFGYR02JHD4H 292 AAAGACAAGTGGTATCAACG...
E3MFGYR02GFKUC 295 CGGCCGGGCCTCTCATCGGT...
E3MFGYR02FTGED 277 TGGTAATGGGGGGAAATTTA...
E3MFGYR02FR9G7 256 CTCCGTAAGAAGGTGCTGCC...
E3MFGYR02GAZMS 271 AAAGAAGTAAGGTAAATAAC...
E3MFGYR02HHZ8O 150 ACTTTCTTCTTTACCGTAAC...
E3MFGYR02GPGB1 221 AAGCAGTGGTATCAACGCAG...
E3MFGYR02F7Z7G 130 AATCATCCACTTTTTAACGT...

Looking at the final record in more detail, note how this differs to the example above:

>>> print("%s %i" % (record.id, len(record)))
E3MFGYR02F7Z7G 130
>>> print("%s..." % record.seq[:10])
AATCATCCAC...
>>> print("%r..." % record.letter_annotations["phred_quality"][:10])
[26, 15, 12, 21, 28, 21, 36, 28, 27, 27]...
>>> len(record.annotations)
4
>>> print(record.annotations["region"])
2
>>> print(record.annotations["coords"])
(2434, 1658)
>>> print(record.annotations["time"])
[2008, 1, 9, 16, 16, 0]
>>> print(record.annotations["molecule_type"])
DNA

You might use the Bio.SeqIO.convert() function to convert the (trimmed) SFF reads into a FASTQ file (or a FASTA file and a QUAL file), e.g.

>>> from Bio import SeqIO
>>> from io import StringIO
>>> out_handle = StringIO()
>>> count = SeqIO.convert("Roche/E3MFGYR02_random_10_reads.sff", "sff",
...                       out_handle, "fastq")
...
>>> print("Converted %i records" % count)
Converted 10 records

The output FASTQ file would start like this:

>>> print("%s..." % out_handle.getvalue()[:50])
@E3MFGYR02JWQ7T
tcagGGTCTACATGTTGGTTAACCCGTACTGATT...

Bio.SeqIO.index() provides memory efficient random access to the reads in an SFF file by name. SFF files can include an index within the file, which can be read in making this very fast. If the index is missing (or in a format not yet supported in Biopython) the file is indexed by scanning all the reads - which is a little slower. For example,

>>> from Bio import SeqIO
>>> reads = SeqIO.index("Roche/E3MFGYR02_random_10_reads.sff", "sff")
>>> record = reads["E3MFGYR02JHD4H"]
>>> print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
E3MFGYR02JHD4H 310 tcagAAAGACAAGTGGTATC...
>>> reads.close()

Or, using the trimmed reads:

>>> from Bio import SeqIO
>>> reads = SeqIO.index("Roche/E3MFGYR02_random_10_reads.sff", "sff-trim")
>>> record = reads["E3MFGYR02JHD4H"]
>>> print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
E3MFGYR02JHD4H 292 AAAGACAAGTGGTATCAACG...
>>> reads.close()

You can also use the Bio.SeqIO.write() function with the “sff” format. Note that this requires all the flow information etc, and thus is probably only useful for SeqRecord objects originally from reading another SFF file (and not the trimmed SeqRecord objects from parsing an SFF file as “sff-trim”).

As an example, let’s pretend this example SFF file represents some DNA which was pre-amplified with a PCR primers AAAGANNNNN. The following script would produce a sub-file containing all those reads whose post-quality clipping region (i.e. the sequence after trimming) starts with AAAGA exactly (the non- degenerate bit of this pretend primer):

>>> from Bio import SeqIO
>>> records = (record for record in
...            SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff")
...            if record.seq[record.annotations["clip_qual_left"]:].startswith("AAAGA"))
...
>>> count = SeqIO.write(records, "temp_filtered.sff", "sff")
>>> print("Selected %i records" % count)
Selected 2 records

Of course, for an assembly you would probably want to remove these primers. If you want FASTA or FASTQ output, you could just slice the SeqRecord. However, if you want SFF output we have to preserve all the flow information - the trick is just to adjust the left clip position!

>>> from Bio import SeqIO
>>> def filter_and_trim(records, primer):
...     for record in records:
...         if record.seq[record.annotations["clip_qual_left"]:].startswith(primer):
...             record.annotations["clip_qual_left"] += len(primer)
...             yield record
...
>>> records = SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff")
>>> count = SeqIO.write(filter_and_trim(records, "AAAGA"),
...                     "temp_filtered.sff", "sff")
...
>>> print("Selected %i records" % count)
Selected 2 records

We can check the results, note the lower case clipped region now includes the “AAAGA” sequence:

>>> for record in SeqIO.parse("temp_filtered.sff", "sff"):
...     print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
...
E3MFGYR02JHD4H 310 tcagaaagaCAAGTGGTATC...
E3MFGYR02GAZMS 278 tcagaaagaAGTAAGGTAAA...
>>> for record in SeqIO.parse("temp_filtered.sff", "sff-trim"):
...     print("%s %i %s..." % (record.id, len(record), record.seq[:20]))
...
E3MFGYR02JHD4H 287 CAAGTGGTATCAACGCAGAG...
E3MFGYR02GAZMS 266 AGTAAGGTAAATAACAAACG...
>>> import os
>>> os.remove("temp_filtered.sff")

For a description of the file format, please see the Roche manuals and: http://www.ncbi.nlm.nih.gov/Traces/trace.cgi?cmd=show&f=formats&m=doc&s=formats

Bio.SeqIO.SffIO.ReadRocheXmlManifest(handle)

Read any Roche style XML manifest data in the SFF “index”.

The SFF file format allows for multiple different index blocks, and Roche took advantage of this to define their own index block which also embeds an XML manifest string. This is not a publicly documented extension to the SFF file format, this was reverse engineered.

The handle should be to an SFF file opened in binary mode. This function will use the handle seek/tell functions and leave the handle in an arbitrary location.

Any XML manifest found is returned as a Python string, which you can then parse as appropriate, or reuse when writing out SFF files with the SffWriter class.

Returns a string, or raises a ValueError if an Roche manifest could not be found.

class Bio.SeqIO.SffIO.SffIterator(source, alphabet=None, trim=False)

Bases: Bio.SeqIO.Interfaces.SequenceIterator

Parser for Standard Flowgram Format (SFF) files.

__init__(source, alphabet=None, trim=False)

Iterate over Standard Flowgram Format (SFF) reads (as SeqRecord objects).

  • source - path to an SFF file, e.g. from Roche 454 sequencing, or a file-like object opened in binary mode.

  • alphabet - optional alphabet, unused. Leave as None.

  • trim - should the sequences be trimmed?

The resulting SeqRecord objects should match those from a paired FASTA and QUAL file converted from the SFF file using the Roche 454 tool ssfinfo. i.e. The sequence will be mixed case, with the trim regions shown in lower case.

This function is used internally via the Bio.SeqIO functions:

>>> from Bio import SeqIO
>>> for record in SeqIO.parse("Roche/E3MFGYR02_random_10_reads.sff", "sff"):
...     print("%s %i" % (record.id, len(record)))
...
E3MFGYR02JWQ7T 265
E3MFGYR02JA6IL 271
E3MFGYR02JHD4H 310
E3MFGYR02GFKUC 299
E3MFGYR02FTGED 281
E3MFGYR02FR9G7 261
E3MFGYR02GAZMS 278
E3MFGYR02HHZ8O 221
E3MFGYR02GPGB1 269
E3MFGYR02F7Z7G 219

You can also call it directly:

>>> with open("Roche/E3MFGYR02_random_10_reads.sff", "rb") as handle:
...     for record in SffIterator(handle):
...         print("%s %i" % (record.id, len(record)))
...
E3MFGYR02JWQ7T 265
E3MFGYR02JA6IL 271
E3MFGYR02JHD4H 310
E3MFGYR02GFKUC 299
E3MFGYR02FTGED 281
E3MFGYR02FR9G7 261
E3MFGYR02GAZMS 278
E3MFGYR02HHZ8O 221
E3MFGYR02GPGB1 269
E3MFGYR02F7Z7G 219

Or, with the trim option:

>>> with open("Roche/E3MFGYR02_random_10_reads.sff", "rb") as handle:
...     for record in SffIterator(handle, trim=True):
...         print("%s %i" % (record.id, len(record)))
...
E3MFGYR02JWQ7T 260
E3MFGYR02JA6IL 265
E3MFGYR02JHD4H 292
E3MFGYR02GFKUC 295
E3MFGYR02FTGED 277
E3MFGYR02FR9G7 256
E3MFGYR02GAZMS 271
E3MFGYR02HHZ8O 150
E3MFGYR02GPGB1 221
E3MFGYR02F7Z7G 130
parse(handle)

Start parsing the file, and return a SeqRecord generator.

iterate(handle)

Parse the file and generate SeqRecord objects.

__abstractmethods__ = frozenset({})
class Bio.SeqIO.SffIO.SffWriter(target, index=True, xml=None)

Bases: Bio.SeqIO.Interfaces.SequenceWriter

SFF file writer.

__init__(target, index=True, xml=None)

Initialize an SFF writer object.

Arguments:
  • target - Output stream opened in binary mode, or a path to a file.

  • index - Boolean argument, should we try and write an index?

  • xml - Optional string argument, xml manifest to be recorded in the index block (see function ReadRocheXmlManifest for reading this data).

write_file(records)

Use this to write an entire file containing the given records.

write_header()

Write the SFF file header.

write_record(record)

Write a single additional record to the output file.

This assumes the header has been done.