Bio.AlignIO package¶
Submodules¶
Module contents¶
Multiple sequence alignment input/output as alignment objects.
The Bio.AlignIO interface is deliberately very similar to Bio.SeqIO, and in fact the two are connected internally. Both modules use the same set of file format names (lower case strings). From the user’s perspective, you can read in a PHYLIP file containing one or more alignments using Bio.AlignIO, or you can read in the sequences within these alignmenta using Bio.SeqIO.
Bio.AlignIO is also documented at http://biopython.org/wiki/AlignIO and by a whole chapter in our tutorial:
Input¶
For the typical special case when your file or handle contains one and only one alignment, use the function Bio.AlignIO.read(). This takes an input file handle (or in recent versions of Biopython a filename as a string), format string and optional number of sequences per alignment. It will return a single MultipleSeqAlignment object (or raise an exception if there isn’t just one alignment):
>>> from Bio import AlignIO
>>> align = AlignIO.read("Phylip/interlaced.phy", "phylip")
>>> print(align)
SingleLetterAlphabet() alignment with 3 rows and 384 columns
-----MKVILLFVLAVFTVFVSS---------------RGIPPE...I-- CYS1_DICDI
MAHARVLLLALAVLATAAVAVASSSSFADSNPIRPVTDRAASTL...VAA ALEU_HORVU
------MWATLPLLCAGAWLLGV--------PVCGAAELSVNSL...PLV CATH_HUMAN
For the general case, when the handle could contain any number of alignments, use the function Bio.AlignIO.parse(…) which takes the same arguments, but returns an iterator giving MultipleSeqAlignment objects (typically used in a for loop). If you want random access to the alignments by number, turn this into a list:
>>> from Bio import AlignIO
>>> alignments = list(AlignIO.parse("Emboss/needle.txt", "emboss"))
>>> print(alignments[2])
SingleLetterAlphabet() alignment with 2 rows and 120 columns
-KILIVDDQYGIRILLNEVFNKEGYQTFQAANGLQALDIVTKER...--- ref_rec
LHIVVVDDDPGTCVYIESVFAELGHTCKSFVRPEAAEEYILTHP...HKE gi|94967506|receiver
Most alignment file formats can be concatenated so as to hold as many different multiple sequence alignments as possible. One common example is the output of the tool seqboot in the PHLYIP suite. Sometimes there can be a file header and footer, as seen in the EMBOSS alignment output.
Output¶
Use the function Bio.AlignIO.write(…), which takes a complete set of Alignment objects (either as a list, or an iterator), an output file handle (or filename in recent versions of Biopython) and of course the file format:
from Bio import AlignIO
alignments = ...
count = SeqIO.write(alignments, "example.faa", "fasta")
If using a handle make sure to close it to flush the data to the disk:
from Bio import AlignIO
alignments = ...
with open("example.faa", "w") as handle:
count = SeqIO.write(alignments, handle, "fasta")
In general, you are expected to call this function once (with all your alignments) and then close the file handle. However, for file formats like PHYLIP where multiple alignments are stored sequentially (with no file header and footer), then multiple calls to the write function should work as expected when using handles.
If you are using a filename, the repeated calls to the write functions will overwrite the existing file each time.
Conversion¶
The Bio.AlignIO.convert(…) function allows an easy interface for simple alignment file format conversions. Additionally, it may use file format specific optimisations so this should be the fastest way too.
In general however, you can combine the Bio.AlignIO.parse(…) function with the Bio.AlignIO.write(…) function for sequence file conversion. Using generator expressions provides a memory efficient way to perform filtering or other extra operations as part of the process.
File Formats¶
When specifying the file format, use lowercase strings. The same format names are also used in Bio.SeqIO and include the following:
clustal - Output from Clustal W or X, see also the module Bio.Clustalw which can be used to run the command line tool from Biopython.
emboss - EMBOSS tools’ “pairs” and “simple” alignment formats.
fasta - The generic sequence file format where each record starts with an identifer line starting with a “>” character, followed by lines of sequence.
fasta-m10 - For the pairswise alignments output by Bill Pearson’s FASTA tools when used with the -m 10 command line option for machine readable output.
ig - The IntelliGenetics file format, apparently the same as the MASE alignment format.
msf - The GCG MSF alignment format, originally from PileUp tool.
nexus - Output from NEXUS, see also the module Bio.Nexus which can also read any phylogenetic trees in these files.
phylip - Interlaced PHYLIP, as used by the PHLIP tools.
phylip-sequential - Sequential PHYLIP.
phylip-relaxed - PHYLIP like format allowing longer names.
stockholm - A richly annotated alignment file format used by PFAM.
mauve - Output from progressiveMauve/Mauve
Note that while Bio.AlignIO can read all the above file formats, it cannot write to all of them.
You can also use any file format supported by Bio.SeqIO, such as “fasta” or “ig” (which are listed above), PROVIDED the sequences in your file are all the same length.
-
Bio.AlignIO.
write
(alignments, handle, format)¶ Write complete set of alignments to a file.
- Arguments:
alignments - A list (or iterator) of MultipleSeqAlignment objects, or a single alignment object.
handle - File handle object to write to, or filename as string (note older versions of Biopython only took a handle).
format - lower case string describing the file format to write.
You should close the handle after calling this function.
Returns the number of alignments written (as an integer).
-
Bio.AlignIO.
parse
(handle, format, seq_count=None, alphabet=None)¶ Iterate over an alignment file as MultipleSeqAlignment objects.
- Arguments:
handle - handle to the file, or the filename as a string (note older versions of Biopython only took a handle).
format - string describing the file format.
alphabet - optional Alphabet object, useful when the sequence type cannot be automatically inferred from the file itself (e.g. fasta, phylip, clustal)
seq_count - Optional integer, number of sequences expected in each alignment. Recommended for fasta format files.
If you have the file name in a string ‘filename’, use:
>>> from Bio import AlignIO >>> filename = "Emboss/needle.txt" >>> format = "emboss" >>> for alignment in AlignIO.parse(filename, format): ... print("Alignment of length %i" % alignment.get_alignment_length()) Alignment of length 124 Alignment of length 119 Alignment of length 120 Alignment of length 118 Alignment of length 125
If you have a string ‘data’ containing the file contents, use:
from Bio import AlignIO from StringIO import StringIO my_iterator = AlignIO.parse(StringIO(data), format)
Use the Bio.AlignIO.read() function when you expect a single record only.
-
Bio.AlignIO.
read
(handle, format, seq_count=None, alphabet=None)¶ Turn an alignment file into a single MultipleSeqAlignment object.
- Arguments:
handle - handle to the file, or the filename as a string (note older versions of Biopython only took a handle).
format - string describing the file format.
alphabet - optional Alphabet object, useful when the sequence type cannot be automatically inferred from the file itself (e.g. fasta, phylip, clustal)
seq_count - Optional integer, number of sequences expected in each alignment. Recommended for fasta format files.
If the handle contains no alignments, or more than one alignment, an exception is raised. For example, using a PFAM/Stockholm file containing one alignment:
>>> from Bio import AlignIO >>> filename = "Clustalw/protein.aln" >>> format = "clustal" >>> alignment = AlignIO.read(filename, format) >>> print("Alignment of length %i" % alignment.get_alignment_length()) Alignment of length 411
If however you want the first alignment from a file containing multiple alignments this function would raise an exception.
>>> from Bio import AlignIO >>> filename = "Emboss/needle.txt" >>> format = "emboss" >>> alignment = AlignIO.read(filename, format) Traceback (most recent call last): ... ValueError: More than one record found in handle
Instead use:
>>> from Bio import AlignIO >>> filename = "Emboss/needle.txt" >>> format = "emboss" >>> alignment = next(AlignIO.parse(filename, format)) >>> print("First alignment has length %i" % alignment.get_alignment_length()) First alignment has length 124
You must use the Bio.AlignIO.parse() function if you want to read multiple records from the handle.
-
Bio.AlignIO.
convert
(in_file, in_format, out_file, out_format, alphabet=None)¶ Convert between two alignment files, returns number of alignments.
- Arguments:
in_file - an input handle or filename
in_format - input file format, lower case string
output - an output handle or filename
out_file - output file format, lower case string
alphabet - optional alphabet to assume
NOTE - If you provide an output filename, it will be opened which will overwrite any existing file without warning. This may happen if even the conversion is aborted (e.g. an invalid out_format name is given).