biopython v1.71.0 Bio.Seq.Seq

Read-only sequence object (essentially a string with an alphabet).

Like normal python strings, our basic sequence object is immutable. This prevents you from doing my_seq[5] = “A” for example, but does allow Seq objects to be used as dictionary keys.

The Seq object provides a number of string like methods (such as count, find, split and strip), which are alphabet aware where appropriate.

In addition to the string like sequence, the Seq object has an alphabet property. This is an instance of an Alphabet class from Bio.Alphabet, for example generic DNA, or IUPAC DNA. This describes the type of molecule (e.g. RNA, DNA, protein) and may also indicate the expected symbols (letters).

The Seq object also provides some biological methods, such as complement, reverse_complement, transcribe, back_transcribe and translate (which are not applicable to sequences with a protein alphabet).

Link to this section Summary

Functions

Add another sequence or string to this sequence

Implement the ‘in’ keyword, like a python string

Compare the sequence to another sequence or a string (README)

Return a subsequence of single letter, use my_seq[index]

Hash for comparison

Create a Seq object

Implement the less-than or equal operand

Return the length of the sequence, use len(my_seq)

Implement the less-than operand

Implement the not-equal operand

Add a sequence on the left

Return (truncated) representation of the sequence for debugging

Return the full sequence as a python string, use str(my_seq)

Convert string/Seq/MutableSeq to string, checking alphabet (PRIVATE)

Return the DNA sequence from an RNA sequence by creating a new Seq object

Return the complement sequence by creating a new Seq object

Return a non-overlapping count, like that of a python string

Return an overlapping count

Return True if the Seq ends with the given suffix, False otherwise

Find method, like that of a python string

Return a lower case copy of the sequence

Return a new Seq object with leading (left) end stripped

Return the reverse complement sequence by creating a new Seq object

Find from right method, like that of a python string

Do a right split method, like that of a python string

Return a new Seq object with trailing (right) end stripped

Split method, like that of a python string

Return True if the Seq starts with the given prefix, False otherwise

Return a new Seq object with leading and trailing ends stripped

Return the full sequence as a MutableSeq object

Return the full sequence as a python string (DEPRECATED)

Return the RNA sequence from a DNA sequence by creating a new Seq object

Turn a nucleotide sequence into a protein sequence by creating a new Seq object

Return a copy of the sequence without the gap character(s)

Return an upper case copy of the sequence

Link to this section Functions

Add another sequence or string to this sequence.

If adding a string to a Seq, the alphabet is preserved:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_protein
 >>> Seq("MELKI", generic_protein) + "LV"
 Seq('MELKILV', ProteinAlphabet())

When adding two Seq (like) objects, the alphabets are important. Consider this example:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet.IUPAC import unambiguous_dna, ambiguous_dna
 >>> unamb_dna_seq = Seq("ACGT", unambiguous_dna)
 >>> ambig_dna_seq = Seq("ACRGT", ambiguous_dna)
 >>> unamb_dna_seq
 Seq('ACGT', IUPACUnambiguousDNA())
 >>> ambig_dna_seq
 Seq('ACRGT', IUPACAmbiguousDNA())

If we add the ambiguous and unambiguous IUPAC DNA alphabets, we get the more general ambiguous IUPAC DNA alphabet:

 >>> unamb_dna_seq + ambig_dna_seq
 Seq('ACGTACRGT', IUPACAmbiguousDNA())

However, if the default generic alphabet is included, the result is a generic alphabet:

 >>> Seq("") + ambig_dna_seq
 Seq('ACRGT', Alphabet())

You can’t add RNA and DNA sequences:

 >>> from Bio.Alphabet import generic_dna, generic_rna
 >>> Seq("ACGT", generic_dna) + Seq("ACGU", generic_rna)
 Traceback (most recent call last):
    ...
 TypeError: Incompatible alphabets DNAAlphabet() and RNAAlphabet()

