Bio.Align package

Subpackages

Submodules

Module contents

Code for dealing with sequence alignments.

One of the most important things in this module is the MultipleSeqAlignment class, used in the Bio.AlignIO module.

class Bio.Align.AlignmentCounts(gaps, identities, mismatches)

Bases: tuple

__getnewargs__()

Return self as a plain tuple. Used by copy and pickle.

__match_args__ = ('gaps', 'identities', 'mismatches')
static __new__(_cls, gaps, identities, mismatches)

Create new instance of AlignmentCounts(gaps, identities, mismatches)

__repr__()

Return a nicely formatted representation string

__slots__ = ()
gaps

Alias for field number 0

identities

Alias for field number 1

mismatches

Alias for field number 2

class Bio.Align.MultipleSeqAlignment(records, alphabet=None, annotations=None, column_annotations=None)

Bases: object

Represents a classical multiple sequence alignment (MSA).

By this we mean a collection of sequences (usually shown as rows) which are all the same length (usually with gap characters for insertions or padding). The data can then be regarded as a matrix of letters, with well defined columns.

You would typically create an MSA by loading an alignment file with the AlignIO module:

>>> from Bio import AlignIO
>>> align = AlignIO.read("Clustalw/opuntia.aln", "clustal")
>>> print(align)
Alignment with 7 rows and 156 columns
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191
TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191

In some respects you can treat these objects as lists of SeqRecord objects, each representing a row of the alignment. Iterating over an alignment gives the SeqRecord object for each row:

>>> len(align)
7
>>> for record in align:
...     print("%s %i" % (record.id, len(record)))
...
gi|6273285|gb|AF191659.1|AF191 156
gi|6273284|gb|AF191658.1|AF191 156
gi|6273287|gb|AF191661.1|AF191 156
gi|6273286|gb|AF191660.1|AF191 156
gi|6273290|gb|AF191664.1|AF191 156
gi|6273289|gb|AF191663.1|AF191 156
gi|6273291|gb|AF191665.1|AF191 156

You can also access individual rows as SeqRecord objects via their index:

>>> print(align[0].id)
gi|6273285|gb|AF191659.1|AF191
>>> print(align[-1].id)
gi|6273291|gb|AF191665.1|AF191

And extract columns as strings:

>>> print(align[:, 1])
AAAAAAA

Or, take just the first ten columns as a sub-alignment:

>>> print(align[:, :10])
Alignment with 7 rows and 10 columns
TATACATTAA gi|6273285|gb|AF191659.1|AF191
TATACATTAA gi|6273284|gb|AF191658.1|AF191
TATACATTAA gi|6273287|gb|AF191661.1|AF191
TATACATAAA gi|6273286|gb|AF191660.1|AF191
TATACATTAA gi|6273290|gb|AF191664.1|AF191
TATACATTAA gi|6273289|gb|AF191663.1|AF191
TATACATTAA gi|6273291|gb|AF191665.1|AF191

Combining this alignment slicing with alignment addition allows you to remove a section of the alignment. For example, taking just the first and last ten columns:

>>> print(align[:, :10] + align[:, -10:])
Alignment with 7 rows and 20 columns
TATACATTAAGTGTACCAGA gi|6273285|gb|AF191659.1|AF191
TATACATTAAGTGTACCAGA gi|6273284|gb|AF191658.1|AF191
TATACATTAAGTGTACCAGA gi|6273287|gb|AF191661.1|AF191
TATACATAAAGTGTACCAGA gi|6273286|gb|AF191660.1|AF191
TATACATTAAGTGTACCAGA gi|6273290|gb|AF191664.1|AF191
TATACATTAAGTATACCAGA gi|6273289|gb|AF191663.1|AF191
TATACATTAAGTGTACCAGA gi|6273291|gb|AF191665.1|AF191

Note - This object does NOT attempt to model the kind of alignments used in next generation sequencing with multiple sequencing reads which are much shorter than the alignment, and where there is usually a consensus or reference sequence with special status.

__init__(records, alphabet=None, annotations=None, column_annotations=None)

Initialize a new MultipleSeqAlignment object.

Arguments:
  • records - A list (or iterator) of SeqRecord objects, whose

    sequences are all the same length. This may be an empty list.

  • alphabet - For backward compatibility only; its value should always

    be None.

  • annotations - Information about the whole alignment (dictionary).

  • column_annotations - Per column annotation (restricted dictionary).

    This holds Python sequences (lists, strings, tuples) whose length matches the number of columns. A typical use would be a secondary structure consensus string.

You would normally load a MSA from a file using Bio.AlignIO, but you can do this from a list of SeqRecord objects too:

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("AAAACGT"), id="Alpha")
>>> b = SeqRecord(Seq("AAA-CGT"), id="Beta")
>>> c = SeqRecord(Seq("AAAAGGT"), id="Gamma")
>>> align = MultipleSeqAlignment([a, b, c],
...                              annotations={"tool": "demo"},
...                              column_annotations={"stats": "CCCXCCC"})
>>> print(align)
Alignment with 3 rows and 7 columns
AAAACGT Alpha
AAA-CGT Beta
AAAAGGT Gamma
>>> align.annotations
{'tool': 'demo'}
>>> align.column_annotations
{'stats': 'CCCXCCC'}
property column_annotations

Dictionary of per-letter-annotation for the sequence.

__str__()

Return a multi-line string summary of the alignment.

This output is intended to be readable, but large alignments are shown truncated. A maximum of 20 rows (sequences) and 50 columns are shown, with the record identifiers. This should fit nicely on a single screen. e.g.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("ACTGCTAGCTAG"), id="Alpha")
>>> b = SeqRecord(Seq("ACT-CTAGCTAG"), id="Beta")
>>> c = SeqRecord(Seq("ACTGCTAGATAG"), id="Gamma")
>>> align = MultipleSeqAlignment([a, b, c])
>>> print(align)
Alignment with 3 rows and 12 columns
ACTGCTAGCTAG Alpha
ACT-CTAGCTAG Beta
ACTGCTAGATAG Gamma

See also the alignment’s format method.

__repr__()

Return a representation of the object for debugging.

The representation cannot be used with eval() to recreate the object, which is usually possible with simple python objects. For example:

<Bio.Align.MultipleSeqAlignment instance (2 records of length 14) at a3c184c>

The hex string is the memory address of the object, see help(id). This provides a simple way to visually distinguish alignments of the same size.

__format__(format_spec)

Return the alignment as a string in the specified file format.

The format should be a lower case string supported as an output format by Bio.AlignIO (such as “fasta”, “clustal”, “phylip”, “stockholm”, etc), which is used to turn the alignment into a string.

e.g.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("ACTGCTAGCTAG"), id="Alpha", description="")
>>> b = SeqRecord(Seq("ACT-CTAGCTAG"), id="Beta", description="")
>>> c = SeqRecord(Seq("ACTGCTAGATAG"), id="Gamma", description="")
>>> align = MultipleSeqAlignment([a, b, c])
>>> print(format(align, "fasta"))
>Alpha
ACTGCTAGCTAG
>Beta
ACT-CTAGCTAG
>Gamma
ACTGCTAGATAG

>>> print(format(align, "phylip"))
 3 12
Alpha      ACTGCTAGCT AG
Beta       ACT-CTAGCT AG
Gamma      ACTGCTAGAT AG
__iter__()

Iterate over alignment rows as SeqRecord objects.

e.g.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("ACTGCTAGCTAG"), id="Alpha")
>>> b = SeqRecord(Seq("ACT-CTAGCTAG"), id="Beta")
>>> c = SeqRecord(Seq("ACTGCTAGATAG"), id="Gamma")
>>> align = MultipleSeqAlignment([a, b, c])
>>> for record in align:
...    print(record.id)
...    print(record.seq)
...
Alpha
ACTGCTAGCTAG
Beta
ACT-CTAGCTAG
Gamma
ACTGCTAGATAG
__len__()

Return the number of sequences in the alignment.

Use len(alignment) to get the number of sequences (i.e. the number of rows), and alignment.get_alignment_length() to get the length of the longest sequence (i.e. the number of columns).

This is easy to remember if you think of the alignment as being like a list of SeqRecord objects.

