Trees | Indices | Help |
---|
|
Sequence input/output as SeqRecord objects.
Bio.SeqIO is also documented at http://biopython.org/wiki/SeqIO and by a whole chapter in our tutorial:
The main function is Bio.SeqIO.parse(...) which takes an input file handle, and format string. This returns an iterator giving SeqRecord objects:
>>> from Bio import SeqIO >>> handle = open("Fasta/f002", "rU") >>> for record in SeqIO.parse(handle, "fasta") : ... print record.id, len(record) gi|1348912|gb|G26680|G26680 633 gi|1348917|gb|G26685|G26685 413 gi|1592936|gb|G29385|G29385 471 >>> handle.close()
Note that the parse() function will all invoke the relevant parser for the format with its default settings. You may want more control, in which case you need to create a format specific sequence iterator directly.
For non-interlaced files (e.g. Fasta, GenBank, EMBL) with multiple records using a sequence iterator can save you a lot of memory (RAM). There is less benefit for interlaced file formats (e.g. most multiple alignment file formats). However, an iterator only lets you access the records one by one.
If you want random access to the records by number, turn this into a list:
>>> from Bio import SeqIO >>> handle = open("Fasta/f002", "rU") >>> records = list(SeqIO.parse(handle, "fasta")) >>> handle.close() >>> print records[1].id gi|1348917|gb|G26685|G26685
If you want random access to the records by a key such as the record id, turn the iterator into a dictionary:
>>> from Bio import SeqIO >>> handle = open("Fasta/f002", "rU") >>> record_dict = SeqIO.to_dict(SeqIO.parse(handle, "fasta")) >>> handle.close() >>> print len(record_dict["gi|1348917|gb|G26685|G26685"]) 413
If you expect your file to contain one-and-only-one record, then we provide the following 'helper' function which will return a single SeqRecord, or raise an exception if there are no records or more than one record:
>>> from Bio import SeqIO >>> handle = open("Fasta/f001", "rU") >>> record = SeqIO.read(handle, "fasta") >>> handle.close() >>> print record.id, len(record) gi|3318709|pdb|1A91| 79
This style is useful when you expect a single record only (and would consider multiple records an error). For example, when dealing with GenBank files for bacterial genomes or chromosomes, there is normally only a single record. Alternatively, use this with a handle when download a single record from the internet.
However, if you just want the first record from a file containing multiple record, use the iterator's next() method:
>>> from Bio import SeqIO >>> handle = open("Fasta/f002", "rU") >>> record = SeqIO.parse(handle, "fasta").next() >>> handle.close() >>> print record.id, len(record) gi|1348912|gb|G26680|G26680 633
The above code will work as long as the file contains at least one record. Note that if there is more than one record, the remaining records will be silently ignored.
You can read in alignment files as Alignment objects using Bio.AlignIO. Alternatively, reading in an alignment file format via Bio.SeqIO will give you a SeqRecord for each row of each alignment:
>>> from Bio import SeqIO >>> handle = open("Clustalw/hedgehog.aln", "rU") >>> for record in SeqIO.parse(handle, "clustal") : ... print record.id, len(record) gi|167877390|gb|EDS40773.1| 447 gi|167234445|ref|NP_001107837. 447 gi|74100009|gb|AAZ99217.1| 447 gi|13990994|dbj|BAA33523.2| 447 gi|56122354|gb|AAV74328.1| 447 >>> handle.close()
Use the function Bio.SeqIO.write(...), which takes a complete set of SeqRecord objects (either as a list, or an iterator), an output file handle and of course the file format:
from Bio import SeqIO records = ... handle = open("example.faa", "w") SeqIO.write(records, handle, "fasta") handle.close()
In general, you are expected to call this function once (with all your records) and then close the file handle.
The effect of calling write() multiple times on a single file will vary depending on the file format, and is best avoided unless you have a strong reason to do so.
Trying this for certain alignment formats (e.g. phylip, clustal, stockholm) would have the effect of concatenating several multiple sequence alignments together. Such files are created by the PHYLIP suite of programs for bootstrap analysis.
For sequential files formats (e.g. fasta, genbank) each "record block" holds a single sequence. For these files it would probably be safe to call write() multiple times.
When specifying the file format, use lowercase strings. The same format names are also used in Bio.AlignIO and include the following:
Note that while Bio.SeqIO can read all the above file formats, it cannot write to all of them.
You can also use any file format supported by Bio.AlignIO, such as "nexus", "phlip" and "stockholm", which gives you access to the individual sequences making up each alignment as SeqRecords.
|
|||
|
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|
|||
_FormatToIterator = {"fasta": FastaIO.FastaIterator, "gb": Ins
|
|||
_FormatToWriter = {"fasta": FastaIO.FastaWriter, "gb": InsdcIO
|
|
Write complete set of sequences to a file.
You should close the handle after calling this function. Returns the number of records written (as an integer). |
Turns a sequence file into an iterator returning SeqRecords.
Typical usage, opening a file to read in, and looping over the record(s): >>> from Bio import SeqIO >>> filename = "Nucleic/sweetpea.nu" >>> for record in SeqIO.parse(open(filename,"rU"), "fasta") : ... print "ID", record.id ... print "Sequence length", len(record) ... print "Sequence alphabet", record.seq.alphabet ID gi|3176602|gb|U78617.1|LOU78617 Sequence length 309 Sequence alphabet SingleLetterAlphabet() For file formats like FASTA where the alphabet cannot be determined, it may be useful to specify the alphabet explicitly: >>> from Bio import SeqIO >>> from Bio.Alphabet import generic_dna >>> filename = "Nucleic/sweetpea.nu" >>> for record in SeqIO.parse(open(filename,"rU"), "fasta", generic_dna) : ... print "ID", record.id ... print "Sequence length", len(record) ... print "Sequence alphabet", record.seq.alphabet ID gi|3176602|gb|U78617.1|LOU78617 Sequence length 309 Sequence alphabet DNAAlphabet() If you have a string 'data' containing the file contents, you must first turn this into a handle in order to parse it: >>> data = ">Alpha\nACCGGATGTA\n>Beta\nAGGCTCGGTTA\n" >>> from Bio import SeqIO >>> from StringIO import StringIO >>> for record in SeqIO.parse(StringIO(data), "fasta") : ... print record.id, record.seq Alpha ACCGGATGTA Beta AGGCTCGGTTA Use the Bio.SeqIO.read(handle, format) function when you expect a single record only. |
Turns a sequence file into a single SeqRecord.
This function is for use parsing sequence files containing exactly one record. For example, reading a GenBank file: >>> from Bio import SeqIO >>> record = SeqIO.read(open("GenBank/arab1.gb", "rU"), "genbank") >>> print "ID", record.id ID AC007323.5 >>> print "Sequence length", len(record) Sequence length 86436 >>> print "Sequence alphabet", record.seq.alphabet Sequence alphabet IUPACAmbiguousDNA() If the handle contains no records, or more than one record, an exception is raised. For example: >>> from Bio import SeqIO >>> record = SeqIO.read(open("GenBank/cor6_6.gb", "rU"), "genbank") Traceback (most recent call last): ... ValueError: More than one record found in handle If however you want the first record from a file containing multiple records this function would raise an exception (as shown in the example above). Instead use: >>> from Bio import SeqIO >>> record = SeqIO.parse(open("GenBank/cor6_6.gb", "rU"), "genbank").next() >>> print "First record's ID", record.id First record's ID X55053.1 Use the Bio.SeqIO.parse(handle, format) function if you want to read multiple records from the handle. |
Turns a sequence iterator or list into a dictionary.
e.g. key_function = lambda rec : rec.name or, key_function = lambda rec : rec.description.split()[0] If key_function is ommitted then record.id is used, on the assumption that the records objects returned are SeqRecords with a unique id field. If there are duplicate keys, an error is raised. Example usage, defaulting to using the record.id as key: >>> from Bio import SeqIO >>> handle = open("GenBank/cor6_6.gb", "rU") >>> format = "genbank" >>> id_dict = SeqIO.to_dict(SeqIO.parse(handle, format)) >>> print id_dict.keys() ['L31939.1', 'AJ237582.1', 'X62281.1', 'AF297471.1', 'X55053.1', 'M81224.1'] >>> print id_dict["L31939.1"].description Brassica rapa (clone bif72) kin mRNA, complete cds. A more complex example, using the key_function argument in order to use a sequence checksum as the dictionary key: >>> from Bio import SeqIO >>> from Bio.SeqUtils.CheckSum import seguid >>> handle = open("GenBank/cor6_6.gb", "rU") >>> format = "genbank" >>> seguid_dict = SeqIO.to_dict(SeqIO.parse(handle, format), ... key_function = lambda rec : seguid(rec.seq)) >>> for key, record in seguid_dict.iteritems() : ... print key, record.id SabZaA4V2eLE9/2Fm5FnyYy07J4 X55053.1 l7gjJFE6W/S1jJn5+1ASrUKW/FA X62281.1 /wQvmrl87QWcm9llO4/efg23Vgg AJ237582.1 TtWsXo45S3ZclIBy4X/WJc39+CY M81224.1 uVEYeAQSV5EDQOnFoeMmVea+Oow AF297471.1 BUg6YxXSKWEcFFH0L08JzaLGhQs L31939.1 |
Returns a multiple sequence alignment (OBSOLETE).
Using this function is now discouraged. Rather doing this: >>> from Bio import SeqIO >>> handle = open("Clustalw/protein.aln") >>> alignment = SeqIO.to_alignment(SeqIO.parse(handle, "clustal")) >>> handle.close() You are now encouraged to use Bio.AlignIO instead, e.g. >>> from Bio import AlignIO >>> handle = open("Clustalw/protein.aln") >>> alignment = AlignIO.read(handle, "clustal") >>> handle.close() |
Run the Bio.SeqIO module's doctests. This will try and locate the unit tests directory, and run the doctests from there in order that the relative paths used in the examples work. |
|
_FormatToIterator
|
_FormatToWriter
|
Trees | Indices | Help |
---|
Generated by Epydoc 3.0.1 on Sun May 3 15:51:23 2009 | http://epydoc.sourceforge.net |