You can’t add nucleotide and protein sequences:

 >>> from Bio.Alphabet import generic_dna, generic_protein
 >>> Seq("ACGT", generic_dna) + Seq("MELKI", generic_protein)
 Traceback (most recent call last):
    ...
 TypeError: Incompatible alphabets DNAAlphabet() and ProteinAlphabet()

Implement the ‘in’ keyword, like a python string.

e.g.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna, generic_rna, generic_protein
 >>> my_dna = Seq("ATATGAAATTTGAAAA", generic_dna)
 >>> "AAA" in my_dna
 True
 >>> Seq("AAA") in my_dna
 True
 >>> Seq("AAA", generic_dna) in my_dna
 True

Like other Seq methods, this will raise a type error if another Seq (or Seq like) object with an incompatible alphabet is used:

 >>> Seq("AAA", generic_rna) in my_dna
 Traceback (most recent call last):
    ...
 TypeError: Incompatible alphabets DNAAlphabet() and RNAAlphabet()
 >>> Seq("AAA", generic_protein) in my_dna
 Traceback (most recent call last):
    ...
 TypeError: Incompatible alphabets DNAAlphabet() and ProteinAlphabet()

Compare the sequence to another sequence or a string (README).

Historically comparing Seq objects has done Python object comparison. After considerable discussion (keeping in mind constraints of the Python language, hashes and dictionary support), Biopython now uses simple string comparison (with a warning about the change).

Note that incompatible alphabets (e.g. DNA to RNA) will trigger a warning.

During this transition period, please just do explicit comparisons:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> seq1 = Seq("ACGT")
 >>> seq2 = Seq("ACGT")
 >>> id(seq1) == id(seq2)
 False
 >>> str(seq1) == str(seq2)
 True

The new behaviour is to use string-like equality:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> seq1 == seq2
 True
 >>> seq1 == "ACGT"
 True
 >>> seq1 == Seq("ACGT", generic_dna)
 True

Return a subsequence of single letter, use my_seq[index].

 >>> my_seq = Seq('ACTCGACGTCG')
 >>> my_seq[5]
 'A'

Hash for comparison.

See the cmp documentation - this has changed from past versions of Biopython!

Create a Seq object.

Arguments:

  • seq - Sequence, required (string)
  • alphabet - Optional argument, an Alphabet object from Bio.Alphabet

You will typically use Bio.SeqIO to read in sequences from files as SeqRecord objects, whose sequence will be exposed as a Seq object via the seq property.

However, will often want to create your own Seq objects directly:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> my_seq = Seq("MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF",
 ...              IUPAC.protein)
 >>> my_seq
 Seq('MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF', IUPACProtein())
 >>> print(my_seq)
 MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF
 >>> my_seq.alphabet
 IUPACProtein()

Implement the less-than or equal operand.

Return the length of the sequence, use len(my_seq).

Implement the less-than operand.

Implement the not-equal operand.

Add a sequence on the left.

If adding a string to a Seq, the alphabet is preserved:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_protein
 >>> "LV" + Seq("MELKI", generic_protein)
 Seq('LVMELKI', ProteinAlphabet())

Adding two Seq (like) objects is handled via the add method.

Return (truncated) representation of the sequence for debugging.

Return the full sequence as a python string, use str(my_seq).

Note that Biopython 1.44 and earlier would give a truncated version of repr(my_seq) for str(my_seq). If you are writing code which need to be backwards compatible with old Biopython, you should continue to use my_seq.tostring() rather than str(my_seq).

Link to this function _get_seq_str_and_check_alphabet()

Convert string/Seq/MutableSeq to string, checking alphabet (PRIVATE).

For a string argument, returns the string.

For a Seq or MutableSeq, it checks the alphabet is compatible (raising an exception if it isn’t), and then returns a string.

Link to this function back_transcribe()

Return the DNA sequence from an RNA sequence by creating a new Seq object.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> messenger_rna = Seq("AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG",
 ...                     IUPAC.unambiguous_rna)
 >>> messenger_rna
 Seq('AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG', IUPACUnambiguousRNA())
 >>> messenger_rna.back_transcribe()
 Seq('ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG', IUPACUnambiguousDNA())