get_alignment_length()

Return the maximum length of the alignment.

All objects in the alignment should (hopefully) have the same length. This function will go through and find this length by finding the maximum length of sequences in the alignment.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("ACTGCTAGCTAG"), id="Alpha")
>>> b = SeqRecord(Seq("ACT-CTAGCTAG"), id="Beta")
>>> c = SeqRecord(Seq("ACTGCTAGATAG"), id="Gamma")
>>> align = MultipleSeqAlignment([a, b, c])
>>> align.get_alignment_length()
12

If you want to know the number of sequences in the alignment, use len(align) instead:

>>> len(align)
3
extend(records)

Add more SeqRecord objects to the alignment as rows.

They must all have the same length as the original alignment. For example,

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("AAAACGT"), id="Alpha")
>>> b = SeqRecord(Seq("AAA-CGT"), id="Beta")
>>> c = SeqRecord(Seq("AAAAGGT"), id="Gamma")
>>> d = SeqRecord(Seq("AAAACGT"), id="Delta")
>>> e = SeqRecord(Seq("AAA-GGT"), id="Epsilon")

First we create a small alignment (three rows):

>>> align = MultipleSeqAlignment([a, b, c])
>>> print(align)
Alignment with 3 rows and 7 columns
AAAACGT Alpha
AAA-CGT Beta
AAAAGGT Gamma

Now we can extend this alignment with another two rows:

>>> align.extend([d, e])
>>> print(align)
Alignment with 5 rows and 7 columns
AAAACGT Alpha
AAA-CGT Beta
AAAAGGT Gamma
AAAACGT Delta
AAA-GGT Epsilon

Because the alignment object allows iteration over the rows as SeqRecords, you can use the extend method with a second alignment (provided its sequences have the same length as the original alignment).

append(record)

Add one more SeqRecord object to the alignment as a new row.

This must have the same length as the original alignment (unless this is the first record).

>>> from Bio import AlignIO
>>> align = AlignIO.read("Clustalw/opuntia.aln", "clustal")
>>> print(align)
Alignment with 7 rows and 156 columns
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191
TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191
>>> len(align)
7

We’ll now construct a dummy record to append as an example:

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> dummy = SeqRecord(Seq("N"*156), id="dummy")

Now append this to the alignment,

>>> align.append(dummy)
>>> print(align)
Alignment with 8 rows and 156 columns
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191
TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191
TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191
TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191
NNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNN...NNN dummy
>>> len(align)
8
__add__(other)

Combine two alignments with the same number of rows by adding them.

If you have two multiple sequence alignments (MSAs), there are two ways to think about adding them - by row or by column. Using the extend method adds by row. Using the addition operator adds by column. For example,

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a1 = SeqRecord(Seq("AAAAC"), id="Alpha")
>>> b1 = SeqRecord(Seq("AAA-C"), id="Beta")
>>> c1 = SeqRecord(Seq("AAAAG"), id="Gamma")
>>> a2 = SeqRecord(Seq("GT"), id="Alpha")
>>> b2 = SeqRecord(Seq("GT"), id="Beta")
>>> c2 = SeqRecord(Seq("GT"), id="Gamma")
>>> left = MultipleSeqAlignment([a1, b1, c1],
...                             annotations={"tool": "demo", "name": "start"},
...                             column_annotations={"stats": "CCCXC"})
>>> right = MultipleSeqAlignment([a2, b2, c2],
...                             annotations={"tool": "demo", "name": "end"},
...                             column_annotations={"stats": "CC"})

Now, let’s look at these two alignments:

>>> print(left)
Alignment with 3 rows and 5 columns
AAAAC Alpha
AAA-C Beta
AAAAG Gamma
>>> print(right)
Alignment with 3 rows and 2 columns
GT Alpha
GT Beta
GT Gamma

And add them:

>>> combined = left + right
>>> print(combined)
Alignment with 3 rows and 7 columns
AAAACGT Alpha
AAA-CGT Beta
AAAAGGT Gamma

For this to work, both alignments must have the same number of records (here they both have 3 rows):

>>> len(left)
3
>>> len(right)
3
>>> len(combined)
3

The individual rows are SeqRecord objects, and these can be added together. Refer to the SeqRecord documentation for details of how the annotation is handled. This example is a special case in that both original alignments shared the same names, meaning when the rows are added they also get the same name.

Any common annotations are preserved, but differing annotation is lost. This is the same behaviour used in the SeqRecord annotations and is designed to prevent accidental propagation of inappropriate values:

>>> combined.annotations
{'tool': 'demo'}

Similarly any common per-column-annotations are combined:

>>> combined.column_annotations
{'stats': 'CCCXCCC'}
__getitem__(index)

Access part of the alignment.

Depending on the indices, you can get a SeqRecord object (representing a single row), a Seq object (for a single column), a string (for a single character) or another alignment (representing some part or all of the alignment).

align[r,c] gives a single character as a string align[r] gives a row as a SeqRecord align[r,:] gives a row as a SeqRecord align[:,c] gives a column as a Seq

align[:] and align[:,:] give a copy of the alignment

Anything else gives a sub alignment, e.g. align[0:2] or align[0:2,:] uses only row 0 and 1 align[:,1:3] uses only columns 1 and 2 align[0:2,1:3] uses only rows 0 & 1 and only cols 1 & 2

We’ll use the following example alignment here for illustration:

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> a = SeqRecord(Seq("AAAACGT"), id="Alpha")
>>> b = SeqRecord(Seq("AAA-CGT"), id="Beta")
>>> c = SeqRecord(Seq("AAAAGGT"), id="Gamma")
>>> d = SeqRecord(Seq("AAAACGT"), id="Delta")
>>> e = SeqRecord(Seq("AAA-GGT"), id="Epsilon")
>>> align = MultipleSeqAlignment([a, b, c, d, e])

You can access a row of the alignment as a SeqRecord using an integer index (think of the alignment as a list of SeqRecord objects here):

>>> first_record = align[0]
>>> print("%s %s" % (first_record.id, first_record.seq))
Alpha AAAACGT
>>> last_record = align[-1]
>>> print("%s %s" % (last_record.id, last_record.seq))
Epsilon AAA-GGT

You can also access use python’s slice notation to create a sub-alignment containing only some of the SeqRecord objects:

>>> sub_alignment = align[2:5]
>>> print(sub_alignment)
Alignment with 3 rows and 7 columns
AAAAGGT Gamma
AAAACGT Delta
AAA-GGT Epsilon

This includes support for a step, i.e. align[start:end:step], which can be used to select every second sequence:

>>> sub_alignment = align[::2]
>>> print(sub_alignment)
Alignment with 3 rows and 7 columns
AAAACGT Alpha
AAAAGGT Gamma
AAA-GGT Epsilon

Or to get a copy of the alignment with the rows in reverse order:

>>> rev_alignment = align[::-1]
>>> print(rev_alignment)
Alignment with 5 rows and 7 columns
AAA-GGT Epsilon
AAAACGT Delta
AAAAGGT Gamma
AAA-CGT Beta
AAAACGT Alpha

You can also use two indices to specify both rows and columns. Using simple integers gives you the entry as a single character string. e.g.

>>> align[3, 4]
'C'

This is equivalent to:

>>> align[3][4]
'C'

or:

>>> align[3].seq[4]
'C'

To get a single column (as a string) use this syntax:

>>> align[:, 4]
'CCGCG'

Or, to get part of a column,

>>> align[1:3, 4]
'CG'

However, in general you get a sub-alignment,

>>> print(align[1:5, 3:6])
Alignment with 4 rows and 3 columns
-CG Beta
AGG Gamma
ACG Delta
-GG Epsilon

This should all seem familiar to anyone who has used the NumPy array or matrix objects.

__delitem__(index)

Delete SeqRecord by index or multiple SeqRecords by slice.

sort(key=None, reverse=False)

Sort the rows (SeqRecord objects) of the alignment in place.

This sorts the rows alphabetically using the SeqRecord object id by default. The sorting can be controlled by supplying a key function which must map each SeqRecord to a sort value.

