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.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 replaced the older Alignment object defined in module Bio.Align.Generic but is not fully backwards compatible with it.
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__(self, 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 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__(self)¶
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__(self)¶
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__(self, 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.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__(self)¶
Iterate over alignment rows as SeqRecord objects.
e.g.
>>> 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__(self)¶
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(self)¶
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.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(self, 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(self, 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__(self, 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__(self, 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 columns), a string (for a single characters) 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.
- sort(self, 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 >>> print(align1) Alignment with 3 rows and 4 columns ACGC Chicken ACGT Human ACGG Mouse >>> align1.sort(key = lambda record: GC(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
- class Bio.Align.PairwiseAlignment(target, query, path, score)¶
Bases:
object
Represents a pairwise sequence alignment.
Internally, the pairwise alignment is stored as the path through the traceback matrix, i.e. a tuple of pairs of indices corresponding to the vertices of the path in the traceback matrix.
- __init__(self, target, query, path, score)¶
Initialize a new PairwiseAlignment object.
- Arguments:
target - The first sequence, as a plain string, without gaps.
query - The second sequence, as a plain string, without gaps.
- path - The path through the traceback matrix, defining an
alignment.
score - The alignment score.
You would normally obtain a PairwiseAlignment object by iterating over a PairwiseAlignments object.
- __eq__(self, other)¶
Check if two PairwiseAlignment objects have the same path.
- __ne__(self, other)¶
Check if two PairwiseAlignment objects have different paths.
- __lt__(self, other)¶
Check if self should come before other.
- __le__(self, other)¶
Check if self should come before or is equal to other.
- __gt__(self, other)¶
Check if self should come after other.
- __ge__(self, other)¶
Check if self should come after or is equal to other.
- __getitem__(self, key)¶
Return self[key].
Currently, this is implemented only for indices of the form
self[:, :]
which returns a copy of the PairwiseAlignment object;
self[:, i:] self[:, :j] self[:, i:j]
which returns a new PairwiseAlignment object spanning the indicated columns; and
self[k, i] self[k, i:] self[k, :j] self[k, i:j]
which returns a string with the aligned sequence (including gaps) for the indicated columns, where k = 0 represents the target and k = 1 represents the query sequence.
>>> from Bio.Align import PairwiseAligner >>> aligner = PairwiseAligner() >>> alignments = aligner.align("ACCGGTTT", "ACGGGTT") >>> alignment = alignments[0] >>> print(alignment) ACCGG-TTT ||-||-||- AC-GGGTT- >>> alignment[0, :] 'ACCGG-TTT' >>> alignment[1, :] 'AC-GGGTT-' >>> alignment[0, 1:-2] 'CCGG-T' >>> alignment[1, 1:-2] 'C-GGGT' >>> alignment[:, 1:] <Bio.Align.PairwiseAlignment object at ...> >>> print(alignment[:, 1:]) ACCGG-TTT |-||-||- AC-GGGTT- >>> print(alignment[:, 2:]) ACCGG-TTT -||-||- AC-GGGTT- >>> print(alignment[:, 3:]) ACCGG-TTT ||-||- ACGGGTT- >>> print(alignment[:, 3:-1]) ACCGG-TTT ||-|| ACGGGTT
- __format__(self, format_spec)¶
Return the alignment as a string in the specified file format.
Wrapper for self.format() .
- format(self, fmt='', **kwargs)¶
Return the alignment as a string in the specified file format.
- Arguments:
- fmt - File format. Acceptable values are
- “”create a human-readable representation of the
alignment (default);
- “BED”: create a line representing the alignment in
the Browser Extensible Data (BED) file format;
- “PSL”: create a line representing the alignment in
the Pattern Space Layout (PSL) file format as generated by BLAT;
- “SAM”: create a line representing the alignment in
the Sequence Alignment/Map (SAM) format.
- 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’.
- __str__(self)¶
Return a string representation of the PairwiseAlignment object.
Wrapper for self.format() .
- __len__(self)¶
Return the number of sequences in the alignment, which is always 2.
- 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) -GACCT-G -||--|-| CGA--TCG >>> len(alignment) 2 >>> alignment.shape (2, 8) >>> aligner.mode = "local" >>> alignments = aligner.align("GACCTG", "CGATCG") >>> alignment = alignments[0] >>> print(alignment) GACCT-G ||--|-| CGA--TCG >>> 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) GAACT ||--| GA--T >>> alignment.aligned (((0, 2), (4, 5)), ((0, 2), (2, 3))) >>> alignment = alignments[1] >>> print(alignment) GAACT |-|-| G-A-T >>> alignment.aligned (((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]) AAAC-AAA |||--||| AAA-GAAA >>> alignments[0].aligned (((0, 3), (4, 7)), ((0, 3), (4, 7))) >>> print(alignments[1]) AAA-CAAA |||--||| AAAG-AAA >>> alignments[1].aligned (((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.
- sort(self, 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) AATAA ||.|| AAGAA >>> alignment.sort() >>> print(alignment) AAGAA ||.|| AATAA
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 >>> alignment.sort(key=GC) >>> print(alignment) AATAA ||.|| AAGAA
You can reverse the sort order by passing reverse=True:
>>> alignment.sort(key=GC, reverse=True) >>> print(alignment) AAGAA ||.|| AATAA
The sequences are now sorted by decreasing GC content value.
- map(self, 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) AAAAAAAACCCCCCCAAAAAAAAAAAGGGGGGAAAAAAAA |||||||-----------|||||| CCCCCCC-----------GGGGGG >>> sequence = "CCCCGGGG" >>> alignments2 = aligner.align(transcript, sequence) >>> len(alignments2) 1 >>> alignment2 = alignments2[0] >>> print(alignment2) CCCCCCCGGGGGG |||||||| CCCCGGGG >>> alignment = alignment1.map(alignment2) >>> print(alignment) AAAAAAAACCCCCCCAAAAAAAAAAAGGGGGGAAAAAAAA ||||-----------|||| CCCC-----------GGGG >>> 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'
- 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) ATACTTACCTGGCAGGGGAGATACCATGATCACGAAGGTGGTTTTCCCAGGGCGAGGCTTATCCATTGCACTCCGGATGTGCTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACTGCATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTTTCCCCTG |||||||||||.||||||||..|||||||||||..|||||||..|||||||||||||||..|||||||||||.|||..|.|.|||||||||||||||||||||||||||||||||||||||.||||||||||||||||||||||||||||||||||.|||||.| ATACTTACCTGACAGGGGAGGCACCATGATCACACAGGTGGTCCTCCCAGGGCGAGGCTCTTCCATTGCACTGCGGGAGGGTTGACCCCTGCGATTTCCCCAAATGTGGGAAACTCGACTGTATAATTTGTGGTAGTGGGGGACTGCGTTCGCGCTATCCCCCG >>> 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.
- __hash__ = None¶
- class Bio.Align.PairwiseAlignments(seqA, seqB, score, paths)¶
Bases:
object
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__(self, 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 an PairwiseAlignments object by calling aligner.align(seqA, seqB), where aligner is a PairwiseAligner object.
- __len__(self)¶
Return the number of alignments.
- __getitem__(self, index)¶
- __iter__(self)¶
- __next__(self)¶
- class Bio.Align.PairwiseAligner(**kwargs)¶
Bases:
_algorithms.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: TACCG -|-|| -A-CG Score = 3.0: TACCG -||-| -AC-G
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: TACCG |-|| A-CG Score = 3.0: TACCG ||-| AC-G
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: TACCG -||-| -AC-G Score = 6.0: TACCG -|-|| -A-CG
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: TACCG -|-|| -A-CG Score = 5.5: TACCG -||-| -AC-G
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) KEVLA -|||- -EVL-
You can also set the value of attributes directly during construction of the PairwiseAligner object by providing them as keyword arguemnts:
>>> 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: TACCG -||-| -AC-G Score = 6.0: TACCG -|-|| -A-CG
- __init__(self, **kwargs)¶
Initialize a new PairwiseAligner with the keyword arguments as attributes.
Loops over the keyword arguments and sets them as attributes on the object.
- __setattr__(self, key, value)¶
Implement setattr(self, name, value).
- align(self, seqA, seqB, strand='+')¶
Return the alignments of two sequences using PairwiseAligner.
- score(self, seqA, seqB, strand='+')¶
Return the alignments score of two sequences using PairwiseAligner.
- __getstate__(self)¶
- __setstate__(self, state)¶