Trying to back-transcribe a protein or DNA sequence raises an exception:

 >>> my_protein = Seq("MAIVMGR", IUPAC.protein)
 >>> my_protein.back_transcribe()
 Traceback (most recent call last):
    ...
 ValueError: Proteins cannot be back transcribed!

Return the complement sequence by creating a new Seq object.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> my_dna = Seq("CCCCCGATAG", IUPAC.unambiguous_dna)
 >>> my_dna
 Seq('CCCCCGATAG', IUPACUnambiguousDNA())
 >>> my_dna.complement()
 Seq('GGGGGCTATC', IUPACUnambiguousDNA())

You can of course used mixed case sequences,

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> my_dna = Seq("CCCCCgatA-GD", generic_dna)
 >>> my_dna
 Seq('CCCCCgatA-GD', DNAAlphabet())
 >>> my_dna.complement()
 Seq('GGGGGctaT-CH', DNAAlphabet())

Note in the above example, ambiguous character D denotes G, A or T so its complement is H (for C, T or A).

Trying to complement a protein sequence raises an exception.

 >>> my_protein = Seq("MAIVMGR", IUPAC.protein)
 >>> my_protein.complement()
 Traceback (most recent call last):
    ...
 ValueError: Proteins do not have complements!

Return a non-overlapping count, like that of a python string.

This behaves like the python string method of the same name, which does a non-overlapping count!

For an overlapping search use the newer count_overlap() method.

Returns an integer, the number of occurrences of substring argument sub in the (sub)sequence given by [start:end]. Optional arguments start and end are interpreted as in slice notation.

Arguments:

  • sub - a string or another Seq object to look for
  • start - optional integer, slice start
  • end - optional integer, slice end

e.g.

 >>> from Bio.Seq import Seq
 >>> my_seq = Seq("AAAATGA")
 >>> print(my_seq.count("A"))
 5
 >>> print(my_seq.count("ATG"))
 1
 >>> print(my_seq.count(Seq("AT")))
 1
 >>> print(my_seq.count("AT", 2, -1))
 1

HOWEVER, please note because python strings and Seq objects (and MutableSeq objects) do a non-overlapping search, this may not give the answer you expect:

 >>> "AAAA".count("AA")
 2
 >>> print(Seq("AAAA").count("AA"))
 2

An overlapping search, as implemented in .count_overlap(), would give the answer as three!

Link to this function count_overlap()

Return an overlapping count.

For a non-overlapping search use the count() method.

Returns an integer, the number of occurrences of substring argument sub in the (sub)sequence given by [start:end]. Optional arguments start and end are interpreted as in slice notation.

Arguments:

  • sub - a string or another Seq object to look for
  • start - optional integer, slice start
  • end - optional integer, slice end

e.g.

 >>> from Bio.Seq import Seq
 >>> print(Seq("AAAA").count_overlap("AA"))
 3
 >>> print(Seq("ATATATATA").count_overlap("ATA"))
 4
 >>> print(Seq("ATATATATA").count_overlap("ATA", 3, -1))
 1

Where substrings do not overlap, should behave the same as the count() method:

 >>> from Bio.Seq import Seq
 >>> my_seq = Seq("AAAATGA")
 >>> print(my_seq.count_overlap("A"))
 5
 >>> my_seq.count_overlap("A") == my_seq.count("A")
 True
 >>> print(my_seq.count_overlap("ATG"))
 1
 >>> my_seq.count_overlap("ATG") == my_seq.count("ATG")
 True
 >>> print(my_seq.count_overlap(Seq("AT")))
 1
 >>> my_seq.count_overlap(Seq("AT")) == my_seq.count(Seq("AT"))
 True
 >>> print(my_seq.count_overlap("AT", 2, -1))
 1
 >>> my_seq.count_overlap("AT", 2, -1) == my_seq.count("AT", 2, -1)
 True

HOWEVER, do not use this method for such cases because the count() method is much for efficient.

Return True if the Seq ends with the given suffix, False otherwise.

This behaves like the python string method of the same name.