This is useful if you want to add two alignments which use the same record identifiers, but in a different order. For example,

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> align1 = MultipleSeqAlignment([
...              SeqRecord(Seq("ACGT"), id="Human"),
...              SeqRecord(Seq("ACGG"), id="Mouse"),
...              SeqRecord(Seq("ACGC"), id="Chicken"),
...          ])
>>> align2 = MultipleSeqAlignment([
...              SeqRecord(Seq("CGGT"), id="Mouse"),
...              SeqRecord(Seq("CGTT"), id="Human"),
...              SeqRecord(Seq("CGCT"), id="Chicken"),
...          ])

If you simple try and add these without sorting, you get this:

>>> print(align1 + align2)
Alignment with 3 rows and 8 columns
ACGTCGGT <unknown id>
ACGGCGTT <unknown id>
ACGCCGCT Chicken

Consult the SeqRecord documentation which explains why you get a default value when annotation like the identifier doesn’t match up. However, if we sort the alignments first, then add them we get the desired result:

>>> align1.sort()
>>> align2.sort()
>>> print(align1 + align2)
Alignment with 3 rows and 8 columns
ACGCCGCT Chicken
ACGTCGTT Human
ACGGCGGT Mouse

As an example using a different sort order, you could sort on the GC content of each sequence.

>>> from Bio.SeqUtils import gc_fraction
>>> print(align1)
Alignment with 3 rows and 4 columns
ACGC Chicken
ACGT Human
ACGG Mouse
>>> align1.sort(key = lambda record: gc_fraction(record.seq))
>>> print(align1)
Alignment with 3 rows and 4 columns
ACGT Human
ACGC Chicken
ACGG Mouse

There is also a reverse argument, so if you wanted to sort by ID but backwards:

>>> align1.sort(reverse=True)
>>> print(align1)
Alignment with 3 rows and 4 columns
ACGG Mouse
ACGT Human
ACGC Chicken
property substitutions

Return an Array with the number of substitutions of letters in the alignment.

As an example, consider a multiple sequence alignment of three DNA sequences:

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> seq1 = SeqRecord(Seq("ACGT"), id="seq1")
>>> seq2 = SeqRecord(Seq("A--A"), id="seq2")
>>> seq3 = SeqRecord(Seq("ACGT"), id="seq3")
>>> seq4 = SeqRecord(Seq("TTTC"), id="seq4")
>>> alignment = MultipleSeqAlignment([seq1, seq2, seq3, seq4])
>>> print(alignment)
Alignment with 4 rows and 4 columns
ACGT seq1
A--A seq2
ACGT seq3
TTTC seq4
>>> m = alignment.substitutions
>>> print(m)
    A   C   G   T
A 3.0 0.5 0.0 2.5
C 0.5 1.0 0.0 2.0
G 0.0 0.0 1.0 1.0
T 2.5 2.0 1.0 1.0

Note that the matrix is symmetric, with counts divided equally on both sides of the diagonal. For example, the total number of substitutions between A and T in the alignment is 3.5 + 3.5 = 7.

Any weights associated with the sequences are taken into account when calculating the substitution matrix. For example, given the following multiple sequence alignment:

GTATC  0.5
AT--C  0.8
CTGTC  1.0

For the first column we have:

('A', 'G') : 0.5 * 0.8 = 0.4
('C', 'G') : 0.5 * 1.0 = 0.5
('A', 'C') : 0.8 * 1.0 = 0.8
property alignment

Return an Alignment object based on the MultipleSeqAlignment object.

This makes a copy of each SeqRecord with a gap-less sequence. Any future changes to the original records in the MultipleSeqAlignment will not affect the new records in the Alignment.

class Bio.Align.Alignment(sequences, coordinates=None)

Bases: object

Represents a sequence alignment.

An Alignment object has a .sequences attribute storing the sequences (Seq, MutableSeq, SeqRecord, or string objects) that were aligned, as well as a .coordinates attribute storing the sequence coordinates defining the alignment as a NumPy array.

Other commonly used attributes (which may or may not be present) are:
  • annotations - A dictionary with annotations describing the

    alignment;

  • column_annotations - A dictionary with annotations describing each

    column in the alignment;

  • score - The alignment score.

classmethod infer_coordinates(lines, skipped_columns=None)

Infer the coordinates from a printed alignment.

This method is primarily employed in Biopython’s alignment parsers, though it may be useful for other purposes.

For an alignment consisting of N sequences, printed as N lines with the same number of columns, where gaps are represented by dashes, this method will calculate the sequence coordinates that define the alignment. The coordinates are returned as a NumPy array of integers, and can be used to create an Alignment object.

The argument skipped columns should be None (the default) or an empty list. If skipped_columns is a list, then the indices of any columns in the alignment with a gap in all lines are appended to skipped_columns.

This is an example for the alignment of three sequences TAGGCATACGTG, AACGTACGT, and ACGCATACTTG, with gaps in the second and third sequence:

>>> from Bio.Align import Alignment
>>> lines = ["TAGGCATACGTG",
...          "AACG--TACGT-",
...          "-ACGCATACTTG",
...         ]
>>> sequences = [line.replace("-", "") for line in lines]
>>> sequences
['TAGGCATACGTG', 'AACGTACGT', 'ACGCATACTTG']
>>> coordinates = Alignment.infer_coordinates(lines)
>>> coordinates
array([[ 0,  1,  4,  6, 11, 12],
       [ 0,  1,  4,  4,  9,  9],
       [ 0,  0,  3,  5, 10, 11]])
>>> alignment = Alignment(sequences, coordinates)
__init__(sequences, coordinates=None)

Initialize a new Alignment object.

Arguments:
  • sequences - A list of the sequences (Seq, MutableSeq, SeqRecord,

    or string objects) that were aligned.

  • coordinates - The sequence coordinates that define the alignment.

    If None (the default value), assume that the sequences align to each other without any gaps.

__array__(dtype=None)
__add__(other)

Combine two alignments by adding them row-wise.

For example,

>>> import numpy as np
>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import Alignment
>>> a1 = SeqRecord(Seq("AAAAC"), id="Alpha")
>>> b1 = SeqRecord(Seq("AAAC"), id="Beta")
>>> c1 = SeqRecord(Seq("AAAAG"), id="Gamma")
>>> a2 = SeqRecord(Seq("GTT"), id="Alpha")
>>> b2 = SeqRecord(Seq("TT"), id="Beta")
>>> c2 = SeqRecord(Seq("GT"), id="Gamma")
>>> left = Alignment([a1, b1, c1],
...                  coordinates=np.array([[0, 3, 4, 5],
...                                        [0, 3, 3, 4],
...                                        [0, 3, 4, 5]]))
>>> left.annotations = {"tool": "demo", "name": "start"}
>>> left.column_annotations = {"stats": "CCCXC"}
>>> right = Alignment([a2, b2, c2],
...                   coordinates=np.array([[0, 1, 2, 3],
...                                         [0, 0, 1, 2],
...                                         [0, 1, 1, 2]]))
>>> right.annotations = {"tool": "demo", "name": "end"}
>>> right.column_annotations = {"stats": "CXC"}

Now, let’s look at these two alignments:

>>> print(left)
Alpha             0 AAAAC 5
Beta              0 AAA-C 4
Gamma             0 AAAAG 5

>>> print(right)
Alpha             0 GTT 3
Beta              0 -TT 2
Gamma             0 G-T 2

And add them:

>>> combined = left + right
>>> print(combined)
Alpha             0 AAAACGTT 8
Beta              0 AAA-C-TT 6
Gamma             0 AAAAGG-T 7

For this to work, both alignments must have the same number of sequences (here they both have 3 rows):

>>> len(left)
3
>>> len(right)
3
>>> len(combined)
3

The sequences are SeqRecord objects, and these can be added together. Refer to the SeqRecord documentation for details of how the annotation is handled. This example is a special case in that both original alignments shared the same names, meaning when the rows are added they also get the same name.

Any common annotations are preserved, but differing annotation is lost. This is the same behavior used in the SeqRecord annotations and is designed to prevent accidental propagation of inappropriate values:

>>> combined.annotations
{'tool': 'demo'}

Similarly any common per-column-annotations are combined:

>>> combined.column_annotations
{'stats': 'CCCXCCXC'}
property frequencies

Return the frequency of each letter in each column of the alignment.

Gaps are represented by a dash (“-”) character. For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = "global"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 -GACCT-G 6
                  0 -||--|-| 8
query             0 CGA--TCG 6

>>> alignment.frequencies
{'-': array([1., 0., 0., 1., 1., 0., 1., 0.]), 'G': array([0., 2., 0., 0., 0., 0., 0., 2.]), 'A': array([0., 0., 2., 0., 0., 0., 0., 0.]), 'C': array([1., 0., 0., 1., 1., 0., 1., 0.]), 'T': array([0., 0., 0., 0., 0., 2., 0., 0.])}
>>> aligner.mode = "local"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 GACCT-G 6
                  0 ||--|-| 7
query             1 GA--TCG 6

>>> alignment.frequencies
{'G': array([2., 0., 0., 0., 0., 0., 2.]), 'A': array([0., 2., 0., 0., 0., 0., 0.]), 'C': array([0., 0., 1., 1., 0., 1., 0.]), 'T': array([0., 0., 0., 0., 2., 0., 0.]), '-': array([0., 0., 1., 1., 0., 1., 0.])}
property target

Return self.sequences[0] for a pairwise alignment.

property query

Return self.sequences[1] for a pairwise alignment.

__eq__(other)

Check if two Alignment objects specify the same alignment.

__ne__(other)

Check if two Alignment objects have different alignments.

__lt__(other)

Check if self should come before other.

__le__(other)

Check if self should come before or is equal to other.

__gt__(other)

Check if self should come after other.

__ge__(other)

Check if self should come after or is equal to other.

__getitem__(key)

Return self[key].

Indices of the form

self[:, :]

return a copy of the Alignment object;

self[:, i:] self[:, :j] self[:, i:j] self[:, iterable] (where iterable returns integers)

return a new Alignment object spanning the selected columns;

self[k, i] self[k, i:] self[k, :j] self[k, i:j] self[k, iterable] (where iterable returns integers) self[k] (equivalent to self[k, :])

return a string with the aligned sequence (including gaps) for the selected columns, where k = 0 represents the target and k = 1 represents the query sequence; and

self[:, i]

returns a string with the selected column in the alignment.

>>> from Bio.Align import PairwiseAligner
>>> aligner = PairwiseAligner()
>>> alignments = aligner.align("ACCGGTTT", "ACGGGTT")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 ACCGG-TTT 8
                  0 ||-||-||- 9
query             0 AC-GGGTT- 7

>>> alignment[0, :]
'ACCGG-TTT'
>>> alignment[1, :]
'AC-GGGTT-'
>>> alignment[0]
'ACCGG-TTT'
>>> alignment[1]
'AC-GGGTT-'
>>> alignment[0, 1:-2]
'CCGG-T'
>>> alignment[1, 1:-2]
'C-GGGT'
>>> alignment[0, (1, 5, 2)]
'C-C'
>>> alignment[1, ::2]
'A-GT-'
>>> alignment[1, range(0, 9, 2)]
'A-GT-'
>>> alignment[:, 0]
'AA'
>>> alignment[:, 5]
'-G'
>>> alignment[:, 1:]  
<Alignment object (2 rows x 8 columns) at 0x...>
>>> print(alignment[:, 1:])
target            1 CCGG-TTT 8
                  0 |-||-||- 8
query             1 C-GGGTT- 7

>>> print(alignment[:, 2:])
target            2 CGG-TTT 8
                  0 -||-||- 7
query             2 -GGGTT- 7

>>> print(alignment[:, 3:])
target            3 GG-TTT 8
                  0 ||-||- 6
query             2 GGGTT- 7

>>> print(alignment[:, 3:-1])
target            3 GG-TT 7
                  0 ||-|| 5
query             2 GGGTT 7

>>> print(alignment[:, ::2])
target            0 ACGTT 5
                  0 |-||- 5
query             0 A-GT- 3

>>> print(alignment[:, range(1, 9, 2)])
target            0 CG-T 3
                  0 ||-| 4
query             0 CGGT 4

>>> print(alignment[:, (2, 7, 3)])
target            0 CTG 3
                  0 -|| 3
query             0 -TG 2
__format__(format_spec)

Return the alignment as a string in the specified file format.

Wrapper for self.format().

format(fmt='', *args, **kwargs)

Return the alignment as a string in the specified file format.

Arguments:
  • fmt - File format. Acceptable values are an empty string to

    create a human-readable representation of the alignment, or any of the alignment file formats supported by Bio.Align (some have not yet been implemented).

All other arguments are passed to the format-specific writer functions:
  • mask - PSL format only. Specify if repeat regions in the target

    sequence are masked and should be reported in the repMatches field of the PSL file instead of in the matches field. Acceptable values are None : no masking (default); “lower”: masking by lower-case characters; “upper”: masking by upper-case characters.

  • wildcard - PSL format only. Report alignments to the wildcard

    character in the target or query sequence in the nCount field of the PSL file instead of in the matches, misMatches, or repMatches fields. Default value is ‘N’.

  • md - SAM format only. If True, calculate the MD tag from

    the alignment and include it in the output. If False (default), do not include the MD tag in the output.

__str__()

Return a human-readable string representation of the alignment.

For sequence alignments, each line has at most 80 columns. The first 10 columns show the (possibly truncated) sequence name, which may be the id attribute of a SeqRecord, or otherwise ‘target’ or ‘query’ for pairwise alignments. The next 10 columns show the sequence coordinate, using zero-based counting as usual in Python. The remaining 60 columns shown the sequence, using dashes to represent gaps. At the end of the alignment, the end coordinates are shown on the right of the sequence, again in zero-based coordinates.

Pairwise alignments have an additional line between the two sequences showing whether the sequences match (‘|’) or mismatch (‘.’), or if there is a gap (‘-‘). The coordinates shown for this line are the column indices, which can be useful when extracting a subalignment.

For example,

>>> from Bio.Align import PairwiseAligner
>>> aligner = PairwiseAligner()
>>> seqA = "TTAACCCCATTTG"
>>> seqB = "AAGCCCCTTT"
>>> seqC = "AAAGGGGCTT"
>>> alignments = aligner.align(seqA, seqB)
>>> len(alignments)
1
>>> alignment = alignments[0]
>>> print(alignment)
target            0 TTAA-CCCCATTTG 13
                  0 --||-||||-|||- 14
query             0 --AAGCCCC-TTT- 10

Note that seqC is the reverse complement of seqB. Aligning it to the reverse strand gives the same alignment, but the query coordinates are switched:

>>> alignments = aligner.align(seqA, seqC, strand="-")
>>> len(alignments)
1
>>> alignment = alignments[0]
>>> print(alignment)
target            0 TTAA-CCCCATTTG 13
                  0 --||-||||-|||- 14
query            10 --AAGCCCC-TTT-  0
__repr__()

Return a representation of the alignment, including its shape.

The representation cannot be used with eval() to recreate the object, which is usually possible with simple python objects. For example:

<Alignment object (2 rows x 14 columns) at 0x10403d850>

The hex string is the memory address of the object and can be used to distinguish different Alignment objects. See help(id) for more information.

>>> import numpy as np
>>> from Bio.Align import Alignment
>>> alignment = Alignment(("ACCGT", "ACGT"),
...                       coordinates = np.array([[0, 2, 3, 5],
...                                               [0, 2, 2, 4],
...                                              ]))
>>> print(alignment)
target            0 ACCGT 5
                  0 ||-|| 5
query             0 AC-GT 4

>>> alignment  
<Alignment object (2 rows x 5 columns) at 0x...>
__len__()

Return the number of sequences in the alignment.

property length

Return the alignment length, i.e. the number of columns when printed..

The alignment length is the number of columns in the alignment when it is printed, and is equal to the sum of the number of matches, number of mismatches, and the total length of gaps in the target and query. Sequence sections beyond the aligned segment are not included in the number of columns.

For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = "global"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 -GACCT-G 6
                  0 -||--|-| 8
query             0 CGA--TCG 6

