Package Bio :: Package Align :: Module Generic :: Class Alignment
[hide private]
[frames] | no frames]

Class Alignment

source code

object --+
         |
        Alignment
Known Subclasses:

Represent a set of alignments (DEPRECATED).

This is a base class to represent alignments, which can be subclassed to deal with an alignment in a specific format.

With the introduction of the MultipleSeqAlignment class in Bio.Align, this base class is deprecated and is likely to be removed in future releases of Biopython.

Instance Methods [hide private]
 
__init__(self, alphabet)
Initialize a new Alignment object.
source code
 
_str_line(self, record, length=50)
Returns a truncated string representation of a SeqRecord (PRIVATE).
source code
 
__str__(self)
Returns a multi-line string summary of the alignment.
source code
 
__repr__(self)
Returns a representation of the object for debugging.
source code
 
format(self, format)
Returns the alignment as a string in the specified file format.
source code
 
__format__(self, format_spec)
Returns the alignment as a string in the specified file format.
source code
 
get_all_seqs(self)
Return all of the sequences involved in the alignment (DEPRECATED).
source code
 
__iter__(self)
Iterate over alignment rows as SeqRecord objects.
source code
 
get_seq_by_num(self, number)
Retrieve a sequence by row number (DEPRECATED).
source code
 
__len__(self)
Returns the number of sequences in the alignment.
source code
 
get_alignment_length(self)
Return the maximum length of the alignment.
source code
 
add_sequence(self, descriptor, sequence, start=None, end=None, weight=1.0)
Add a sequence to the alignment.
source code
 
get_column(self, col)
Returns a string containing a given column.
source code
 
__getitem__(self, index)
Access part of the alignment.
source code

Inherited from object: __delattr__, __getattribute__, __hash__, __new__, __reduce__, __reduce_ex__, __setattr__, __sizeof__, __subclasshook__

Properties [hide private]

Inherited from object: __class__

Method Details [hide private]

__init__(self, alphabet)
(Constructor)

source code 

Initialize a new Alignment object.

Arguments:

  • alphabet - The alphabet to use for the sequence objects that are created. This alphabet must be a gapped type.

e.g.

>>> from Bio.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> print(align)
Gapped(IUPACUnambiguousDNA(), '-') alignment with 3 rows and 12 columns
ACTGCTAGCTAG Alpha
ACT-CTAGCTAG Beta
ACTGCTAGATAG Gamma
Overrides: object.__init__

_str_line(self, record, length=50)

source code 

Returns a truncated string representation of a SeqRecord (PRIVATE).

This is a PRIVATE function used by the __str__ method.

__str__(self)
(Informal representation operator)

source code 

Returns 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.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> print(align)
Gapped(IUPACUnambiguousDNA(), '-') alignment with 3 rows and 12 columns
ACTGCTAGCTAG Alpha
ACT-CTAGCTAG Beta
ACTGCTAGATAG Gamma

See also the alignment's format method.

Overrides: object.__str__

__repr__(self)
(Representation operator)

source code 

Returns 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 ojects. For example:

<Bio.Align.Generic.Alignment instance (2 records of length 14, SingleLetterAlphabet()) 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.

Overrides: object.__repr__

format(self, format)

source code 

Returns 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.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> print(align.format("fasta"))
>Alpha
ACTGCTAGCTAG
>Beta
ACT-CTAGCTAG
>Gamma
ACTGCTAGATAG
<BLANKLINE>
>>> print(align.format("phylip"))
 3 12
Alpha      ACTGCTAGCT AG
Beta       ACT-CTAGCT AG
Gamma      ACTGCTAGAT AG
<BLANKLINE>

For Python 2.6, 3.0 or later see also the built in format() function.

__format__(self, format_spec)

source code 

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

This method supports the python format() function added in Python 2.6/3.0. The format_spec should be a lower case string supported by Bio.AlignIO as an output file format. See also the alignment's format() method.

Overrides: object.__format__

get_all_seqs(self)

