Bio.AlignIO.FastaIO module

Bio.AlignIO support for “fasta-m10” output from Bill Pearson’s FASTA tools.

You are expected to use this module via the Bio.AlignIO functions (or the Bio.SeqIO functions if you want to work directly with the gapped sequences).

This module contains a parser for the pairwise alignments produced by Bill Pearson’s FASTA tools, for use from the Bio.AlignIO interface where it is referred to as the “fasta-m10” file format (as we only support the machine readable output format selected with the -m 10 command line option).

This module does NOT cover the generic “fasta” file format originally developed as an input format to the FASTA tools. The Bio.AlignIO and Bio.SeqIO both use the Bio.SeqIO.FastaIO module to deal with these files, which can also be used to store a multiple sequence alignments.

Bio.AlignIO.FastaIO.FastaM10Iterator(handle, alphabet=SingleLetterAlphabet())

Alignment iterator for the FASTA tool’s pairwise alignment output.

This is for reading the pairwise alignments output by Bill Pearson’s FASTA program when called with the -m 10 command line option for machine readable output. For more details about the FASTA tools, see the website http://fasta.bioch.virginia.edu/ and the paper:

W.R. Pearson & D.J. Lipman PNAS (1988) 85:2444-2448

This class is intended to be used via the Bio.AlignIO.parse() function by specifying the format as “fasta-m10” as shown in the following code:

from Bio import AlignIO
handle = ...
for a in AlignIO.parse(handle, "fasta-m10"):
    assert len(a) == 2, "Should be pairwise!"
    print("Alignment length %i" % a.get_alignment_length())
    for record in a:
        print("%s %s %s" % (record.seq, record.name, record.id))

Note that this is not a full blown parser for all the information in the FASTA output - for example, most of the header and all of the footer is ignored. Also, the alignments are not batched according to the input queries.

Also note that there can be up to about 30 letters of flanking region included in the raw FASTA output as contextual information. This is NOT part of the alignment itself, and is not included in the resulting MultipleSeqAlignment objects returned.