>>> alignment.length
8
>>> aligner.mode = "local"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 GACCT-G 6
                  0 ||--|-| 7
query             1 GA--TCG 6

>>> len(alignment)
2
>>> alignment.length
7
property shape

Return the shape of the alignment as a tuple of two integer values.

The first integer value is the number of sequences in the alignment as returned by len(alignment), which is always 2 for pairwise alignments.

The second integer value is the number of columns in the alignment when it is printed, and is equal to the sum of the number of matches, number of mismatches, and the total length of gaps in the target and query. Sequence sections beyond the aligned segment are not included in the number of columns.

For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = "global"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 -GACCT-G 6
                  0 -||--|-| 8
query             0 CGA--TCG 6

>>> len(alignment)
2
>>> alignment.shape
(2, 8)
>>> aligner.mode = "local"
>>> alignments = aligner.align("GACCTG", "CGATCG")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 GACCT-G 6
                  0 ||--|-| 7
query             1 GA--TCG 6

>>> len(alignment)
2
>>> alignment.shape
(2, 7)
property aligned

Return the indices of subsequences aligned to each other.

This property returns the start and end indices of subsequences in the target and query sequence that were aligned to each other. If the alignment between target (t) and query (q) consists of N chunks, you get two tuples of length N:

(((t_start1, t_end1), (t_start2, t_end2), …, (t_startN, t_endN)),

((q_start1, q_end1), (q_start2, q_end2), …, (q_startN, q_endN)))

For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> alignments = aligner.align("GAACT", "GAT")
>>> alignment = alignments[0]
>>> print(alignment)
target            0 GAACT 5
                  0 ||--| 5
query             0 GA--T 3

>>> alignment.aligned
array([[[0, 2],
        [4, 5]],

       [[0, 2],
        [2, 3]]])
>>> alignment = alignments[1]
>>> print(alignment)
target            0 GAACT 5
                  0 |-|-| 5
query             0 G-A-T 3

>>> alignment.aligned
array([[[0, 1],
        [2, 3],
        [4, 5]],

       [[0, 1],
        [1, 2],
        [2, 3]]])

Note that different alignments may have the same subsequences aligned to each other. In particular, this may occur if alignments differ from each other in terms of their gap placement only:

>>> aligner.mismatch_score = -10
>>> alignments = aligner.align("AAACAAA", "AAAGAAA")
>>> len(alignments)
2
>>> print(alignments[0])
target            0 AAAC-AAA 7
                  0 |||--||| 8
query             0 AAA-GAAA 7

>>> alignments[0].aligned
array([[[0, 3],
        [4, 7]],

       [[0, 3],
        [4, 7]]])
>>> print(alignments[1])
target            0 AAA-CAAA 7
                  0 |||--||| 8
query             0 AAAG-AAA 7

>>> alignments[1].aligned
array([[[0, 3],
        [4, 7]],

       [[0, 3],
        [4, 7]]])

The property can be used to identify alignments that are identical to each other in terms of their aligned sequences.

property indices

Return the sequence index of each lettter in the alignment.

This property returns a 2D NumPy array with the sequence index of each letter in the alignment. Gaps are indicated by -1. The array has the same number of rows and columns as the alignment, as given by self.shape.

For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = "local"
>>> alignments = aligner.align("GAACTGG", "AATG")
>>> alignment = alignments[0]
>>> print(alignment)
target            1 AACTG 6
                  0 ||-|| 5
query             0 AA-TG 4

>>> alignment.indices
array([[ 1,  2,  3,  4,  5],
       [ 0,  1, -1,  2,  3]])
>>> alignment = alignments[1]
>>> print(alignment)
target            1 AACTGG 7
                  0 ||-|-| 6
query             0 AA-T-G 4

>>> alignment.indices
array([[ 1,  2,  3,  4,  5,  6],
       [ 0,  1, -1,  2, -1,  3]])
>>> alignments = aligner.align("GAACTGG", "CATT", strand="-")
>>> alignment = alignments[0]
>>> print(alignment)
target            1 AACTG 6
                  0 ||-|| 5
query             4 AA-TG 0

>>> alignment.indices
array([[ 1,  2,  3,  4,  5],
       [ 3,  2, -1,  1,  0]])
>>> alignment = alignments[1]
>>> print(alignment)
target            1 AACTGG 7
                  0 ||-|-| 6
query             4 AA-T-G 0

>>> alignment.indices
array([[ 1,  2,  3,  4,  5,  6],
       [ 3,  2, -1,  1, -1,  0]])
property inverse_indices

Return the alignment column index for each letter in each sequence.

This property returns a list of 1D NumPy arrays; the number of arrays is equal to the number of aligned sequences, and the length of each array is equal to the length of the corresponding sequence. For each letter in each sequence, the array contains the corresponding column index in the alignment. Letters not included in the alignment are indicated by -1.

For example,

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = "local"
>>> alignments = aligner.align("GAACTGG", "AATG")
>>> alignment = alignments[0]
>>> print(alignment)
target            1 AACTG 6
                  0 ||-|| 5
query             0 AA-TG 4

>>> alignment.inverse_indices
[array([-1,  0,  1,  2,  3,  4, -1]), array([0, 1, 3, 4])]
>>> alignment = alignments[1]
>>> print(alignment)
target            1 AACTGG 7
                  0 ||-|-| 6
query             0 AA-T-G 4

>>> alignment.inverse_indices
[array([-1,  0,  1,  2,  3,  4,  5]), array([0, 1, 3, 5])]
>>> alignments = aligner.align("GAACTGG", "CATT", strand="-")
>>> alignment = alignments[0]
>>> print(alignment)
target            1 AACTG 6
                  0 ||-|| 5
query             4 AA-TG 0

>>> alignment.inverse_indices
[array([-1,  0,  1,  2,  3,  4, -1]), array([4, 3, 1, 0])]
>>> alignment = alignments[1]
>>> print(alignment)
target            1 AACTGG 7
                  0 ||-|-| 6
query             4 AA-T-G 0

>>> alignment.inverse_indices
[array([-1,  0,  1,  2,  3,  4,  5]), array([5, 3, 1, 0])]
sort(key=None, reverse=False)

Sort the sequences of the alignment in place.

By default, this sorts the sequences alphabetically using their id attribute if available, or by their sequence contents otherwise. For example,

>>> from Bio.Align import PairwiseAligner
>>> aligner = PairwiseAligner()
>>> aligner.gap_score = -1
>>> alignments = aligner.align("AATAA", "AAGAA")
>>> len(alignments)
1
>>> alignment = alignments[0]
>>> print(alignment)
target            0 AATAA 5
                  0 ||.|| 5
query             0 AAGAA 5

>>> alignment.sort()
>>> print(alignment)
target            0 AAGAA 5
                  0 ||.|| 5
query             0 AATAA 5

Alternatively, a key function can be supplied that maps each sequence to a sort value. For example, you could sort on the GC content of each sequence.

>>> from Bio.SeqUtils import gc_fraction
>>> alignment.sort(key=gc_fraction)
>>> print(alignment)
target            0 AATAA 5
                  0 ||.|| 5
query             0 AAGAA 5

You can reverse the sort order by passing reverse=True:

>>> alignment.sort(key=gc_fraction, reverse=True)
>>> print(alignment)
target            0 AAGAA 5
                  0 ||.|| 5
query             0 AATAA 5

The sequences are now sorted by decreasing GC content value.

map(alignment)

Map the alignment to self.target and return the resulting alignment.

Here, self.query and alignment.target are the same sequence.

A typical example is where self is the pairwise alignment between a chromosome and a transcript, the argument is the pairwise alignment between the transcript and a sequence (e.g., as obtained by RNA-seq), and we want to find the alignment of the sequence to the chromosome:

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> aligner.mode = 'local'
>>> aligner.open_gap_score = -1
>>> aligner.extend_gap_score = 0
>>> chromosome = "AAAAAAAACCCCCCCAAAAAAAAAAAGGGGGGAAAAAAAA"
>>> transcript = "CCCCCCCGGGGGG"
>>> alignments1 = aligner.align(chromosome, transcript)
>>> len(alignments1)
1
>>> alignment1 = alignments1[0]
>>> print(alignment1)
target            8 CCCCCCCAAAAAAAAAAAGGGGGG 32
                  0 |||||||-----------|||||| 24