source code 

Return all of the sequences involved in the alignment (DEPRECATED).

The return value is a list of SeqRecord objects.

This method is deprecated, as the Alignment object itself now offers much of the functionality of a list of SeqRecord objects (e.g. iteration or slicing to create a sub-alignment). Instead use the Python builtin function list, i.e. my_list = list(my_align)

__iter__(self)

source code 

Iterate over alignment rows as SeqRecord objects.

e.g.

>>> from Bio.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> for record in align:
...    print(record.id)
...    print(record.seq)
Alpha
ACTGCTAGCTAG
Beta
ACT-CTAGCTAG
Gamma
ACTGCTAGATAG

get_seq_by_num(self, number)

source code 

Retrieve a sequence by row number (DEPRECATED).

Returns:

  • A Seq object for the requested sequence.

Raises:

  • IndexError - If the specified number is out of range.

NOTE: This is a legacy method. In new code where you need to access the rows of the alignment (i.e. the sequences) consider iterating over them or accessing them as SeqRecord objects.

__len__(self)
(Length operator)

source code 

Returns 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)

source code 

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.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> align.get_alignment_length()
12

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

>>> len(align)
3

add_sequence(self, descriptor, sequence, start=None, end=None, weight=1.0)

source code 

Add a sequence to the alignment.

This doesn't do any kind of alignment, it just adds in the sequence object, which is assumed to be prealigned with the existing sequences.

Arguments:

  • descriptor - The descriptive id of the sequence being added. This will be used as the resulting SeqRecord's .id property (and, for historical compatibility, also the .description property)
  • sequence - A string with sequence info.
  • start - You can explicitly set the start point of the sequence. This is useful (at least) for BLAST alignments, which can just be partial alignments of sequences.
  • end - Specify the end of the sequence, which is important for the same reason as the start.
  • weight - The weight to place on the sequence in the alignment. By default, all sequences have the same weight. (0.0 => no weight, 1.0 => highest weight)

get_column(self, col)

source code 

Returns a string containing a given column.

e.g.

>>> from Bio.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha", "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",  "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma", "ACTGCTAGATAG")
>>> align.get_column(0)
'AAA'
>>> align.get_column(3)
'G-G'

__getitem__(self, index)
(Indexing operator)

source code 

Access part of the alignment.

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

>>> from Bio.Alphabet import IUPAC, Gapped
>>> align = Alignment(Gapped(IUPAC.unambiguous_dna, "-"))
>>> align.add_sequence("Alpha",  "ACTGCTAGCTAG")
>>> align.add_sequence("Beta",   "ACT-CTAGCTAG")
>>> align.add_sequence("Gamma",  "ACTGCTAGATAG")
>>> align.add_sequence("Delta",  "ACTGCTTGCTAG")
>>> align.add_sequence("Epsilon", "ACTGCTTGATAG")

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 ACTGCTAGCTAG
>>> last_record = align[-1]
>>> print("%s %s" % (last_record.id, last_record.seq))
Epsilon ACTGCTTGATAG

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)
Gapped(IUPACUnambiguousDNA(), '-') alignment with 3 rows and 12 columns
ACTGCTAGATAG Gamma
ACTGCTTGCTAG Delta
ACTGCTTGATAG 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)
Gapped(IUPACUnambiguousDNA(), '-') alignment with 3 rows and 12 columns
ACTGCTAGCTAG Alpha
ACTGCTAGATAG Gamma
ACTGCTTGATAG Epsilon

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

>>> rev_alignment = align[::-1]
>>> print(rev_alignment)
Gapped(IUPACUnambiguousDNA(), '-') alignment with 5 rows and 12 columns
ACTGCTTGATAG Epsilon
ACTGCTTGCTAG Delta
ACTGCTAGATAG Gamma
ACT-CTAGCTAG Beta
ACTGCTAGCTAG Alpha

Right now, these are the ONLY indexing operations supported. The use of a second column based index is under discussion for a future update.