Bio.Seq module

Provide objects to represent biological sequences.

See also the Seq wiki and the chapter in our tutorial:
class Bio.Seq.Seq(data)

Bases: object

Read-only sequence object (essentially a string with biological methods).

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

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

__init__(self, data)

Create a Seq object.

Arguments:
  • data - Sequence, required (string)

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
>>> my_seq = Seq("MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF")
>>> my_seq
Seq('MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF')
>>> print(my_seq)
MKQHKAMIVALIVICITAVVAALVTRKDLCEVHIRTGQTEVAVF
__repr__(self)

Return (truncated) representation of the sequence for debugging.

__str__(self)

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 really old Biopython, you should continue to use my_seq.tostring() as follows:

try:
    # The old way, removed in Biopython 1.73
    as_string = seq_obj.tostring()
except AttributeError:
    # The new way, needs Biopython 1.45 or later.
    # Don't use this on Biopython 1.44 or older as truncates
    as_string = str(seq_obj)
__hash__(self)

Hash of the sequence as a string for comparison.

See Seq object comparison documentation (method __eq__ in particular) as this has changed in Biopython 1.65. Older versions would hash on object identity.

__eq__(self, other)

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

If you still need to support releases prior to Biopython 1.65, please just do explicit comparisons:

>>> from Bio.Seq import Seq
>>> 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
>>> seq1 == seq2
True
>>> seq1 == "ACGT"
True
__lt__(self, other)

Implement the less-than operand.

__le__(self, other)

Implement the less-than or equal operand.

__gt__(self, other)

Implement the greater-than operand.

__ge__(self, other)

Implement the greater-than or equal operand.

__len__(self)

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

__getitem__(self, index)

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

>>> my_seq = Seq('ACTCGACGTCG')
>>> my_seq[5]
'A'
__add__(self, other)

Add another sequence or string to this sequence.

>>> from Bio.Seq import Seq
>>> Seq("MELKI") + "LV"
Seq('MELKILV')
__radd__(self, other)

Add a sequence on the left.

>>> from Bio.Seq import Seq
>>> "LV" + Seq("MELKI")
Seq('LVMELKI')

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

__mul__(self, other)

Multiply Seq by integer.

>>> from Bio.Seq import Seq
>>> Seq('ATG') * 2
Seq('ATGATG')
__rmul__(self, other)

Multiply integer by Seq.

>>> from Bio.Seq import Seq
>>> 2 * Seq('ATG')
Seq('ATGATG')
__imul__(self, other)

Multiply Seq in-place.

>>> from Bio.Seq import Seq
>>> seq = Seq('ATG')
>>> seq *= 2
>>> seq
Seq('ATGATG')
tomutable(self)

Return the full sequence as a MutableSeq object.

>>> from Bio.Seq import Seq
>>> my_seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> my_seq
Seq('MKQHKAMIVALIVICITAVVAAL')
>>> my_seq.tomutable()
MutableSeq('MKQHKAMIVALIVICITAVVAAL')
count(self, sub, start=0, end=9223372036854775807)

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!

count_overlap(self, sub, start=0, end=9223372036854775807)

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.

__contains__(self, char)

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

e.g.

>>> from Bio.Seq import Seq
>>> my_dna = Seq("ATATGAAATTTGAAAA")
>>> "AAA" in my_dna
True
>>> Seq("AAA") in my_dna
True
find(self, sub, start=0, end=9223372036854775807)

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
rfind(self, sub, start=0, end=9223372036854775807)

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
index(self, sub, start=0, end=9223372036854775807)

Like find() but raise ValueError when the substring is not found.

>>> from Bio.Seq import Seq
>>> my_rna = Seq("GUCAUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAGUUG")
>>> my_rna.find("T")
-1
>>> my_rna.index("T")
Traceback (most recent call last):
           ...
ValueError: substring not found...
rindex(self, sub, start=0, end=9223372036854775807)

Like rfind() but raise ValueError when the substring is not found.

startswith(self, prefix, start=0, end=9223372036854775807)

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
endswith(self, suffix, start=0, end=9223372036854775807)

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
split(self, sep=None, maxsplit=- 1)

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')
>>> for pep in my_aa.split("*"):
...     pep
Seq('VMAIVMGR')
Seq('KGAR')
Seq('L')
>>> for pep in my_aa.split("*", 1):
...     pep
Seq('VMAIVMGR')
Seq('KGAR*L')