query             0 CCCCCCC-----------GGGGGG 13

>>> sequence = "CCCCGGGG"
>>> alignments2 = aligner.align(transcript, sequence)
>>> len(alignments2)
1
>>> alignment2 = alignments2[0]
>>> print(alignment2)
target            3 CCCCGGGG 11
                  0 ||||||||  8
query             0 CCCCGGGG  8

>>> alignment = alignment1.map(alignment2)
>>> print(alignment)
target           11 CCCCAAAAAAAAAAAGGGG 30
                  0 ||||-----------|||| 19
query             0 CCCC-----------GGGG  8

>>> format(alignment, "psl")
'8\t0\t0\t0\t0\t0\t1\t11\t+\tquery\t8\t0\t8\ttarget\t40\t11\t30\t2\t4,4,\t0,4,\t11,26,\n'

Mapping the alignment does not depend on the sequence contents. If we delete the sequence contents, the same alignment is found in PSL format (though we obviously lose the ability to print the sequence alignment):

>>> alignment1.target = Seq(None, len(alignment1.target))
>>> alignment1.query = Seq(None, len(alignment1.query))
>>> alignment2.target = Seq(None, len(alignment2.target))
>>> alignment2.query = Seq(None, len(alignment2.query))
>>> alignment = alignment1.map(alignment2)
>>> format(alignment, "psl")
'8\t0\t0\t0\t0\t0\t1\t11\t+\tquery\t8\t0\t8\ttarget\t40\t11\t30\t2\t4,4,\t0,4,\t11,26,\n'

The map method can also be used to lift over an alignment between different genome assemblies. In this case, self is a DNA alignment between two genome assemblies, and the argument is an alignment of a transcript against one of the genome assemblies:

>>> np.set_printoptions(threshold=5)  # print 5 array elements per row
>>> chain = Align.read("Blat/panTro5ToPanTro6.over.chain", "chain")
>>> chain.sequences[0].id
'chr1'
>>> len(chain.sequences[0].seq)
228573443
>>> chain.sequences[1].id
'chr1'
>>> len(chain.sequences[1].seq)
224244399
>>> print(chain.coordinates)
[[122250000 122250400 122250400 ... 122909818 122909819 122909835]
 [111776384 111776784 111776785 ... 112019962 112019962 112019978]]

showing that the range 122250000:122909835 of chr1 on chimpanzee genome assembly panTro5 aligns to range 111776384:112019978 of chr1 of chimpanzee genome assembly panTro6.

>>> alignment = Align.read("Blat/est.panTro5.psl", "psl")
>>> alignment.sequences[0].id
'chr1'
>>> len(alignment.sequences[0].seq)
228573443
>>> alignment.sequences[1].id
'DC525629'
>>> len(alignment.sequences[1].seq)
407
>>> print(alignment.coordinates)
[[122835789 122835847 122840993 122841145 122907212 122907314]
 [       32        90        90       242       242       344]]

This shows that nucleotide range 32:344 of expressed sequence tag DC525629 aligns to range 122835789:122907314 of chr1 of chimpanzee genome assembly panTro5.

Note that the target sequence chain.sequences[0].seq and the target sequence alignment.sequences[0] have the same length:

>>> len(chain.sequences[0].seq) == len(alignment.sequences[0].seq)
True

We swap the target and query of the chain such that the query of the chain corresponds to the target of alignment:

>>> chain = chain[::-1]
>>> chain.sequences[0].id
'chr1'
>>> len(chain.sequences[0].seq)
224244399
>>> chain.sequences[1].id
'chr1'
>>> len(chain.sequences[1].seq)
228573443
>>> print(chain.coordinates)
[[111776384 111776784 111776785 ... 112019962 112019962 112019978]
 [122250000 122250400 122250400 ... 122909818 122909819 122909835]]

Now we can get the coordinates of DC525629 against chimpanzee genome assembly panTro6 by calling map on the chain, with alignment as the argument:

>>> lifted_alignment = chain.map(alignment)
>>> lifted_alignment.sequences[0].id
'chr1'
>>> len(lifted_alignment.sequences[0].seq)
224244399
>>> lifted_alignment.sequences[1].id
'DC525629'
>>> len(lifted_alignment.sequences[1].seq)
407
>>> print(lifted_alignment.coordinates)
[[111982717 111982775 111987921 111988073 112009200 112009302]
 [       32        90        90       242       242       344]]

This shows that nucleotide range 32:344 of expressed sequence tag DC525629 aligns to range 111982717:112009302 of chr1 of chimpanzee genome assembly panTro6. Note that the genome span of DC525629 on chimpanzee genome assembly panTro5 is 122907314 - 122835789 = 71525 bp, while on panTro6 the genome span is 112009302 - 111982717 = 26585 bp.

mapall(alignments)

Map each of the alignments to self, and return the mapped alignment.

property substitutions

Return an Array with the number of substitutions of letters in the alignment.

As an example, consider a sequence alignment of two RNA sequences:

>>> from Bio.Align import PairwiseAligner
>>> target = "ATACTTACCTGGCAGGGGAGATACCATGATCACGAAGGTGGTTTTCCCAGGGCGAGGCTTATCCATTGCACTCCGGATGTGCTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACTGCATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTTTCCCCTG"  # human spliceosomal small nuclear RNA U1
>>> query = "ATACTTACCTGACAGGGGAGGCACCATGATCACACAGGTGGTCCTCCCAGGGCGAGGCTCTTCCATTGCACTGCGGGAGGGTTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACTGTATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTATCCCCCG"  # sea lamprey spliceosomal small RNA U1
>>> aligner = PairwiseAligner()
>>> aligner.gap_score = -10
>>> alignments = aligner.align(target, query)
>>> len(alignments)
1
>>> alignment = alignments[0]
>>> print(alignment)
target            0 ATACTTACCTGGCAGGGGAGATACCATGATCACGAAGGTGGTTTTCCCAGGGCGAGGCTT
                  0 |||||||||||.||||||||..|||||||||||..|||||||..|||||||||||||||.
query             0 ATACTTACCTGACAGGGGAGGCACCATGATCACACAGGTGGTCCTCCCAGGGCGAGGCTC

target           60 ATCCATTGCACTCCGGATGTGCTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACT
                 60 .|||||||||||.|||..|.|.||||||||||||||||||||||||||||||||||||||
query            60 TTCCATTGCACTGCGGGAGGGTTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACT

target          120 GCATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTTTCCCCTG 164
                120 |.||||||||||||||||||||||||||||||||||.|||||.| 164
query           120 GTATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTATCCCCCG 164

>>> m = alignment.substitutions
>>> print(m)
     A    C    G    T
A 28.0  1.0  2.0  1.0
C  0.0 39.0  1.0  2.0
G  2.0  0.0 45.0  0.0
T  2.0  5.0  1.0 35.0

Note that the matrix is not symmetric: rows correspond to the target sequence, and columns to the query sequence. For example, the number of T’s in the target sequence that are aligned to a C in the query sequence is

>>> m['T', 'C']
5.0

and the number of C’s in the query sequence tat are aligned to a T in the query sequence is

>>> m['C', 'T']
2.0

For some applications (for example, to define a scoring matrix from the substitution matrix), a symmetric matrix may be preferred, which can be calculated as follows:

>>> m += m.transpose()
>>> m /= 2.0
>>> print(m)
     A    C    G    T
A 28.0  0.5  2.0  1.5
C  0.5 39.0  0.5  3.5
G  2.0  0.5 45.0  0.5
T  1.5  3.5  0.5 35.0

The matrix is now symmetric, with counts divided equally on both sides of the diagonal:

>>> m['C', 'T']
3.5
>>> m['T', 'C']
3.5

The total number of substitutions between T’s and C’s in the alignment is 3.5 + 3.5 = 7.

counts()

Return number of identities, mismatches, and gaps of a pairwise alignment.

