Bio.SeqUtils package

Module contents

Miscellaneous functions for dealing with sequences.

Bio.SeqUtils.GC(seq)

Calculate G+C content, return percentage (as float between 0 and 100).

Copes mixed case sequences, and with the ambiguous nucleotide S (G or C) when counting the G and C content. The percentage is calculated against the full length, e.g.:

>>> from Bio.SeqUtils import GC
>>> GC("ACTGN")
40.0

Note that this will return zero for an empty sequence.

Bio.SeqUtils.GC123(seq)

Calculate G+C content: total, for first, second and third positions.

Returns a tuple of four floats (percentages between 0 and 100) for the entire sequence, and the three codon positions. e.g.

>>> from Bio.SeqUtils import GC123
>>> GC123("ACTGTN")
(40.0, 50.0, 50.0, 0.0)

Copes with mixed case sequences, but does NOT deal with ambiguous nucleotides.

Bio.SeqUtils.GC_skew(seq, window=100)

Calculate GC skew (G-C)/(G+C) for multiple windows along the sequence.

Returns a list of ratios (floats), controlled by the length of the sequence and the size of the window.

Returns 0 for windows without any G/C by handling zero division errors.

Does NOT look at any ambiguous nucleotides.

Bio.SeqUtils.xGC_skew(seq, window=1000, zoom=100, r=300, px=100, py=100)

Calculate and plot normal and accumulated GC skew (GRAPHICS !!!).

Search for a DNA subseq in sequence, return list of [subseq, positions].

Use ambiguous values (like N = A or T or C or G, R = A or G etc.), searches only on forward strand.

Bio.SeqUtils.seq3(seq, custom_map=None, undef_code='Xaa')

Convert protein sequence from one-letter to three-letter code.

The single required input argument ‘seq’ should be a protein sequence using single letter codes, either as a Python string or as a Seq or MutableSeq object.

This function returns the amino acid sequence as a string using the three letter amino acid codes. Output follows the IUPAC standard (including ambiguous characters B for “Asx”, J for “Xle” and X for “Xaa”, and also U for “Sel” and O for “Pyl”) plus “Ter” for a terminator given as an asterisk. Any unknown character (including possible gap characters), is changed into ‘Xaa’ by default.

e.g.

>>> from Bio.SeqUtils import seq3
>>> seq3("MAIVMGRWKGAR*")
'MetAlaIleValMetGlyArgTrpLysGlyAlaArgTer'

You can set a custom translation of the codon termination code using the dictionary “custom_map” argument (which defaults to {‘*’: ‘Ter’}), e.g.

>>> seq3("MAIVMGRWKGAR*", custom_map={"*": "***"})
'MetAlaIleValMetGlyArgTrpLysGlyAlaArg***'

You can also set a custom translation for non-amino acid characters, such as ‘-‘, using the “undef_code” argument, e.g.

>>> seq3("MAIVMGRWKGA--R*", undef_code='---')
'MetAlaIleValMetGlyArgTrpLysGlyAla------ArgTer'

If not given, “undef_code” defaults to “Xaa”, e.g.

>>> seq3("MAIVMGRWKGA--R*")
'MetAlaIleValMetGlyArgTrpLysGlyAlaXaaXaaArgTer'

This function was inspired by BioPerl’s seq3.

Bio.SeqUtils.seq1(seq, custom_map=None, undef_code='X')

Convert protein sequence from three-letter to one-letter code.

The single required input argument ‘seq’ should be a protein sequence using three-letter codes, either as a Python string or as a Seq or MutableSeq object.

This function returns the amino acid sequence as a string using the one letter amino acid codes. Output follows the IUPAC standard (including ambiguous characters “B” for “Asx”, “J” for “Xle”, “X” for “Xaa”, “U” for “Sel”, and “O” for “Pyl”) plus “*” for a terminator given the “Ter” code. Any unknown character (including possible gap characters), is changed into ‘-‘ by default.

e.g.

>>> from Bio.SeqUtils import seq1
>>> seq1("MetAlaIleValMetGlyArgTrpLysGlyAlaArgTer")
'MAIVMGRWKGAR*'