Return True if the sequence ends with the specified suffix (a string or another Seq object), False otherwise. With optional start, test sequence beginning at that position. With optional end, stop comparing sequence at that position. suffix can also be a tuple of strings to try. e.g.

 >>> from Bio.Seq import Seq
 >>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
 >>> my_rna.endswith("UUG")
 True
 >>> my_rna.endswith("AUG")
 False
 >>> my_rna.endswith("AUG", 0, 18)
 True
 >>> my_rna.endswith(("UCC", "UCA", "UUG"))
 True

Find method, like that of a python string.

This behaves like the python string method of the same name.

Returns an integer, the index of the first occurrence of substring argument sub in the (sub)sequence given by [start:end].

Arguments:

  • sub - a string or another Seq object to look for
  • start - optional integer, slice start
  • end - optional integer, slice end

Returns -1 if the subsequence is NOT found.

e.g. Locating the first typical start codon, AUG, in an RNA sequence:

 >>> from Bio.Seq import Seq
 >>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
 >>> my_rna.find("AUG")
 3

Return a lower case copy of the sequence.

This will adjust the alphabet if required. Note that the IUPAC alphabets are upper case only, and thus a generic alphabet must be substituted.

 >>> from Bio.Alphabet import Gapped, generic_dna
 >>> from Bio.Alphabet import IUPAC
 >>> from Bio.Seq import Seq
 >>> my_seq = Seq("CGGTACGCTTATGTCACGTAG*AAAAAA",
 ...          Gapped(IUPAC.unambiguous_dna, "*"))
 >>> my_seq
 Seq('CGGTACGCTTATGTCACGTAG*AAAAAA', Gapped(IUPACUnambiguousDNA(), '*'))
 >>> my_seq.lower()
 Seq('cggtacgcttatgtcacgtag*aaaaaa', Gapped(DNAAlphabet(), '*'))

See also the upper method.

Return a new Seq object with leading (left) end stripped.

This behaves like the python string method of the same name.

Optional argument chars defines which characters to remove. If omitted or None (default) then as for the python string method, this defaults to removing any white space.

e.g. print(my_seq.lstrip(“-“))

See also the strip and rstrip methods.

Link to this function reverse_complement()

Return the reverse complement sequence by creating a new Seq object.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> my_dna = Seq("CCCCCGATAGNR", IUPAC.ambiguous_dna)
 >>> my_dna
 Seq('CCCCCGATAGNR', IUPACAmbiguousDNA())
 >>> my_dna.reverse_complement()
 Seq('YNCTATCGGGGG', IUPACAmbiguousDNA())

Note in the above example, since R = G or A, its complement is Y (which denotes C or T).

You can of course used mixed case sequences,

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> my_dna = Seq("CCCCCgatA-G", generic_dna)
 >>> my_dna
 Seq('CCCCCgatA-G', DNAAlphabet())
 >>> my_dna.reverse_complement()
 Seq('C-TatcGGGGG', DNAAlphabet())

Trying to complement a protein sequence raises an exception:

 >>> my_protein = Seq("MAIVMGR", IUPAC.protein)
 >>> my_protein.reverse_complement()
 Traceback (most recent call last):
    ...
 ValueError: Proteins do not have complements!

Find from right method, like that of a python string.

This behaves like the python string method of the same name.

Returns an integer, the index of the last (right most) occurrence of substring argument sub in the (sub)sequence given by [start:end].

Arguments:

  • sub - a string or another Seq object to look for
  • start - optional integer, slice start
  • end - optional integer, slice end

Returns -1 if the subsequence is NOT found.

e.g. Locating the last typical start codon, AUG, in an RNA sequence:

 >>> from Bio.Seq import Seq
 >>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
 >>> my_rna.rfind("AUG")
 15

Do a right split method, like that of a python string.

This behaves like the python string method of the same name.

Return a list of the ‘words’ in the string (as Seq objects), using sep as the delimiter string. If maxsplit is given, at most maxsplit splits are done COUNTING FROM THE RIGHT. If maxsplit is omitted, all splits are made.