See also the rsplit method:

>>> for pep in my_aa.rsplit("*", 1):
...     pep
Seq('VMAIVMGR*KGAR')
Seq('L')
rsplit(self, sep=None, maxsplit=- 1)

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.

strip(self, chars=None)

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.

lstrip(self, chars=None)

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.

rstrip(self, chars=None)

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.Seq import Seq
>>> my_seq = Seq("CGGTACGCTTATGTCACGTAGAAAAAA")
>>> my_seq
Seq('CGGTACGCTTATGTCACGTAGAAAAAA')
>>> my_seq.rstrip("A")
Seq('CGGTACGCTTATGTCACGTAG')

See also the strip and lstrip methods.

upper(self)

Return an upper case copy of the sequence.

>>> from Bio.Seq import Seq
>>> my_seq = Seq("VHLTPeeK*")
>>> my_seq
Seq('VHLTPeeK*')
>>> my_seq.lower()
Seq('vhltpeek*')
>>> my_seq.upper()
Seq('VHLTPEEK*')
lower(self)

Return a lower case copy of the sequence.

>>> from Bio.Seq import Seq
>>> my_seq = Seq("CGGTACGCTTATGTCACGTAGAAAAAA")
>>> my_seq
Seq('CGGTACGCTTATGTCACGTAGAAAAAA')
>>> my_seq.lower()
Seq('cggtacgcttatgtcacgtagaaaaaa')

See also the upper method.

encode(self, encoding='utf-8', errors='strict')

Return an encoded version of the sequence as a bytes object.

The Seq object aims to match the interface of a Python string.

This is essentially to save you doing str(my_seq).encode() when you need a bytes string, for example for computing a hash:

>>> from Bio.Seq import Seq
>>> Seq("ACGT").encode("ascii")
b'ACGT'
complement(self)

Return the complement sequence by creating a new Seq object.

This method is intended for use with DNA sequences:

>>> from Bio.Seq import Seq
>>> my_dna = Seq("CCCCCGATAG")
>>> my_dna
Seq('CCCCCGATAG')
>>> my_dna.complement()
Seq('GGGGGCTATC')

You can of course used mixed case sequences,

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

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

Note that if the sequence contains neither T nor U, we assume it is DNA and map any A character to T:

>>> Seq("CGA").complement()
Seq('GCT')
>>> Seq("CGAT").complement()
Seq('GCTA')

If you actually have RNA, this currently works but we may deprecate this later. We recommend using the new complement_rna method instead:

>>> Seq("CGAU").complement()
Seq('GCUA')
>>> Seq("CGAU").complement_rna()
Seq('GCUA')

If the sequence contains both T and U, an exception is raised:

>>> Seq("CGAUT").complement()
Traceback (most recent call last):
   ...
ValueError: Mixed RNA/DNA found

Trying to complement a protein sequence gives a meaningless sequence:

>>> my_protein = Seq("MAIVMGR")
>>> my_protein.complement()
Seq('KTIBKCY')

Here “M” was interpreted as the IUPAC ambiguity code for “A” or “C”, with complement “K” for “T” or “G”. Likewise “A” has complement “T”. The letter “I” has no defined meaning under the IUPAC convention, and is unchanged.

reverse_complement(self)

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

This method is intended for use with DNA sequences:

>>> from Bio.Seq import Seq
>>> my_dna = Seq("CCCCCGATAGNR")
>>> my_dna
Seq('CCCCCGATAGNR')
>>> my_dna.reverse_complement()
Seq('YNCTATCGGGGG')

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
>>> my_dna = Seq("CCCCCgatA-G")
>>> my_dna
Seq('CCCCCgatA-G')
>>> my_dna.reverse_complement()
Seq('C-TatcGGGGG')

As discussed for the complement method, if the sequence contains neither T nor U, is is assumed to be DNA and will map any letter A to T.

If you are dealing with RNA you should use the new reverse_complement_rna method instead

>>> Seq("CGA").reverse_complement()  # defaults to DNA
Seq('TCG')
>>> Seq("CGA").reverse_complement_rna()
Seq('UCG')

If the sequence contains both T and U, an exception is raised:

>>> Seq("CGAUT").reverse_complement()
Traceback (most recent call last):
   ...