The input is case insensitive, e.g.

>>> from Bio.SeqUtils import seq1
>>> seq1("METalaIlEValMetGLYArgtRplysGlyAlaARGTer")
'MAIVMGRWKGAR*'

You can set a custom translation of the codon termination code using the dictionary “custom_map” argument (defaulting to {‘Ter’: ‘*’}), e.g.

>>> seq1("MetAlaIleValMetGlyArgTrpLysGlyAla***", custom_map={"***": "*"})
'MAIVMGRWKGA*'

You can also set a custom translation for non-amino acid characters, such as ‘-‘, using the “undef_code” argument, e.g.

>>> seq1("MetAlaIleValMetGlyArgTrpLysGlyAla------ArgTer", undef_code='?')
'MAIVMGRWKGA??R*'

If not given, “undef_code” defaults to “X”, e.g.

>>> seq1("MetAlaIleValMetGlyArgTrpLysGlyAla------ArgTer")
'MAIVMGRWKGAXXR*'
Bio.SeqUtils.molecular_weight(seq, seq_type=None, double_stranded=False, circular=False, monoisotopic=False)

Calculate the molecular mass of DNA, RNA or protein sequences as float.

Only unambiguous letters are allowed. Nucleotide sequences are assumed to have a 5’ phosphate.

Arguments:
  • seq: String or Biopython sequence object.

  • seq_type: The default (None) is to take the alphabet from the seq argument, or assume DNA if the seq argument is a string. Override this with a string ‘DNA’, ‘RNA’, or ‘protein’.

  • double_stranded: Calculate the mass for the double stranded molecule?

  • circular: Is the molecule circular (has no ends)?

  • monoisotopic: Use the monoisotopic mass tables?

Note that for backwards compatibility, if the seq argument is a string, or Seq object with a generic alphabet, and no seq_type is specified (i.e. left as None), then DNA is assumed.

>>> print("%0.2f" % molecular_weight("AGC"))
949.61
>>> print("%0.2f" % molecular_weight(Seq("AGC")))
949.61

However, it is better to be explicit - for example with strings:

>>> print("%0.2f" % molecular_weight("AGC", "DNA"))
949.61
>>> print("%0.2f" % molecular_weight("AGC", "RNA"))
997.61
>>> print("%0.2f" % molecular_weight("AGC", "protein"))
249.29

Or, with the sequence alphabet:

>>> from Bio.Seq import Seq
>>> from Bio.Alphabet import generic_dna, generic_rna, generic_protein
>>> print("%0.2f" % molecular_weight(Seq("AGC", generic_dna)))
949.61
>>> print("%0.2f" % molecular_weight(Seq("AGC", generic_rna)))
997.61
>>> print("%0.2f" % molecular_weight(Seq("AGC", generic_protein)))
249.29

Also note that contradictory sequence alphabets and seq_type will also give an exception:

>>> from Bio.Seq import Seq
>>> from Bio.Alphabet import generic_dna
>>> print("%0.2f" % molecular_weight(Seq("AGC", generic_dna), "RNA"))
Traceback (most recent call last):
  ...
ValueError: seq_type='RNA' contradicts DNA from seq alphabet
Bio.SeqUtils.six_frame_translations(seq, genetic_code=1)

Return pretty string showing the 6 frame translations and GC content.

Nice looking 6 frame translation with GC content - code from xbbtools similar to DNA Striders six-frame translation

>>> from Bio.SeqUtils import six_frame_translations
>>> print(six_frame_translations("AUGGCCAUUGUAAUGGGCCGCUGA"))
GC_Frame: a:5 t:0 g:8 c:5 
Sequence: auggccauug ... gggccgcuga, 24 nt, 54.17 %GC


1/1
  G  H  C  N  G  P  L
 W  P  L  *  W  A  A
M  A  I  V  M  G  R  *
auggccauuguaaugggccgcuga   54 %
uaccgguaacauuacccggcgacu
A  M  T  I  P  R  Q 
 H  G  N  Y  H  A  A  S
  P  W  Q  L  P  G  S