Following the python string method, sep will by default be any white space (tabs, spaces, newlines) but this is unlikely to apply to biological sequences.

e.g. print(my_seq.rsplit(“*”,1))

See also the split method.

Return a new Seq object with trailing (right) end stripped.

This behaves like the python string method of the same name.

Optional argument chars defines which characters to remove. If omitted or None (default) then as for the python string method, this defaults to removing any white space.

e.g. Removing a nucleotide sequence’s polyadenylation (poly-A tail):

 >>> from Bio.Alphabet import IUPAC
 >>> from Bio.Seq import Seq
 >>> my_seq = Seq("CGGTACGCTTATGTCACGTAGAAAAAA", IUPAC.unambiguous_dna)
 >>> my_seq
 Seq('CGGTACGCTTATGTCACGTAGAAAAAA', IUPACUnambiguousDNA())
 >>> my_seq.rstrip("A")
 Seq('CGGTACGCTTATGTCACGTAG', IUPACUnambiguousDNA())

See also the strip and lstrip methods.

Split method, like that of a python string.

This behaves like the python string method of the same name.

Return a list of the ‘words’ in the string (as Seq objects), using sep as the delimiter string. If maxsplit is given, at most maxsplit splits are done. If maxsplit is omitted, all splits are made.

Following the python string method, sep will by default be any white space (tabs, spaces, newlines) but this is unlikely to apply to biological sequences.

e.g.

 >>> from Bio.Seq import Seq
 >>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
 >>> my_aa = my_rna.translate()
 >>> my_aa
 Seq('VMAIVMGR*KGAR*L', HasStopCodon(ExtendedIUPACProtein(), '*'))
 >>> for pep in my_aa.split("*"):
 ...     pep
 Seq('VMAIVMGR', HasStopCodon(ExtendedIUPACProtein(), '*'))
 Seq('KGAR', HasStopCodon(ExtendedIUPACProtein(), '*'))
 Seq('L', HasStopCodon(ExtendedIUPACProtein(), '*'))
 >>> for pep in my_aa.split("*", 1):
 ...     pep
 Seq('VMAIVMGR', HasStopCodon(ExtendedIUPACProtein(), '*'))
 Seq('KGAR*L', HasStopCodon(ExtendedIUPACProtein(), '*'))

See also the rsplit method:

 >>> for pep in my_aa.rsplit("*", 1):
 ...     pep
 Seq('VMAIVMGR*KGAR', HasStopCodon(ExtendedIUPACProtein(), '*'))
 Seq('L', HasStopCodon(ExtendedIUPACProtein(), '*'))

Return True if the Seq starts with the given prefix, False otherwise.

This behaves like the python string method of the same name.

Return True if the sequence starts with the specified prefix (a string or another Seq object), False otherwise. With optional start, test sequence beginning at that position. With optional end, stop comparing sequence at that position. prefix can also be a tuple of strings to try. e.g.

 >>> from Bio.Seq import Seq
 >>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
 >>> my_rna.startswith("GUC")
 True
 >>> my_rna.startswith("AUG")
 False
 >>> my_rna.startswith("AUG", 3)
 True
 >>> my_rna.startswith(("UCC", "UCA", "UCG"), 1)
 True

Return a new Seq object with leading and trailing ends stripped.

This behaves like the python string method of the same name.

Optional argument chars defines which characters to remove. If omitted or None (default) then as for the python string method, this defaults to removing any white space.

e.g. print(my_seq.strip(“-“))

See also the lstrip and rstrip methods.

Return the full sequence as a MutableSeq object.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> my_seq = Seq("MKQHKAMIVALIVICITAVVAAL",
 ...              IUPAC.protein)
 >>> my_seq
 Seq('MKQHKAMIVALIVICITAVVAAL', IUPACProtein())
 >>> my_seq.tomutable()
 MutableSeq('MKQHKAMIVALIVICITAVVAAL', IUPACProtein())

Note that the alphabet is preserved.

Return the full sequence as a python string (DEPRECATED).