ValueError: Mixed RNA/DNA found

Trying to reverse complement a protein sequence will give a meaningless sequence:

>>> from Bio.Seq import Seq
>>> my_protein = Seq("MAIVMGR")
>>> my_protein.reverse_complement()
Seq('YCKBITK')

Here “M” was interpretted as the IUPAC ambiguity code for “A” or “C”, with complement “K” for “T” or “G” - and so on.

complement_rna(self)

Complement of an RNA sequence.

>>> Seq("CGA").complement()  # defaults to DNA
Seq('GCT')
>>> Seq("CGA").complement_rna()
Seq('GCU')

If the sequence contains both T and U, an exception is raised:

>>> Seq("CGAUT").complement_rna()
Traceback (most recent call last):
   ...
ValueError: Mixed RNA/DNA found

If the sequence contains T, an exception is raised:

>>> Seq("ACGT").complement_rna()
Traceback (most recent call last):
   ...
ValueError: DNA found, RNA expected
reverse_complement_rna(self)

Reverse complement of an RNA sequence.

>>> from Bio.Seq import Seq
>>> Seq("ACG").reverse_complement_rna()
Seq('CGU')
transcribe(self)

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

>>> from Bio.Seq import Seq
>>> coding_dna = Seq("ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG")
>>> coding_dna
Seq('ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG')
>>> coding_dna.transcribe()
Seq('AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG')

Trying to transcribe an RNA sequence should have no effect. If you have a nucleotide sequence which might be DNA or RNA (or even a mixture), calling the transcribe method will ensure any T becomes U.

Trying to transcribe a protein sequence will replace any T for Threonine with U for Selenocysteine, which has no biologically plausible rational. Older versions of Biopython would throw an exception.

>>> from Bio.Seq import Seq
>>> my_protein = Seq("MAIVMGRT")
>>> my_protein.transcribe()
Seq('MAIVMGRU')
back_transcribe(self)

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

>>> from Bio.Seq import Seq
>>> messenger_rna = Seq("AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG")
>>> messenger_rna
Seq('AUGGCCAUUGUAAUGGGCCGCUGAAAGGGUGCCCGAUAG')
>>> messenger_rna.back_transcribe()
Seq('ATGGCCATTGTAATGGGCCGCTGAAAGGGTGCCCGATAG')

Trying to back-transcribe DNA has no effect, If you have a nucleotide sequence which might be DNA or RNA (or even a mixture), calling the back-transcribe method will ensure any T becomes U.

Trying to back-transcribe a protein sequence will replace any U for Selenocysteine with T for Threonine, which is biologically meaningless. Older versions of Biopython would raise an exception here:

>>> from Bio.Seq import Seq
>>> my_protein = Seq("MAIVMGRU")
>>> my_protein.back_transcribe()
Seq('MAIVMGRT')
translate(self, table='Standard', stop_symbol='*', to_stop=False, cds=False, gap='-')

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

This method will translate DNA or RNA sequences. It should not be used on protein sequences as any result will be biologically meaningless.

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. Defaults to the minus sign.

e.g. Using the standard table:

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

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

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

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

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):
    ...
Bio.Data.CodonTable.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')
>>> coding_dna2.translate(to_stop=True)
Seq('LAIVMGR')

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

ungap(self, gap='-')

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

The gap character now defaults to the minus sign, and can only be specified via the method argument. This is no longer possible via the sequence’s alphabet (as was possible up to Biopython 1.77):

>>> from Bio.Seq import Seq
>>> my_dna = Seq("-ATA--TGAAAT-TTGAAAA")
>>> my_dna
Seq('-ATA--TGAAAT-TTGAAAA')
>>> my_dna.ungap("-")
Seq('ATATGAAATTTGAAAA')
join(self, other)

Return a merge of the sequences in other, spaced by the sequence from self.

Accepts either a Seq or string (and iterates over the letters), or an iterable containing Seq or string objects. These arguments will be concatenated with the calling sequence as the spacer:

>>> concatenated = Seq('NNNNN').join([Seq("AAA"), Seq("TTT"), Seq("PPP")])
>>> concatenated
Seq('AAANNNNNTTTNNNNNPPP')

Joining the letters of a single sequence:

>>> Seq('NNNNN').join(Seq("ACGT"))
Seq('ANNNNNCNNNNNGNNNNNT')
>>> Seq('NNNNN').join("ACGT")
Seq('ANNNNNCNNNNNGNNNNNT')
class Bio.Seq.UnknownSeq(length, alphabet=None, character='?')

Bases: Bio.Seq.Seq

Read-only sequence object of known length but unknown contents.

If you have an unknown sequence, you can represent this with a normal Seq object, for example:

>>> my_seq = Seq("N"*5)
>>> my_seq
Seq('NNNNN')
>>> len(my_seq)
5
>>> print(my_seq)
NNNNN

However, this is rather wasteful of memory (especially for large sequences), which is where this class is most useful:

>>> unk_five = UnknownSeq(5)
>>> unk_five
UnknownSeq(5, character='?')
>>> len(unk_five)
5
>>> print(unk_five)
?????

You can add unknown sequence together. Provided the characters are the same, you get another memory saving UnknownSeq:

>>> unk_four = UnknownSeq(4)
>>> unk_four
UnknownSeq(4, character='?')
>>> unk_four + unk_five
UnknownSeq(9, character='?')

If the characters are different, addition gives an ordinary Seq object:

>>> unk_nnnn = UnknownSeq(4, character="N")
>>> unk_nnnn
UnknownSeq(4, character='N')
>>> unk_nnnn + unk_four
Seq('NNNN????')

Combining with a real Seq gives a new Seq object:

>>> known_seq = Seq("ACGT")
>>> unk_four + known_seq
Seq('????ACGT')
>>> known_seq + unk_four
Seq('ACGT????')
__init__(self, length, alphabet=None, character='?')

Create a new UnknownSeq object.

Arguments:
  • length - Integer, required.

  • alphabet - no longer used, must be None.

  • character - single letter string, default “?”. Typically “N” for nucleotides, “X” for proteins, and “?” otherwise.

__len__(self)

Return the stated length of the unknown sequence.

__str__(self)

Return the unknown sequence as full string of the given length.

__repr__(self)

Return (truncated) representation of the sequence for debugging.

__add__(self, other)

Add another sequence or string to this sequence.

Adding two UnknownSeq objects returns another UnknownSeq object provided the character is the same.

>>> from Bio.Seq import UnknownSeq
>>> UnknownSeq(10, character='X') + UnknownSeq(5, character='X')
UnknownSeq(15, character='X')

If the characters differ, an UnknownSeq object cannot be used, so a Seq object is returned:

>>> from Bio.Seq import UnknownSeq
>>> UnknownSeq(10, character='X') + UnknownSeq(5, character="x")
Seq('XXXXXXXXXXxxxxx')

If adding a string to an UnknownSeq, a new Seq is returned:

>>> from Bio.Seq import UnknownSeq
>>> UnknownSeq(5, character='X') + "LV"
Seq('XXXXXLV')
__radd__(self, other)

Add a sequence on the left.

__mul__(self, other)

Multiply UnknownSeq by integer.

>>> from Bio.Seq import UnknownSeq
>>> UnknownSeq(3) * 2
UnknownSeq(6, character='?')
>>> UnknownSeq(3, character="N") * 2
UnknownSeq(6, character='N')
__rmul__(self, other)

Multiply integer by UnknownSeq.

>>> from Bio.Seq import UnknownSeq
>>> 2 * UnknownSeq(3)
UnknownSeq(6, character='?')
>>> 2 * UnknownSeq(3, character="N")
UnknownSeq(6, character='N')
__imul__(self, other)

Multiply UnknownSeq in-place.

>>> from Bio.Seq import UnknownSeq
>>> seq = UnknownSeq(3, character="N")
>>> seq *= 2
>>> seq
UnknownSeq(6, character='N')
__getitem__(self, index)

Get a subsequence from the UnknownSeq object.

>>> unk = UnknownSeq(8, character="N")
>>> print(unk[:])
NNNNNNNN
>>> print(unk[5:3])

>>> print(unk[1:-1])
NNNNNN
>>> print(unk[1:-1:2])
NNN
count(self, sub, start=0, end=9223372036854775807)

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

This behaves like the python string (and Seq object) 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

>>> "NNNN".count("N")
4
>>> Seq("NNNN").count("N")
4
>>> UnknownSeq(4, character="N").count("N")
4
>>> UnknownSeq(4, character="N").count("A")
0
>>> UnknownSeq(4, character="N").count("AA")
0

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

