Bio.codonalign package

Module contents

Code for dealing with Codon Alignments.

Bio.codonalign.build(pro_align, nucl_seqs, corr_dict=None, gap_char='-', unknown='X', codon_table=NCBICodonTable(id=1, names=['Standard', 'SGC0'], ...), alphabet=None, complete_protein=False, anchor_len=10, max_score=10)

Build a codon alignment from protein alignment and corresponding nucleotides.

Arguments:
  • pro_align - a protein MultipleSeqAlignment object

  • nucl_seqs - an object returned by SeqIO.parse or SeqIO.index or a collection of SeqRecord.

  • alphabet - alphabet for the returned codon alignment

  • corr_dict - a dict that maps protein id to nucleotide id

  • complete_protein - whether the sequence begins with a start codon

Return a CodonAlignment object.

The example below answers this Biostars question: https://www.biostars.org/p/89741/

>>> from Bio.Alphabet import generic_dna, generic_protein
>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> from Bio.Align import MultipleSeqAlignment
>>> from Bio.codonalign import build
>>> seq1 = SeqRecord(Seq('ATGTCTCGT', alphabet=generic_dna), id='pro1')
>>> seq2 = SeqRecord(Seq('ATGCGT', alphabet=generic_dna), id='pro2')
>>> pro1 = SeqRecord(Seq('MSR', alphabet=generic_protein), id='pro1')
>>> pro2 = SeqRecord(Seq('M-R', alphabet=generic_protein), id='pro2')
>>> aln = MultipleSeqAlignment([pro1, pro2])
>>> codon_aln = build(aln, [seq1, seq2])
>>> print(codon_aln)
CodonAlphabet(Standard) CodonAlignment with 2 rows and 9 columns (3 codons)
ATGTCTCGT pro1
ATG---CGT pro2