You are now encouraged to use str(my_seq) instead of my_seq.tostring().

Return the RNA sequence from a DNA sequence by creating a new Seq object.

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC
 >>> coding_dna = Seq("ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG",
 ...                  IUPAC.unambiguous_dna)
 >>> coding_dna
 Seq('ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG', IUPACUnambiguousDNA())
 >>> coding_dna.transcribe()
 Seq('AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG', IUPACUnambiguousRNA())

Trying to transcribe a protein or RNA sequence raises an exception:

 >>> my_protein = Seq("MAIVMGR", IUPAC.protein)
 >>> my_protein.transcribe()
 Traceback (most recent call last):
    ...
 ValueError: Proteins cannot be transcribed!

Turn a nucleotide sequence into a protein sequence by creating a new Seq object.

This method will translate DNA or RNA sequences, and those with a nucleotide or generic alphabet. Trying to translate a protein sequence raises an exception.

Arguments:

  • table - Which codon table to use? This can be either a name (string), an NCBI identifier (integer), or a CodonTable object (useful for non-standard genetic codes). This defaults to the “Standard” table.
  • stop_symbol - Single character string, what to use for terminators. This defaults to the asterisk, “*”.
  • to_stop - Boolean, defaults to False meaning do a full translation continuing on past any stop codons (translated as the specified stop_symbol). If True, translation is terminated at the first in frame stop codon (and the stop_symbol is not appended to the returned protein sequence).
  • cds - Boolean, indicates this is a complete CDS. If True, this checks the sequence starts with a valid alternative start codon (which will be translated as methionine, M), that the sequence length is a multiple of three, and that there is a single in frame stop codon at the end (this will be excluded from the protein sequence, regardless of the to_stop option). If these tests fail, an exception is raised.
  • gap - Single character string to denote symbol used for gaps. It will try to guess the gap character from the alphabet.

e.g. Using the standard table:

 >>> coding_dna = Seq("GTGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG")
 >>> coding_dna.translate()
 Seq('VAIVMGR*KGAR*', HasStopCodon(ExtendedIUPACProtein(), '*'))
 >>> coding_dna.translate(stop_symbol="@")
 Seq('VAIVMGR@KGAR@', HasStopCodon(ExtendedIUPACProtein(), '@'))
 >>> coding_dna.translate(to_stop=True)
 Seq('VAIVMGR', ExtendedIUPACProtein())

Now using NCBI table 2, where TGA is not a stop codon:

 >>> coding_dna.translate(table=2)
 Seq('VAIVMGRWKGAR*', HasStopCodon(ExtendedIUPACProtein(), '*'))
 >>> coding_dna.translate(table=2, to_stop=True)
 Seq('VAIVMGRWKGAR', ExtendedIUPACProtein())

In fact, GTG is an alternative start codon under NCBI table 2, meaning this sequence could be a complete CDS:

 >>> coding_dna.translate(table=2, cds=True)
 Seq('MAIVMGRWKGAR', ExtendedIUPACProtein())

It isn’t a valid CDS under NCBI table 1, due to both the start codon and also the in frame stop codons:

 >>> coding_dna.translate(table=1, cds=True)
 Traceback (most recent call last):
     ...
 TranslationError: First codon 'GTG' is not a start codon

If the sequence has no in-frame stop codon, then the to_stop argument has no effect:

 >>> coding_dna2 = Seq("TTGGCCATTGTAATGGGCCGC")
 >>> coding_dna2.translate()
 Seq('LAIVMGR', ExtendedIUPACProtein())
 >>> coding_dna2.translate(to_stop=True)
 Seq('LAIVMGR', ExtendedIUPACProtein())

When translating gapped sequences, the gap character is inferred from the alphabet:

 >>> from Bio.Alphabet import Gapped
 >>> coding_dna3 = Seq("GTG---GCCATT", Gapped(IUPAC.unambiguous_dna))
 >>> coding_dna3.translate()
 Seq('V-AI', Gapped(ExtendedIUPACProtein(), '-'))