>>> aligner = PairwiseAligner(mode='global', match_score=2, mismatch_score=-1)
>>> for alignment in aligner.align("TACCG", "ACG"):
...     print("Score = %.1f:" % alignment.score)
...     c = alignment.counts()  # namedtuple
...     print(f"{c.gaps} gaps, {c.identities} identities, {c.mismatches} mismatches")
...     print(alignment)
...
Score = 6.0:
2 gaps, 3 identities, 0 mismatches
target            0 TACCG 5
                  0 -||-| 5
query             0 -AC-G 3

Score = 6.0:
2 gaps, 3 identities, 0 mismatches
target            0 TACCG 5
                  0 -|-|| 5
query             0 -A-CG 3

This classifies each pair of letters in a pairwise alignment into gaps, perfect matches, or mismatches. It has been defined as a method (not a property) so that it may in future take optional argument(s) allowing the behavior to be customized. These three values are returned as a namedtuple. This is calculated for all the pairs of sequences in the alignment.

reverse_complement()

Reverse-complement the alignment and return it.

>>> sequences = ["ATCG", "AAG", "ATC"]
>>> coordinates = np.array([[0, 2, 3, 4], [0, 2, 2, 3], [0, 2, 3, 3]])
>>> alignment = Alignment(sequences, coordinates)
>>> print(alignment)
                  0 ATCG 4
                  0 AA-G 3
                  0 ATC- 3

>>> rc_alignment = alignment.reverse_complement()
>>> print(rc_alignment)
                  0 CGAT 4
                  0 C-TT 3
                  0 -GAT 3

The attribute column_annotations, if present, is associated with the reverse-complemented alignment, with its values in reverse order.

>>> alignment.column_annotations = {"score": [3, 2, 2, 2]}
>>> rc_alignment = alignment.reverse_complement()
>>> print(rc_alignment.column_annotations)
{'score': [2, 2, 2, 3]}
__hash__ = None
class Bio.Align.AlignmentsAbstractBaseClass

Bases: abc.ABC

Abstract base class for sequence alignments.

Most users will not need to use this class. It is used internally as a base class for the list-like Alignments class, and for the AlignmentIterator class in Bio.Align.interfaces, which itself is the abstract base class for the alignment parsers in Bio/Align/.

__iter__()

Iterate over the alignments as Alignment objects.

This method SHOULD NOT be overridden by any subclass.

abstract __next__()

Return the next alignment.

abstract rewind()

Rewind the iterator to let it loop over the alignments from the beginning.

abstract __len__()

Return the number of alignments.

__abstractmethods__ = frozenset({'__len__', '__next__', 'rewind'})
class Bio.Align.Alignments(alignments=())

Bases: Bio.Align.AlignmentsAbstractBaseClass, list

__init__(alignments=())
__next__()

Return the next alignment.

rewind()

Rewind the iterator to let it loop over the alignments from the beginning.

__len__()

Return the number of alignments.

__abstractmethods__ = frozenset({})
__annotations__ = {}
class Bio.Align.PairwiseAlignments(seqA, seqB, score, paths)

Bases: Bio.Align.AlignmentsAbstractBaseClass

Implements an iterator over pairwise alignments returned by the aligner.

This class also supports indexing, which is fast for increasing indices, but may be slow for random access of a large number of alignments.

Note that pairwise aligners can return an astronomical number of alignments, even for relatively short sequences, if they align poorly to each other. We therefore recommend to first check the number of alignments, accessible as len(alignments), which can be calculated quickly even if the number of alignments is very large.

__init__(seqA, seqB, score, paths)

Initialize a new PairwiseAlignments object.

Arguments:
  • seqA - The first sequence, as a plain string, without gaps.

  • seqB - The second sequence, as a plain string, without gaps.

  • score - The alignment score.

  • paths - An iterator over the paths in the traceback matrix;

    each path defines one alignment.

You would normally obtain a PairwiseAlignments object by calling aligner.align(seqA, seqB), where aligner is a PairwiseAligner object or a CodonAligner object.

__len__()

Return the number of alignments.

__iter__()

Iterate over the alignments as Alignment objects.

This method SHOULD NOT be overridden by any subclass.

__getitem__(index)
__next__()

Return the next alignment.

rewind()

Rewind the iterator to let it loop over the alignments from the beginning.

__abstractmethods__ = frozenset({})
__annotations__ = {}
class Bio.Align.PairwiseAligner(scoring=None, **kwargs)

Bases: _pairwisealigner.PairwiseAligner

Performs pairwise sequence alignment using dynamic programming.

This provides functions to get global and local alignments between two sequences. A global alignment finds the best concordance between all characters in two sequences. A local alignment finds just the subsequences that align the best.

To perform a pairwise sequence alignment, first create a PairwiseAligner object. This object stores the match and mismatch scores, as well as the gap scores. Typically, match scores are positive, while mismatch scores and gap scores are negative or zero. By default, the match score is 1, and the mismatch and gap scores are zero. Based on the values of the gap scores, a PairwiseAligner object automatically chooses the appropriate alignment algorithm (the Needleman-Wunsch, Smith-Waterman, Gotoh, or Waterman-Smith-Beyer global or local alignment algorithm).

Calling the “score” method on the aligner with two sequences as arguments will calculate the alignment score between the two sequences. Calling the “align” method on the aligner with two sequences as arguments will return a generator yielding the alignments between the two sequences.

Some examples:

>>> from Bio import Align
>>> aligner = Align.PairwiseAligner()
>>> alignments = aligner.align("TACCG", "ACG")
>>> for alignment in sorted(alignments):
...     print("Score = %.1f:" % alignment.score)
...     print(alignment)
...
Score = 3.0:
target            0 TACCG 5
                  0 -|-|| 5
query             0 -A-CG 3

Score = 3.0:
target            0 TACCG 5
                  0 -||-| 5
query             0 -AC-G 3

Specify the aligner mode as local to generate local alignments:

>>> aligner.mode = 'local'
>>> alignments = aligner.align("TACCG", "ACG")
>>> for alignment in sorted(alignments):
...     print("Score = %.1f:" % alignment.score)
...     print(alignment)
...
Score = 3.0:
target            1 ACCG 5
                  0 |-|| 4
query             0 A-CG 3

Score = 3.0:
target            1 ACCG 5
                  0 ||-| 4
query             0 AC-G 3

Do a global alignment. Identical characters are given 2 points, 1 point is deducted for each non-identical character.

>>> aligner.mode = 'global'
>>> aligner.match_score = 2
>>> aligner.mismatch_score = -1
>>> for alignment in aligner.align("TACCG", "ACG"):
...     print("Score = %.1f:" % alignment.score)
...     print(alignment)
...
Score = 6.0:
target            0 TACCG 5
                  0 -||-| 5
query             0 -AC-G 3

Score = 6.0:
target            0 TACCG 5
                  0 -|-|| 5
query             0 -A-CG 3

Same as above, except now 0.5 points are deducted when opening a gap, and 0.1 points are deducted when extending it.

>>> aligner.open_gap_score = -0.5
>>> aligner.extend_gap_score = -0.1
>>> aligner.target_end_gap_score = 0.0
>>> aligner.query_end_gap_score = 0.0
>>> for alignment in aligner.align("TACCG", "ACG"):
...     print("Score = %.1f:" % alignment.score)
...     print(alignment)
...
Score = 5.5:
target            0 TACCG 5
                  0 -|-|| 5
query             0 -A-CG 3

Score = 5.5:
target            0 TACCG 5
                  0 -||-| 5
query             0 -AC-G 3

The alignment function can also use known matrices already included in Biopython:

>>> from Bio.Align import substitution_matrices
>>> aligner = Align.PairwiseAligner()
>>> aligner.substitution_matrix = substitution_matrices.load("BLOSUM62")
>>> alignments = aligner.align("KEVLA", "EVL")
>>> alignments = list(alignments)
>>> print("Number of alignments: %d" % len(alignments))
Number of alignments: 1
>>> alignment = alignments[0]
>>> print("Score = %.1f" % alignment.score)
Score = 13.0
>>> print(alignment)
target            0 KEVLA 5
                  0 -|||- 5
query             0 -EVL- 3

You can also set the value of attributes directly during construction of the PairwiseAligner object by providing them as keyword arguments:

>>> aligner = Align.PairwiseAligner(mode='global', match_score=2, mismatch_score=-1)
>>> for alignment in aligner.align("TACCG", "ACG"):
...     print("Score = %.1f:" % alignment.score)
...     print(alignment)
...
Score = 6.0:
target            0 TACCG 5
                  0 -||-| 5
query             0 -AC-G 3

Score = 6.0:
target            0 TACCG 5
                  0 -|-|| 5
query             0 -A-CG 3
__init__(scoring=None, **kwargs)

Initialize a PairwiseAligner as specified by the keyword arguments.

If scoring is None, use the default scoring scheme match = 1.0, mismatch = 0.0, gap score = 0.0 If scoring is “blastn”, “megablast”, or “blastp”, use the default substitution matrix and gap scores for BLASTN, MEGABLAST, or BLASTP, respectively.

Loops over the remaining keyword arguments and sets them as attributes on the object.

__setattr__(key, value)

Implement setattr(self, name, value).

align(seqA, seqB, strand='+')

Return the alignments of two sequences using PairwiseAligner.

score(seqA, seqB, strand='+')

Return the alignment score of two sequences using PairwiseAligner.

__getstate__()
__setstate__(state)
class Bio.Align.CodonAligner(codon_table=None, anchor_len=10)

Bases: _codonaligner.CodonAligner

Aligns a nucleotide sequence to an amino acid sequence.

This class implements a dynamic programming algorithm to align a nucleotide sequence to an amino acid sequence.

__init__(codon_table=None, anchor_len=10)

Initialize a CodonAligner for a specific genetic code.

Arguments:
  • codon_table - a CodonTable object representing the genetic code. If codon_table is None, the standard genetic code is used.

score(seqA, seqB)

Return the alignment score of a protein sequence and nucleotide sequence.

Arguments:
  • seqA - the protein sequence of amino acids (plain string, Seq, MutableSeq, or SeqRecord).

  • seqB - the nucleotide sequence (plain string, Seq, MutableSeq, or SeqRecord); both DNA and RNA sequences are accepted.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> aligner = CodonAligner()
>>> dna = SeqRecord(Seq('ATGTCTCGT'), id='dna')
>>> pro = SeqRecord(Seq('MSR'), id='pro')
>>> score = aligner.score(pro, dna)
>>> print(score)
3.0
>>> rna = SeqRecord(Seq('AUGUCUCGU'), id='rna')
>>> score = aligner.score(pro, rna)
>>> print(score)
3.0

This is an example with a frame shift in the DNA sequence:

>>> dna = "ATGCTGGGCTCGAACGAGTCCGTGTATGCCCTAAGCTGAGCCCGTCG"
>>> pro = "MLGSNESRVCPKLSPS"
>>> len(pro)
16
>>> aligner.frameshift_score = -3.0
>>> score = aligner.score(pro, dna)
>>> print(score)
13.0

In the following example, the position of the frame shift is ambiguous:

>>> dna = 'TTTAAAAAAAAAAATTT'
>>> pro = 'FKKKKF'
>>> len(pro)
6
>>> aligner.frameshift_score = -1.0
>>> alignments = aligner.align(pro, dna)
>>> print(alignments.score)
5.0
>>> len(alignments)
3
>>> print(next(alignments))
target            0 F  K  K  K   4
query             0 TTTAAAAAAAAA 12

target            4 K  F    6
query            11 AAATTT 17

>>> print(next(alignments))
target            0 F  K  K   3
query             0 TTTAAAAAA 9

target            3 K  K  F    6
query             8 AAAAAATTT 17

>>> print(next(alignments))
target            0 F  K   2
query             0 TTTAAA 6

target            2 K  K  K  F    6
query             5 AAAAAAAAATTT 17

>>> print(next(alignments))
Traceback (most recent call last):
...
StopIteration
align(seqA, seqB)

Align a nucleotide sequence to its corresponding protein sequence.

Arguments:
  • seqA - the protein sequence of amino acids (plain string, Seq, MutableSeq, or SeqRecord).

  • seqB - the nucleotide sequence (plain string, Seq, MutableSeq, or SeqRecord); both DNA and RNA sequences are accepted.

Returns an iterator of Alignment objects.

>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> aligner = CodonAligner()
>>> dna = SeqRecord(Seq('ATGTCTCGT'), id='dna')
>>> pro = SeqRecord(Seq('MSR'), id='pro')
>>> alignments = aligner.align(pro, dna)
>>> alignment = alignments[0]
>>> print(alignment)
pro               0 M  S  R   3
dna               0 ATGTCTCGT 9

>>> rna = SeqRecord(Seq('AUGUCUCGU'), id='rna')
>>> alignments = aligner.align(pro, rna)
>>> alignment = alignments[0]
>>> print(alignment)
pro               0 M  S  R   3
rna               0 AUGUCUCGU 9

This is an example with a frame shift in the DNA sequence:

>>> dna = "ATGCTGGGCTCGAACGAGTCCGTGTATGCCCTAAGCTGAGCCCGTCG"
>>> pro = "MLGSNESRVCPKLSPS"
>>> alignments = aligner.align(pro, dna)
>>> alignment = alignments[0]
>>> print(alignment)
target            0 M  L  G  S  N  E  S   7
query             0 ATGCTGGGCTCGAACGAGTCC 21

target            7 R  V  C  P  K  L  S  P  S   16
query            20 CGTGTATGCCCTAAGCTGAGCCCGTCG 47
Bio.Align.write(alignments, target, fmt, *args, **kwargs)

Write alignments to a file.

Arguments:
  • alignments - An Alignments object, an iterator of Alignment objects, or a single Alignment.

  • target - File or file-like object to write to, or filename as string.

  • fmt - String describing the file format (case-insensitive).

Note if providing a file or file-like object, your code should close the target after calling this function, or call .flush(), to ensure the data gets flushed to disk.

Returns the number of alignments written (as an integer).

Bio.Align.parse(source, fmt)

Parse an alignment file and return an iterator over alignments.

Arguments:
  • source - File or file-like object to read from, or filename as string.

  • fmt - String describing the file format (case-insensitive).

Typical usage, opening a file to read in, and looping over the aligments:

>>> from Bio import Align
>>> filename = "Exonerate/exn_22_m_ner_cigar.exn"
>>> for alignment in Align.parse(filename, "exonerate"):
...    print("Number of sequences in alignment", len(alignment))
...    print("Alignment score:", alignment.score)
Number of sequences in alignment 2
Alignment score: 6150.0
Number of sequences in alignment 2
Alignment score: 502.0
Number of sequences in alignment 2
Alignment score: 440.0

For lazy-loading file formats such as bigMaf, for which the file contents is read on demand only, ensure that the file remains open while extracting alignment data.

You can use the Bio.Align.read(…) function when the file contains only one alignment.

Bio.Align.read(handle, fmt)

Parse a file containing one alignment, and return it.

Arguments:
  • source - File or file-like object to read from, or filename as string.

  • fmt - String describing the file format (case-insensitive).

This function is for use parsing alignment files containing exactly one alignment. For example, reading a Clustal file:

>>> from Bio import Align
>>> alignment = Align.read("Clustalw/opuntia.aln", "clustal")
>>> print("Alignment shape:", alignment.shape)
Alignment shape: (7, 156)
>>> for sequence in alignment.sequences:
...     print(sequence.id, len(sequence))
gi|6273285|gb|AF191659.1|AF191 146
gi|6273284|gb|AF191658.1|AF191 148
gi|6273287|gb|AF191661.1|AF191 146
gi|6273286|gb|AF191660.1|AF191 146
gi|6273290|gb|AF191664.1|AF191 150
gi|6273289|gb|AF191663.1|AF191 150
gi|6273291|gb|AF191665.1|AF191 156

If the file contains no records, or more than one record, an exception is raised. For example:

>>> from Bio import Align
>>> filename = "Exonerate/exn_22_m_ner_cigar.exn"
>>> alignment = Align.read(filename, "exonerate")
Traceback (most recent call last):
    ...
ValueError: More than one alignment found in file

Use the Bio.Align.parse function if you want to read a file containing more than one alignment.