>>> UnknownSeq(4, character="N").count("NN")
2
>>> UnknownSeq(4, character="N").count("NNN")
1
count_overlap(self, sub, start=0, end=9223372036854775807)

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 UnknownSeq
>>> UnknownSeq(4, character="N").count_overlap("NN")
3
>>> UnknownSeq(4, character="N").count_overlap("NNN")
2

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

>>> UnknownSeq(4, character="N").count_overlap("N")
4
>>> UnknownSeq(4, character="N").count_overlap("N") == UnknownSeq(4, character="N").count("N")
True
>>> UnknownSeq(4, character="N").count_overlap("A")
0
>>> UnknownSeq(4, character="N").count_overlap("A") == UnknownSeq(4, character="N").count("A")
True
>>> UnknownSeq(4, character="N").count_overlap("AA")
0
>>> UnknownSeq(4, character="N").count_overlap("AA") == UnknownSeq(4, character="N").count("AA")
True
complement(self)

Return the complement of an unknown nucleotide equals itself.

>>> my_nuc = UnknownSeq(8)
>>> my_nuc
UnknownSeq(8, character='?')
>>> print(my_nuc)
????????
>>> my_nuc.complement()
UnknownSeq(8, character='?')
>>> print(my_nuc.complement())
????????
complement_rna(self)

Return the complement assuming it is RNA.

reverse_complement(self)

Return the reverse complement of an unknown sequence.

The reverse complement of an unknown nucleotide equals itself:

>>> from Bio.Seq import UnknownSeq
>>> example = UnknownSeq(6, character="N")
>>> print(example)
NNNNNN
>>> print(example.reverse_complement())
NNNNNN
reverse_complement_rna(self)

Return the reverse complement assuming it is RNA.

transcribe(self)

Return an unknown RNA sequence from an unknown DNA sequence.

>>> my_dna = UnknownSeq(10, character="N")
>>> my_dna
UnknownSeq(10, character='N')
>>> print(my_dna)
NNNNNNNNNN
>>> my_rna = my_dna.transcribe()
>>> my_rna
UnknownSeq(10, character='N')
>>> print(my_rna)
NNNNNNNNNN
back_transcribe(self)

Return an unknown DNA sequence from an unknown RNA sequence.

>>> my_rna = UnknownSeq(20, character="N")
>>> my_rna
UnknownSeq(20, character='N')
>>> print(my_rna)
NNNNNNNNNNNNNNNNNNNN
>>> my_dna = my_rna.back_transcribe()
>>> my_dna
UnknownSeq(20, character='N')
>>> print(my_dna)
NNNNNNNNNNNNNNNNNNNN
upper(self)

Return an upper case copy of the sequence.

>>> from Bio.Seq import UnknownSeq
>>> my_seq = UnknownSeq(20, character="n")
>>> my_seq
UnknownSeq(20, character='n')
>>> print(my_seq)
nnnnnnnnnnnnnnnnnnnn
>>> my_seq.upper()
UnknownSeq(20, character='N')
>>> print(my_seq.upper())
NNNNNNNNNNNNNNNNNNNN

See also the lower method.

lower(self)

Return a lower case copy of the sequence.

>>> from Bio.Seq import UnknownSeq
>>> my_seq = UnknownSeq(20, character="X")
>>> my_seq
UnknownSeq(20, character='X')
>>> print(my_seq)
XXXXXXXXXXXXXXXXXXXX
>>> my_seq.lower()
UnknownSeq(20, character='x')
>>> print(my_seq.lower())
xxxxxxxxxxxxxxxxxxxx

See also the upper method.

translate(self, table='Standard', stop_symbol='*', to_stop=False, cds=False, gap='-')

Translate an unknown nucleotide sequence into an unknown protein.

e.g.

>>> my_seq = UnknownSeq(9, character="N")
>>> print(my_seq)
NNNNNNNNN
>>> my_protein = my_seq.translate()
>>> my_protein
UnknownSeq(3, character='X')
>>> print(my_protein)
XXX

In comparison, using a normal Seq object:

>>> my_seq = Seq("NNNNNNNNN")
>>> print(my_seq)
NNNNNNNNN
>>> my_protein = my_seq.translate()
>>> my_protein
Seq('XXX')
>>> print(my_protein)
XXX
ungap(self, gap='-')

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