It is possible to pass the gap character when the alphabet is missing:

 >>> coding_dna4 = Seq("GTG---GCCATT")
 >>> coding_dna4.translate(gap='-')
 Seq('V-AI', Gapped(ExtendedIUPACProtein(), '-'))

NOTE - Ambiguous codons like “TAN” or “NNN” could be an amino acid or a stop codon. These are translated as “X”. Any invalid codon (e.g. “TA?” or “T-A”) will throw a TranslationError.

NOTE - This does NOT behave like the python string’s translate method. For that use str(my_seq).translate(…) instead.

Return a copy of the sequence without the gap character(s).

The gap character can be specified in two ways - either as an explicit argument, or via the sequence’s alphabet. For example:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> my_dna = Seq("-ATA--TGAAAT-TTGAAAA", generic_dna)
 >>> my_dna
 Seq('-ATA--TGAAAT-TTGAAAA', DNAAlphabet())
 >>> my_dna.ungap("-")
 Seq('ATATGAAATTTGAAAA', DNAAlphabet())

If the gap character is not given as an argument, it will be taken from the sequence’s alphabet (if defined). Notice that the returned sequence’s alphabet is adjusted since it no longer requires a gapped alphabet:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC, Gapped, HasStopCodon
 >>> my_pro = Seq("MVVLE=AD*", HasStopCodon(Gapped(IUPAC.protein, "=")))
 >>> my_pro
 Seq('MVVLE=AD*', HasStopCodon(Gapped(IUPACProtein(), '='), '*'))
 >>> my_pro.ungap()
 Seq('MVVLEAD*', HasStopCodon(IUPACProtein(), '*'))

Or, with a simpler gapped DNA example:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import IUPAC, Gapped
 >>> my_seq = Seq("CGGGTAG=AAAAAA", Gapped(IUPAC.unambiguous_dna, "="))
 >>> my_seq
 Seq('CGGGTAG=AAAAAA', Gapped(IUPACUnambiguousDNA(), '='))
 >>> my_seq.ungap()
 Seq('CGGGTAGAAAAAA', IUPACUnambiguousDNA())

As long as it is consistent with the alphabet, although it is redundant, you can still supply the gap character as an argument to this method:

 >>> my_seq
 Seq('CGGGTAG=AAAAAA', Gapped(IUPACUnambiguousDNA(), '='))
 >>> my_seq.ungap("=")
 Seq('CGGGTAGAAAAAA', IUPACUnambiguousDNA())

However, if the gap character given as the argument disagrees with that declared in the alphabet, an exception is raised:

 >>> my_seq
 Seq('CGGGTAG=AAAAAA', Gapped(IUPACUnambiguousDNA(), '='))
 >>> my_seq.ungap("-")
 Traceback (most recent call last):
    ...
 ValueError: Gap '-' does not match '=' from alphabet

Finally, if a gap character is not supplied, and the alphabet does not define one, an exception is raised:

 >>> from Bio.Seq import Seq
 >>> from Bio.Alphabet import generic_dna
 >>> my_dna = Seq("ATA--TGAAAT-TTGAAAA", generic_dna)
 >>> my_dna
 Seq('ATA--TGAAAT-TTGAAAA', DNAAlphabet())
 >>> my_dna.ungap()
 Traceback (most recent call last):
    ...
 ValueError: Gap character not given and not defined in alphabet

Return an upper case copy of the sequence.

 >>> from Bio.Alphabet import HasStopCodon, generic_protein
 >>> from Bio.Seq import Seq
 >>> my_seq = Seq("VHLTPeeK*", HasStopCodon(generic_protein))
 >>> my_seq
 Seq('VHLTPeeK*', HasStopCodon(ProteinAlphabet(), '*'))
 >>> my_seq.lower()
 Seq('vhltpeek*', HasStopCodon(ProteinAlphabet(), '*'))
 >>> my_seq.upper()
 Seq('VHLTPEEK*', HasStopCodon(ProteinAlphabet(), '*'))

This will adjust the alphabet if required. See also the lower method.