The gap character now defaults to the minus sign, and can only be specified via the method argument. This is no longer possible via the sequence’s alphabet (as was possible up to Biopython 1.77):

>>> from Bio.Seq import UnknownSeq
>>> my_dna = UnknownSeq(20, character='N')
>>> my_dna
UnknownSeq(20, character='N')
>>> my_dna.ungap()  # using default
UnknownSeq(20, character='N')
>>> my_dna.ungap("-")
UnknownSeq(20, character='N')

If the UnknownSeq is using the gap character, then an empty Seq is returned:

>>> my_gap = UnknownSeq(20, character="-")
>>> my_gap
UnknownSeq(20, character='-')
>>> my_gap.ungap()  # using default
Seq('')
>>> my_gap.ungap("-")
Seq('')
join(self, other)

Return a merge of the sequences in other, spaced by the sequence from self.

Accepts either a Seq or string (and iterates over the letters), or an iterable containing Seq or string objects. These arguments will be concatenated with the calling sequence as the spacer:

>>> concatenated = UnknownSeq(5).join([Seq("AAA"), Seq("TTT"), Seq("PPP")])
>>> concatenated
Seq('AAA?????TTT?????PPP')

If all the inputs are also UnknownSeq using the same character, then it returns a new UnknownSeq:

>>> UnknownSeq(5).join([UnknownSeq(3), UnknownSeq(3), UnknownSeq(3)])
UnknownSeq(19, character='?')

Examples taking a single sequence and joining the letters:

>>> UnknownSeq(3).join("ACGT")
Seq('A???C???G???T')
>>> UnknownSeq(3).join(UnknownSeq(4))
UnknownSeq(13, character='?')

Will only return an UnknownSeq object if all of the objects to be joined are also UnknownSeqs with the same character as the spacer, similar to how the addition of an UnknownSeq and another UnknownSeq would work.

class Bio.Seq.MutableSeq(data)

Bases: object

An editable sequence object.

Unlike normal python strings and our basic sequence object (the Seq class) which are immutable, the MutableSeq lets you edit the sequence in place. However, this means you cannot use a MutableSeq object as a dictionary key.

>>> from Bio.Seq import MutableSeq
>>> my_seq = MutableSeq("ACTCGTCGTCG")
>>> my_seq
MutableSeq('ACTCGTCGTCG')
>>> my_seq[5]
'T'
>>> my_seq[5] = "A"
>>> my_seq
MutableSeq('ACTCGACGTCG')
>>> my_seq[5]
'A'
>>> my_seq[5:8] = "NNN"
>>> my_seq
MutableSeq('ACTCGNNNTCG')
>>> len(my_seq)
11

Note that the MutableSeq object does not support as many string-like or biological methods as the Seq object.

__init__(self, data)

Initialize the class.

__repr__(self)

Return (truncated) representation of the sequence for debugging.

__str__(self)

Return the full sequence as a python string.

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 needs to be backwards compatible with old Biopython, you should continue to use my_seq.tostring() rather than str(my_seq).

__eq__(self, other)

Compare the sequence to another sequence or a string.

Historically comparing DNA to RNA, or Nucleotide to Protein would raise an exception. This was later downgraded to a warning, but since Biopython 1.78 the alphabet is ignored for comparisons.

If you need to support older Biopython versions, please just do explicit comparisons:

>>> seq1 = MutableSeq("ACGT")
>>> seq2 = MutableSeq("ACGT")
>>> id(seq1) == id(seq2)
False
>>> str(seq1) == str(seq2)
True

Biopython now does:

>>> seq1 == seq2
True
>>> seq1 == Seq("ACGT")
True
>>> seq1 == "ACGT"
True
__lt__(self, other)

Implement the less-than operand.

__le__(self, other)

Implement the less-than or equal operand.

__gt__(self, other)

Implement the greater-than operand.

__ge__(self, other)

Implement the greater-than or equal operand.

__len__(self)

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

__getitem__(self, index)

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

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq[5]
'A'
__setitem__(self, index, value)

Set a subsequence of single letter via value parameter.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq[0] = 'T'
>>> my_seq
MutableSeq('TCTCGACGTCG')
__delitem__(self, index)

Delete a subsequence of single letter.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> del my_seq[0]
>>> my_seq
MutableSeq('CTCGACGTCG')
__add__(self, other)

Add another sequence or string to this sequence.

Returns a new MutableSeq object.

__radd__(self, other)

Add a sequence on the left.

>>> from Bio.Seq import MutableSeq
>>> "LV" + MutableSeq("MELKI")
MutableSeq('LVMELKI')
__mul__(self, other)

Multiply MutableSeq by integer.

Note this is not in-place and returns a new object, matching native Python list multiplication.

>>> from Bio.Seq import MutableSeq
>>> MutableSeq('ATG') * 2
MutableSeq('ATGATG')
__rmul__(self, other)

Multiply integer by MutableSeq.

Note this is not in-place and returns a new object, matching native Python list multiplication.

>>> from Bio.Seq import MutableSeq
>>> 2 * MutableSeq('ATG')
MutableSeq('ATGATG')
__imul__(self, other)

Multiply MutableSeq in-place.

>>> from Bio.Seq import MutableSeq
>>> seq = MutableSeq('ATG')
>>> seq *= 2
>>> seq
MutableSeq('ATGATG')
append(self, c)

Add a subsequence to the mutable sequence object.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq.append('A')
>>> my_seq
MutableSeq('ACTCGACGTCGA')

No return value.

insert(self, i, c)

Add a subsequence to the mutable sequence object at a given index.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq.insert(0,'A')
>>> my_seq
MutableSeq('AACTCGACGTCG')
>>> my_seq.insert(8,'G')
>>> my_seq
MutableSeq('AACTCGACGGTCG')

No return value.

pop(self, i=- 1)

Remove a subsequence of a single letter at given index.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq.pop()
'G'
>>> my_seq
MutableSeq('ACTCGACGTC')
>>> my_seq.pop()
'C'
>>> my_seq
MutableSeq('ACTCGACGT')

Returns the last character of the sequence.

remove(self, item)

Remove a subsequence of a single letter from mutable sequence.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq.remove('C')
>>> my_seq
MutableSeq('ATCGACGTCG')
>>> my_seq.remove('A')
>>> my_seq
MutableSeq('TCGACGTCG')

No return value.

count(self, sub, start=0, end=9223372036854775807)

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 MutableSeq
>>> my_mseq = MutableSeq("AAAATGA")
>>> print(my_mseq.count("A"))
5
>>> print(my_mseq.count("ATG"))
1
>>> print(my_mseq.count(Seq("AT")))
1
>>> print(my_mseq.count("AT", 2, -1))
1

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

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

An overlapping search would give the answer as three!

count_overlap(self, sub, start=0, end=9223372036854775807)

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 MutableSeq
>>> print(MutableSeq("AAAA").count_overlap("AA"))
3
>>> print(MutableSeq("ATATATATA").count_overlap("ATA"))
4
>>> print(MutableSeq("ATATATATA").count_overlap("ATA", 3, -1))
1

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

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

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

index(self, item)

Return first occurrence position of a single entry (i.e. letter).

>>> my_seq = MutableSeq("ACTCGACGTCG")
>>> my_seq.index("A")
0
>>> my_seq.index("T")
2
>>> my_seq.index(Seq("T"))
2

Note unlike a Biopython Seq object, or Python string, multi-letter subsequences are not supported. Instead this acts like an array or a list of the entries. There is therefore no .rindex() method.

reverse(self)

Modify the mutable sequence to reverse itself.

No return value.

complement(self)

Modify the mutable sequence to take on its complement.

No return value.

If the sequence contains neither T nor U, DNA is assumed and any A will be mapped to T.

If the sequence contains both T and U, an exception is raised.

reverse_complement(self)

Modify the mutable sequence to take on its reverse complement.

No return value.

extend(self, other)

Add a sequence to the original mutable sequence object.

>>> my_seq = MutableSeq('ACTCGACGTCG')
>>> my_seq.extend('A')
>>> my_seq
MutableSeq('ACTCGACGTCGA')
>>> my_seq.extend('TTT')
>>> my_seq
MutableSeq('ACTCGACGTCGATTT')

No return value.

toseq(self)

Return the full sequence as a new immutable Seq object.

>>> from Bio.Seq import MutableSeq
>>> my_mseq = MutableSeq("MKQHKAMIVALIVICITAVVAAL")
>>> my_mseq
MutableSeq('MKQHKAMIVALIVICITAVVAAL')
>>> my_mseq.toseq()
Seq('MKQHKAMIVALIVICITAVVAAL')
join(self, other)

Return a merge of the sequences in other, spaced by the sequence from self.

Accepts all Seq objects and Strings as objects to be concatenated with the spacer

>>> concatenated = MutableSeq('NNNNN').join([Seq("AAA"), Seq("TTT"), Seq("PPP")])
>>> concatenated
Seq('AAANNNNNTTTNNNNNPPP')

Throws error if other is not an iterable and if objects inside of the iterable are not Seq or String objects

__hash__ = None
Bio.Seq.transcribe(dna)

Transcribe a DNA sequence into RNA.

If given a string, returns a new string object.

Given a Seq or MutableSeq, returns a new Seq object.

e.g.

>>> transcribe("ACTGN")
'ACUGN'
Bio.Seq.back_transcribe(rna)

Return the RNA sequence back-transcribed into DNA.

If given a string, returns a new string object.

Given a Seq or MutableSeq, returns a new Seq object.

e.g.

>>> back_transcribe("ACUGN")
'ACTGN'
Bio.Seq.translate(sequence, table='Standard', stop_symbol='*', to_stop=False, cds=False, gap=None)

Translate a nucleotide sequence into amino acids.

If given a string, returns a new string object. Given a Seq or MutableSeq, returns a Seq object.

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). Defaults to the “Standard” table.

  • stop_symbol - Single character string, what to use for any terminators, 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. Defaults to None.

A simple string example using the default (standard) genetic code:

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

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

>>> translate(coding_dna, table=2)
'VAIVMGRWKGAR*'
>>> translate(coding_dna, table=2, to_stop=True)
'VAIVMGRWKGAR'

In fact this example uses an alternative start codon valid under NCBI table 2, GTG, which means this example is a complete valid CDS which when translated should really start with methionine (not valine):

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

Note that if the sequence has no in-frame stop codon, then the to_stop argument has no effect:

>>> coding_dna2 = "GTGGCCATTGTAATGGGCCGC"
>>> translate(coding_dna2)
'VAIVMGR'
>>> translate(coding_dna2, to_stop=True)
'VAIVMGR'

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.

It will however translate either DNA or RNA.

NOTE - Since version 1.71 Biopython contains codon tables with ‘ambiguous stop codons’. These are stop codons with unambiguous sequence but which have a context dependent coding as STOP or as amino acid. With these tables ‘to_stop’ must be False (otherwise a ValueError is raised). The dual coding codons will always be translated as amino acid, except for ‘cds=True’, where the last codon will be translated as STOP.

>>> coding_dna3 = "ATGGCACGGAAGTGA"
>>> translate(coding_dna3)
'MARK*'
>>> translate(coding_dna3, table=27)  # Table 27: TGA -> STOP or W
'MARKW'

It will however raise a BiopythonWarning (not shown).

>>> translate(coding_dna3, table=27, cds=True)
'MARK'
>>> translate(coding_dna3, table=27, to_stop=True)
Traceback (most recent call last):
   ...
ValueError: You cannot use 'to_stop=True' with this table ...
Bio.Seq.reverse_complement(sequence)

Return the reverse complement sequence of a nucleotide string.

If given a string, returns a new string object. Given a Seq or a MutableSeq, returns a new Seq object.

Supports unambiguous and ambiguous nucleotide sequences.

e.g.

>>> reverse_complement("ACTG-NH")
'DN-CAGT'

If neither T nor U is present, DNA is assumed and A is mapped to T:

>>> reverse_complement("A")
'T'
Bio.Seq.complement(sequence)

Return the complement sequence of a DNA string.

If given a string, returns a new string object.

Given a Seq or a MutableSeq, returns a new Seq object.

Supports unambiguous and ambiguous nucleotide sequences.

e.g.

>>> complement("ACTG-NH")
'TGAC-ND'

If neither T nor U is present, DNA is assumed and A is mapped to T:

>>> complement("A")
'T'

However, this may not be supported in future. Please use the complement_rna function if you have RNA.

Bio.Seq.complement_rna(sequence)

Return the complement sequence of an RNA string.

>>> complement("ACG")  # assumed DNA
'TGC'
>>> complement_rna("ACG")
'UGC'

If the sequence contains a T, and error is raised.