Bio.SeqFeature module

Represent a Sequence Feature holding info about a part of a sequence.

This is heavily modeled after the Biocorba SeqFeature objects, and may be pretty biased towards GenBank stuff since I’m writing it for the GenBank parser output…

What’s here:

Base class to hold a Feature

Classes:
  • SeqFeature

Hold information about a Reference

This is an attempt to create a General class to hold Reference type information.

Classes:
  • Reference

Specify locations of a feature on a Sequence

This aims to handle, in Ewan Birney’s words, ‘the dreaded fuzziness issue’. This has the advantages of allowing us to handle fuzzy stuff in case anyone needs it, and also be compatible with BioPerl etc and BioSQL.

Classes:
  • Location - abstract base class of SimpleLocation and CompoundLocation.

  • SimpleLocation - Specify the start and end location of a feature.

  • CompoundLocation - Collection of SimpleLocation objects (for joins etc).

  • Position - abstract base class of ExactPosition, WithinPosition, BetweenPosition, AfterPosition, OneOfPosition, UncertainPosition, and UnknownPosition.

  • ExactPosition - Specify the position as being exact.

  • WithinPosition - Specify a position occurring within some range.

  • BetweenPosition - Specify a position occurring between a range (OBSOLETE?).

  • BeforePosition - Specify the position as being found before some base.

  • AfterPosition - Specify the position as being found after some base.

  • OneOfPosition - Specify a position consisting of multiple alternative positions.

  • UncertainPosition - Specify a specific position which is uncertain.

  • UnknownPosition - Represents missing information like ‘?’ in UniProt.

Exceptions:
  • LocationParserError - Exception indicating a failure to parse a location string.

exception Bio.SeqFeature.LocationParserError

Bases: ValueError

Could not parse a feature location string.

class Bio.SeqFeature.SeqFeature(location=None, type='', location_operator='', strand=None, id='<unknown id>', qualifiers=None, sub_features=None, ref=None, ref_db=None)

Bases: object

Represent a Sequence Feature on an object.

Attributes:
  • location - the location of the feature on the sequence (SimpleLocation)

  • type - the specified type of the feature (ie. CDS, exon, repeat…)

  • location_operator - a string specifying how this SeqFeature may be related to others. For example, in the example GenBank feature shown below, the location_operator would be “join”. This is a proxy for feature.location.operator and only applies to compound locations.

  • strand - A value specifying on which strand (of a DNA sequence, for instance) the feature deals with. 1 indicates the plus strand, -1 indicates the minus strand, 0 indicates stranded but unknown (? in GFF3), while the default of None indicates that strand doesn’t apply (dot in GFF3, e.g. features on proteins). Note this is a shortcut for accessing the strand property of the feature’s location.

  • id - A string identifier for the feature.

  • ref - A reference to another sequence. This could be an accession number for some different sequence. Note this is a shortcut for the reference property of the feature’s location.

  • ref_db - A different database for the reference accession number. Note this is a shortcut for the reference property of the location

  • qualifiers - A dictionary of qualifiers on the feature. These are analogous to the qualifiers from a GenBank feature table. The keys of the dictionary are qualifier names, the values are the qualifier values.

__init__(location=None, type='', location_operator='', strand=None, id='<unknown id>', qualifiers=None, sub_features=None, ref=None, ref_db=None)

Initialize a SeqFeature on a sequence.

location can either be a SimpleLocation (with strand argument also given if required), or None.

e.g. With no strand, on the forward strand, and on the reverse strand:

>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> f1 = SeqFeature(SimpleLocation(5, 10), type="domain")
>>> f1.strand == f1.location.strand == None
True
>>> f2 = SeqFeature(SimpleLocation(7, 110, strand=1), type="CDS")
>>> f2.strand == f2.location.strand == +1
True
>>> f3 = SeqFeature(SimpleLocation(9, 108, strand=-1), type="CDS")
>>> f3.strand == f3.location.strand == -1
True

For exact start/end positions, an integer can be used (as shown above) as shorthand for the ExactPosition object. For non-exact locations, the SimpleLocation must be specified via the appropriate position objects.

Note that the strand, ref and ref_db arguments to the SeqFeature are now deprecated and will later be removed. Set them via the location object instead.

Note that location_operator and sub_features arguments can no longer be used, instead do this via the CompoundLocation object.

property strand

Feature’s strand

This is a shortcut for feature.location.strand

property ref

Feature location reference (e.g. accession).

This is a shortcut for feature.location.ref

property ref_db

Feature location reference’s database.

This is a shortcut for feature.location.ref_db

property location_operator

Location operator for compound locations (e.g. join).

__eq__(other)

Check if two SeqFeature objects should be considered equal.

__repr__()

Represent the feature as a string for debugging.

__str__()

Return the full feature as a python string.

extract(parent_sequence, references=None)

Extract the feature’s sequence from supplied parent sequence.

The parent_sequence can be a Seq like object or a string, and will generally return an object of the same type. The exception to this is a MutableSeq as the parent sequence will return a Seq object.

This should cope with complex locations including complements, joins and fuzzy positions. Even mixed strand features should work! This also covers features on protein sequences (e.g. domains), although here reverse strand features are not permitted. If the location refers to other records, they must be supplied in the optional dictionary references.

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> f = SeqFeature(SimpleLocation(8, 15), type="domain")
>>> f.extract(seq)
Seq('VALIVIC')

If the SimpleLocation is None, e.g. when parsing invalid locus locations in the GenBank parser, extract() will raise a ValueError.

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SeqFeature
>>> seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> f = SeqFeature(None, type="domain")
>>> f.extract(seq)
Traceback (most recent call last):
   ...
ValueError: The feature's .location is None. Check the sequence file for a valid location.

Note - currently only compound features of type “join” are supported.

translate(parent_sequence, table='Standard', start_offset=None, stop_symbol='*', to_stop=False, cds=None, gap=None)

Get a translation of the feature’s sequence.

This method is intended for CDS or other features that code proteins and is a shortcut that will both extract the feature and translate it, taking into account the codon_start and transl_table qualifiers, if they are present. If they are not present the value of the arguments “table” and “start_offset” are used.

The “cds” parameter is set to “True” if the feature is of type “CDS” but can be overridden by giving an explicit argument.

The arguments stop_symbol, to_stop and gap have the same meaning as Seq.translate, refer to that documentation for further information.

Arguments:
  • parent_sequence - A DNA or RNA sequence.

  • table - Which codon table to use if there is no transl_table qualifier for this feature. 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.

  • start_offset - offset at which the first complete codon of a coding feature can be found, relative to the first base of that feature. Has a valid value of 0, 1 or 2. NOTE: this uses python’s 0-based numbering whereas the codon_start qualifier in files from NCBI use 1-based numbering. Will override a codon_start qualifier

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> seq = Seq("GGTTACACTTACCGATAATGTCTCTGATGA")
>>> f = SeqFeature(SimpleLocation(0, 30), type="CDS")
>>> f.qualifiers['transl_table'] = [11]

Note that features of type CDS are subject to the usual checks at translation. But you can override this behavior by giving explicit arguments:

>>> f.translate(seq, cds=False)
Seq('GYTYR*CL**')

Now use the start_offset argument to change the frame. Note this uses python 0-based numbering.

>>> f.translate(seq, start_offset=1, cds=False)
Seq('VTLTDNVSD')

Alternatively use the codon_start qualifier to do the same thing. Note: this uses 1-based numbering, which is found in files from NCBI.

>>> f.qualifiers['codon_start'] = [2]
>>> f.translate(seq, cds=False)
Seq('VTLTDNVSD')
__bool__()

Boolean value of an instance of this class (True).

This behavior is for backwards compatibility, since until the __len__ method was added, a SeqFeature always evaluated as True.

Note that in comparison, Seq objects, strings, lists, etc, will all evaluate to False if they have length zero.

WARNING: The SeqFeature may in future evaluate to False when its length is zero (in order to better match normal python behavior)!

__len__()

Return the length of the region where the feature is located.

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> f = SeqFeature(SimpleLocation(8, 15), type="domain")
>>> len(f)
7
>>> f.extract(seq)
Seq('VALIVIC')
>>> len(f.extract(seq))
7

This is a proxy for taking the length of the feature’s location:

>>> len(f.location)
7

For simple features this is the same as the region spanned (end position minus start position using Pythonic counting). However, for a compound location (e.g. a CDS as the join of several exons) the gaps are not counted (e.g. introns). This ensures that len(f) matches len(f.extract(parent_seq)), and also makes sure things work properly with features wrapping the origin etc.

__iter__()

Iterate over the parent positions within the feature.

The iteration order is strand aware, and can be thought of as moving along the feature using the parent sequence coordinates:

>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> f = SeqFeature(SimpleLocation(5, 10, strand=-1), type="domain")
>>> len(f)
5
>>> for i in f: print(i)
9
8
7
6
5
>>> list(f)
[9, 8, 7, 6, 5]

This is a proxy for iterating over the location,

>>> list(f.location)
[9, 8, 7, 6, 5]
__contains__(value)

Check if an integer position is within the feature.

>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> f = SeqFeature(SimpleLocation(5, 10, strand=-1), type="domain")
>>> len(f)
5
>>> [i for i in range(15) if i in f]
[5, 6, 7, 8, 9]

For example, to see which features include a SNP position, you could use this:

>>> from Bio import SeqIO
>>> record = SeqIO.read("GenBank/NC_000932.gb", "gb")
>>> for f in record.features:
...     if 1750 in f:
...         print("%s %s" % (f.type, f.location))
source [0:154478](+)
gene [1716:4347](-)
tRNA join{[4310:4347](-), [1716:1751](-)}

Note that for a feature defined as a join of several subfeatures (e.g. the union of several exons) the gaps are not checked (e.g. introns). In this example, the tRNA location is defined in the GenBank file as complement(join(1717..1751,4311..4347)), so that position 1760 falls in the gap:

>>> for f in record.features:
...     if 1760 in f:
...         print("%s %s" % (f.type, f.location))
source [0:154478](+)
gene [1716:4347](-)

Note that additional care may be required with fuzzy locations, for example just before a BeforePosition:

>>> from Bio.SeqFeature import SeqFeature, SimpleLocation
>>> from Bio.SeqFeature import BeforePosition
>>> f = SeqFeature(SimpleLocation(BeforePosition(3), 8), type="domain")
>>> len(f)
5
>>> [i for i in range(10) if i in f]
[3, 4, 5, 6, 7]

Note that is is a proxy for testing membership on the location.

>>> [i for i in range(10) if i in f.location]
[3, 4, 5, 6, 7]
__hash__ = None
class Bio.SeqFeature.Reference

Bases: object

Represent a Generic Reference object.

Attributes:
  • location - A list of Location objects specifying regions of the sequence that the references correspond to. If no locations are specified, the entire sequence is assumed.

  • authors - A big old string, or a list split by author, of authors for the reference.

  • title - The title of the reference.

  • journal - Journal the reference was published in.

  • medline_id - A medline reference for the article.

  • pubmed_id - A pubmed reference for the article.

  • comment - A place to stick any comments about the reference.

__init__()

Initialize the class.

__str__()

Return the full Reference object as a python string.

__repr__()

Represent the Reference object as a string for debugging.

__eq__(other)

Check if two Reference objects should be considered equal.

Note prior to Biopython 1.70 the location was not compared, as until then __eq__ for the SimpleLocation class was not defined.

__hash__ = None
class Bio.SeqFeature.Location

Bases: abc.ABC

Abstract base class representing a location.

abstract __repr__()

Represent the Location object as a string for debugging.

fromstring(length=None, circular=False, stranded=True)

Create a Location object from a string.

This should accept any valid location string in the INSDC Feature Table format (https://www.insdc.org/submitting-standards/feature-table/) as used in GenBank, DDBJ and EMBL files.

Simple examples:

>>> Location.fromstring("123..456", 1000)
SimpleLocation(ExactPosition(122), ExactPosition(456), strand=1)
>>> Location.fromstring("complement(<123..>456)", 1000)
SimpleLocation(BeforePosition(122), AfterPosition(456), strand=-1)

A more complex location using within positions,

>>> Location.fromstring("(9.10)..(20.25)", 1000)
SimpleLocation(WithinPosition(8, left=8, right=9), WithinPosition(25, left=20, right=25), strand=1)

Notice how that will act as though it has overall start 8 and end 25.

Zero length between feature,

>>> Location.fromstring("123^124", 1000)
SimpleLocation(ExactPosition(123), ExactPosition(123), strand=1)

The expected sequence length is needed for a special case, a between position at the start/end of a circular genome:

>>> Location.fromstring("1000^1", 1000)
SimpleLocation(ExactPosition(1000), ExactPosition(1000), strand=1)

Apart from this special case, between positions P^Q must have P+1==Q,

>>> Location.fromstring("123^456", 1000)
Traceback (most recent call last):
   ...
Bio.SeqFeature.LocationParserError: invalid feature location '123^456'

You can optionally provide a reference name:

>>> Location.fromstring("AL391218.9:105173..108462", 2000000)
SimpleLocation(ExactPosition(105172), ExactPosition(108462), strand=1, ref='AL391218.9')
>>> Location.fromstring("<2644..159", 2868, "circular")
CompoundLocation([SimpleLocation(BeforePosition(2643), ExactPosition(2868), strand=1), SimpleLocation(ExactPosition(0), ExactPosition(159), strand=1)], 'join')
__abstractmethods__ = frozenset({'__repr__'})
class Bio.SeqFeature.SimpleLocation(start, end, strand=None, ref=None, ref_db=None)

Bases: Bio.SeqFeature.Location

Specify the location of a feature along a sequence.

The SimpleLocation is used for simple continuous features, which can be described as running from a start position to and end position (optionally with a strand and reference information). More complex locations made up from several non-continuous parts (e.g. a coding sequence made up of several exons) are described using a SeqFeature with a CompoundLocation.

Note that the start and end location numbering follow Python’s scheme, thus a GenBank entry of 123..150 (one based counting) becomes a location of [122:150] (zero based counting).

>>> from Bio.SeqFeature import SimpleLocation
>>> f = SimpleLocation(122, 150)
>>> print(f)
[122:150]
>>> print(f.start)
122
>>> print(f.end)
150
>>> print(f.strand)
None

Note the strand defaults to None. If you are working with nucleotide sequences you’d want to be explicit if it is the forward strand:

>>> from Bio.SeqFeature import SimpleLocation
>>> f = SimpleLocation(122, 150, strand=+1)
>>> print(f)
[122:150](+)
>>> print(f.strand)
1

Note that for a parent sequence of length n, the SimpleLocation start and end must satisfy the inequality 0 <= start <= end <= n. This means even for features on the reverse strand of a nucleotide sequence, we expect the ‘start’ coordinate to be less than the ‘end’.

>>> from Bio.SeqFeature import SimpleLocation
>>> r = SimpleLocation(122, 150, strand=-1)
>>> print(r)
[122:150](-)
>>> print(r.start)
122
>>> print(r.end)
150
>>> print(r.strand)
-1

i.e. Rather than thinking of the ‘start’ and ‘end’ biologically in a strand aware manner, think of them as the ‘left most’ or ‘minimum’ boundary, and the ‘right most’ or ‘maximum’ boundary of the region being described. This is particularly important with compound locations describing non-continuous regions.

In the example above we have used standard exact positions, but there are also specialised position objects used to represent fuzzy positions as well, for example a GenBank location like complement(<123..150) would use a BeforePosition object for the start.

__init__(start, end, strand=None, ref=None, ref_db=None)

Initialize the class.

start and end arguments specify the values where the feature begins and ends. These can either by any of the *Position objects that inherit from Position, or can just be integers specifying the position. In the case of integers, the values are assumed to be exact and are converted in ExactPosition arguments. This is meant to make it easy to deal with non-fuzzy ends.

i.e. Short form:

>>> from Bio.SeqFeature import SimpleLocation
>>> loc = SimpleLocation(5, 10, strand=-1)
>>> print(loc)
[5:10](-)

Explicit form:

>>> from Bio.SeqFeature import SimpleLocation, ExactPosition
>>> loc = SimpleLocation(ExactPosition(5), ExactPosition(10), strand=-1)
>>> print(loc)
[5:10](-)

Other fuzzy positions are used similarly,

>>> from Bio.SeqFeature import SimpleLocation
>>> from Bio.SeqFeature import BeforePosition, AfterPosition
>>> loc2 = SimpleLocation(BeforePosition(5), AfterPosition(10), strand=-1)
>>> print(loc2)
[<5:>10](-)

For nucleotide features you will also want to specify the strand, use 1 for the forward (plus) strand, -1 for the reverse (negative) strand, 0 for stranded but strand unknown (? in GFF3), or None for when the strand does not apply (dot in GFF3), e.g. features on proteins.

>>> loc = SimpleLocation(5, 10, strand=+1)
>>> print(loc)
[5:10](+)
>>> print(loc.strand)
1

Normally feature locations are given relative to the parent sequence you are working with, but an explicit accession can be given with the optional ref and db_ref strings:

>>> loc = SimpleLocation(105172, 108462, ref="AL391218.9", strand=1)
>>> print(loc)
AL391218.9[105172:108462](+)
>>> print(loc.ref)
AL391218.9
static fromstring(text, length=None, circular=False)

Create a SimpleLocation object from a string.

property strand

Strand of the location (+1, -1, 0 or None).

__str__()

Return a representation of the SimpleLocation object (with python counting).

For the simple case this uses the python splicing syntax, [122:150] (zero based counting) which GenBank would call 123..150 (one based counting).

__repr__()

Represent the SimpleLocation object as a string for debugging.

__add__(other)

Combine location with another SimpleLocation object, or shift it.

You can add two feature locations to make a join CompoundLocation:

>>> from Bio.SeqFeature import SimpleLocation
>>> f1 = SimpleLocation(5, 10)
>>> f2 = SimpleLocation(20, 30)
>>> combined = f1 + f2
>>> print(combined)
join{[5:10], [20:30]}

This is thus equivalent to:

>>> from Bio.SeqFeature import CompoundLocation
>>> join = CompoundLocation([f1, f2])
>>> print(join)
join{[5:10], [20:30]}

You can also use sum(…) in this way:

>>> join = sum([f1, f2])
>>> print(join)
join{[5:10], [20:30]}

Furthermore, you can combine a SimpleLocation with a CompoundLocation in this way.

Separately, adding an integer will give a new SimpleLocation with its start and end offset by that amount. For example:

>>> print(f1)
[5:10]
>>> print(f1 + 100)
[105:110]
>>> print(200 + f1)
[205:210]

This can be useful when editing annotation.

__radd__(other)

Return a SimpleLocation object by shifting the location by an integer amount.

__sub__(other)

Subtracting an integer will shift the start and end by that amount.

>>> from Bio.SeqFeature import SimpleLocation
>>> f1 = SimpleLocation(105, 150)
>>> print(f1)
[105:150]
>>> print(f1 - 100)
[5:50]

This can be useful when editing annotation. You can also add an integer to a feature location (which shifts in the opposite direction).

__nonzero__()

Return True regardless of the length of the feature.

This behavior is for backwards compatibility, since until the __len__ method was added, a SimpleLocation always evaluated as True.

Note that in comparison, Seq objects, strings, lists, etc, will all evaluate to False if they have length zero.

WARNING: The SimpleLocation may in future evaluate to False when its length is zero (in order to better match normal python behavior)!

__len__()

Return the length of the region described by the SimpleLocation object.

Note that extra care may be needed for fuzzy locations, e.g.

>>> from Bio.SeqFeature import SimpleLocation
>>> from Bio.SeqFeature import BeforePosition, AfterPosition
>>> loc = SimpleLocation(BeforePosition(5), AfterPosition(10))
>>> len(loc)
5
__contains__(value)

Check if an integer position is within the SimpleLocation object.

Note that extra care may be needed for fuzzy locations, e.g.

>>> from Bio.SeqFeature import SimpleLocation
>>> from Bio.SeqFeature import BeforePosition, AfterPosition
>>> loc = SimpleLocation(BeforePosition(5), AfterPosition(10))
>>> len(loc)
5
>>> [i for i in range(15) if i in loc]
[5, 6, 7, 8, 9]
__iter__()

Iterate over the parent positions within the SimpleLocation object.

>>> from Bio.SeqFeature import SimpleLocation
>>> from Bio.SeqFeature import BeforePosition, AfterPosition
>>> loc = SimpleLocation(BeforePosition(5), AfterPosition(10))
>>> len(loc)
5
>>> for i in loc: print(i)
5
6
7
8
9
>>> list(loc)
[5, 6, 7, 8, 9]
>>> [i for i in range(15) if i in loc]
[5, 6, 7, 8, 9]

Note this is strand aware:

>>> loc = SimpleLocation(BeforePosition(5), AfterPosition(10), strand = -1)
>>> list(loc)
[9, 8, 7, 6, 5]
__eq__(other)

Implement equality by comparing all the location attributes.

property parts

Read only list of sections (always one, the SimpleLocation object).

This is a convenience property allowing you to write code handling both SimpleLocation objects (with one part) and more complex CompoundLocation objects (with multiple parts) interchangeably.

property start

Start location - left most (minimum) value, regardless of strand.

Read only, returns an integer like position object, possibly a fuzzy position.

property end

End location - right most (maximum) value, regardless of strand.

Read only, returns an integer like position object, possibly a fuzzy position.

property nofuzzy_start

Start position (integer, approximated if fuzzy, read only) (DEPRECATED).

This is now an alias for int(feature.start), which should be used in preference – unless you are trying to support old versions of Biopython.

property nofuzzy_end

End position (integer, approximated if fuzzy, read only) (DEPRECATED).

This is now an alias for int(feature.end), which should be used in preference – unless you are trying to support old versions of Biopython.

extract(parent_sequence, references=None)

Extract the sequence from supplied parent sequence using the SimpleLocation object.

The parent_sequence can be a Seq like object or a string, and will generally return an object of the same type. The exception to this is a MutableSeq as the parent sequence will return a Seq object. If the location refers to other records, they must be supplied in the optional dictionary references.

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SimpleLocation
>>> seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> feature_loc = SimpleLocation(8, 15)
>>> feature_loc.extract(seq)
Seq('VALIVIC')
__abstractmethods__ = frozenset({})
__hash__ = None
Bio.SeqFeature.FeatureLocation

alias of Bio.SeqFeature.SimpleLocation

class Bio.SeqFeature.CompoundLocation(parts, operator='join')

Bases: Bio.SeqFeature.Location

For handling joins etc where a feature location has several parts.

__init__(parts, operator='join')

Initialize the class.

>>> from Bio.SeqFeature import SimpleLocation, CompoundLocation
>>> f1 = SimpleLocation(10, 40, strand=+1)
>>> f2 = SimpleLocation(50, 59, strand=+1)
>>> f = CompoundLocation([f1, f2])
>>> len(f) == len(f1) + len(f2) == 39 == len(list(f))
True
>>> print(f.operator)
join
>>> 5 in f
False
>>> 15 in f
True
>>> f.strand
1

Notice that the strand of the compound location is computed automatically - in the case of mixed strands on the sub-locations the overall strand is set to None.

>>> f = CompoundLocation([SimpleLocation(3, 6, strand=+1),
...                       SimpleLocation(10, 13, strand=-1)])
>>> print(f.strand)
None
>>> len(f)
6
>>> list(f)
[3, 4, 5, 12, 11, 10]

The example above doing list(f) iterates over the coordinates within the feature. This allows you to use max and min on the location, to find the range covered:

>>> min(f)
3
>>> max(f)
12

More generally, you can use the compound location’s start and end which give the full span covered, 0 <= start <= end <= full sequence length.

>>> f.start == min(f)
True
>>> f.end == max(f) + 1
True

This is consistent with the behavior of the SimpleLocation for a single region, where again the ‘start’ and ‘end’ do not necessarily give the biological start and end, but rather the ‘minimal’ and ‘maximal’ coordinate boundaries.

Note that adding locations provides a more intuitive method of construction:

>>> f = SimpleLocation(3, 6, strand=+1) + SimpleLocation(10, 13, strand=-1)
>>> len(f)
6
>>> list(f)
[3, 4, 5, 12, 11, 10]
__str__()

Return a representation of the CompoundLocation object (with python counting).

__repr__()

Represent the CompoundLocation object as string for debugging.

property strand

Overall strand of the compound location.

If all the parts have the same strand, that is returned. Otherwise for mixed strands, this returns None.

>>> from Bio.SeqFeature import SimpleLocation, CompoundLocation
>>> f1 = SimpleLocation(15, 17, strand=1)
>>> f2 = SimpleLocation(20, 30, strand=-1)
>>> f = f1 + f2
>>> f1.strand
1
>>> f2.strand
-1
>>> f.strand
>>> f.strand is None
True

If you set the strand of a CompoundLocation, this is applied to all the parts - use with caution:

>>> f.strand = 1
>>> f1.strand
1
>>> f2.strand
1
>>> f.strand
1
__add__(other)

Combine locations, or shift the location by an integer offset.

>>> from Bio.SeqFeature import SimpleLocation
>>> f1 = SimpleLocation(15, 17) + SimpleLocation(20, 30)
>>> print(f1)
join{[15:17], [20:30]}

You can add another SimpleLocation:

>>> print(f1 + SimpleLocation(40, 50))
join{[15:17], [20:30], [40:50]}
>>> print(SimpleLocation(5, 10) + f1)
join{[5:10], [15:17], [20:30]}

You can also add another CompoundLocation:

>>> f2 = SimpleLocation(40, 50) + SimpleLocation(60, 70)
>>> print(f2)
join{[40:50], [60:70]}
>>> print(f1 + f2)
join{[15:17], [20:30], [40:50], [60:70]}

Also, as with the SimpleLocation, adding an integer shifts the location’s coordinates by that offset:

>>> print(f1 + 100)
join{[115:117], [120:130]}
>>> print(200 + f1)
join{[215:217], [220:230]}
>>> print(f1 + (-5))
join{[10:12], [15:25]}
__radd__(other)

Add a feature to the left.

__contains__(value)

Check if an integer position is within the CompoundLocation object.

__nonzero__()

Return True regardless of the length of the feature.

This behavior is for backwards compatibility, since until the __len__ method was added, a SimpleLocation always evaluated as True.

Note that in comparison, Seq objects, strings, lists, etc, will all evaluate to False if they have length zero.

WARNING: The SimpleLocation may in future evaluate to False when its length is zero (in order to better match normal python behavior)!

__len__()

Return the length of the CompoundLocation object.

__iter__()

Iterate over the parent positions within the CompoundLocation object.

__eq__(other)

Check if all parts of CompoundLocation are equal to all parts of other CompoundLocation.

property start

Start location - left most (minimum) value, regardless of strand.

Read only, returns an integer like position object, possibly a fuzzy position.

For the special case of a CompoundLocation wrapping the origin of a circular genome, this will return zero.

property end

End location - right most (maximum) value, regardless of strand.

Read only, returns an integer like position object, possibly a fuzzy position.

For the special case of a CompoundLocation wrapping the origin of a circular genome this will match the genome length (minus one given how Python counts from zero).

property nofuzzy_start

Start position (integer, approximated if fuzzy, read only) (DEPRECATED).

This is an alias for int(feature.start), which should be used in preference – unless you are trying to support old versions of Biopython.

property nofuzzy_end

End position (integer, approximated if fuzzy, read only) (DEPRECATED).

This is an alias for int(feature.end), which should be used in preference – unless you are trying to support old versions of Biopython.

property ref

Not present in CompoundLocation, dummy method for API compatibility.

property ref_db

Not present in CompoundLocation, dummy method for API compatibility.

extract(parent_sequence, references=None)

Extract the sequence from supplied parent sequence using the CompoundLocation object.

The parent_sequence can be a Seq like object or a string, and will generally return an object of the same type. The exception to this is a MutableSeq as the parent sequence will return a Seq object. If the location refers to other records, they must be supplied in the optional dictionary references.

>>> from Bio.Seq import Seq
>>> from Bio.SeqFeature import SimpleLocation, CompoundLocation
>>> seq = Seq("MKQHKAMIVALIVICITAVVAAL")
>>> fl1 = SimpleLocation(2, 8)
>>> fl2 = SimpleLocation(10, 15)
>>> fl3 = CompoundLocation([fl1,fl2])
>>> fl3.extract(seq)
Seq('QHKAMILIVIC')
__abstractmethods__ = frozenset({})
__hash__ = None
class Bio.SeqFeature.Position

Bases: abc.ABC

Abstract base class representing a position.

abstract __repr__()

Represent the Position object as a string for debugging.

property position

Legacy attribute to get (left-most) position as an integer (DEPRECATED).

property extension

Legacy attribute to get the position’s ‘width’ as an integer, typically zero (DEPRECATED).

static fromstring(text, offset=0)

Build a Position object from the text string.

For an end position, leave offset as zero (default):

>>> Position.fromstring("5")
ExactPosition(5)

For a start position, set offset to minus one (for Python counting):

>>> Position.fromstring("5", -1)
ExactPosition(4)

This also covers fuzzy positions:

>>> p = Position.fromstring("<5")
>>> p
BeforePosition(5)
>>> print(p)
<5
>>> int(p)
5
>>> Position.fromstring(">5")
AfterPosition(5)

By default assumes an end position, so note the integer behavior:

>>> p = Position.fromstring("one-of(5,8,11)")
>>> p
OneOfPosition(11, choices=[ExactPosition(5), ExactPosition(8), ExactPosition(11)])
>>> print(p)
one-of(5,8,11)
>>> int(p)
11
>>> Position.fromstring("(8.10)")
WithinPosition(10, left=8, right=10)

Fuzzy start positions:

>>> p = Position.fromstring("<5", -1)
>>> p
BeforePosition(4)
>>> print(p)
<4
>>> int(p)
4

Notice how the integer behavior changes too!

>>> p = Position.fromstring("one-of(5,8,11)", -1)
>>> p
OneOfPosition(4, choices=[ExactPosition(4), ExactPosition(7), ExactPosition(10)])
>>> print(p)
one-of(4,7,10)
>>> int(p)
4
__abstractmethods__ = frozenset({'__repr__'})
class Bio.SeqFeature.ExactPosition(position, extension=0)

Bases: int, Bio.SeqFeature.Position

Specify the specific position of a boundary.

Arguments:
  • position - The position of the boundary.

  • extension - An optional argument which must be zero since we don’t have an extension. The argument is provided so that the same number of arguments can be passed to all position types.

In this case, there is no fuzziness associated with the position.

>>> p = ExactPosition(5)
>>> p
ExactPosition(5)
>>> print(p)
5
>>> isinstance(p, Position)
True
>>> isinstance(p, int)
True

Integer comparisons and operations should work as expected:

>>> p == 5
True
>>> p < 6
True
>>> p <= 5
True
>>> p + 10
ExactPosition(15)
static __new__(cls, position, extension=0)

Create an ExactPosition object.

__str__()

Return a representation of the ExactPosition object (with python counting).

__repr__()

Represent the ExactPosition object as a string for debugging.

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.UncertainPosition(position, extension=0)

Bases: Bio.SeqFeature.ExactPosition

Specify a specific position which is uncertain.

This is used in UniProt, e.g. ?222 for uncertain position 222, or in the XML format explicitly marked as uncertain. Does not apply to GenBank/EMBL.

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.UnknownPosition

Bases: Bio.SeqFeature.Position

Specify a specific position which is unknown (has no position).

This is used in UniProt, e.g. ? or in the XML as unknown.

__repr__()

Represent the UnknownPosition object as a string for debugging.

__hash__()

Return the hash value of the UnknownPosition object.

property position

Legacy attribute to get location (None) (DEPRECATED).

In general you can use the location directly as with the exception of UnknownPosition it subclasses int, or use int(location), rather than this location.position legacy attribute.

However, the UnknownPosition cannot be cast to an integer, and thus does not subclass int, and int(…) will fail. The legacy attribute would return None instead.

Note that while None == None, UnknownPosition() != UnknownPosition() which is like the behavour for NaN.

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.WithinPosition(position, left, right)

Bases: int, Bio.SeqFeature.Position

Specify the position of a boundary within some coordinates.

Arguments: - position - The default integer position - left - The start (left) position of the boundary - right - The end (right) position of the boundary

This allows dealing with a location like ((11.14)..100). This indicates that the start of the sequence is somewhere between 11 and 14. Since this is a start coordinate, it should act like it is at position 11 (or in Python counting, 10).

>>> p = WithinPosition(10, 10, 13)
>>> p
WithinPosition(10, left=10, right=13)
>>> print(p)
(10.13)
>>> int(p)
10

Basic integer comparisons and operations should work as though this were a plain integer:

>>> p == 10
True
>>> p in [9, 10, 11]
True
>>> p < 11
True
>>> p + 10
WithinPosition(20, left=20, right=23)
>>> isinstance(p, WithinPosition)
True
>>> isinstance(p, Position)
True
>>> isinstance(p, int)
True

Note this also applies for comparison to other position objects, where again the integer behavior is used:

>>> p == 10
True
>>> p == ExactPosition(10)
True
>>> p == BeforePosition(10)
True
>>> p == AfterPosition(10)
True

If this were an end point, you would want the position to be 13 (the right/larger value, not the left/smaller value as above):

>>> p2 = WithinPosition(13, 10, 13)
>>> p2
WithinPosition(13, left=10, right=13)
>>> print(p2)
(10.13)
>>> int(p2)
13
>>> p2 == 13
True
>>> p2 == ExactPosition(13)
True
static __new__(cls, position, left, right)

Create a WithinPosition object.

__getnewargs__()

Return the arguments accepted by __new__.

Necessary to allow pickling and unpickling of class instances.

__repr__()

Represent the WithinPosition object as a string for debugging.

__str__()

Return a representation of the WithinPosition object (with python counting).

property position

Legacy attribute to get (left) position as integer (DEPRECATED).

property extension

Legacy attribute to get the within-position’s ‘width’ as an integer (DEPRECATED).

__add__(offset)

Return a copy of the position object with its location shifted.

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.BetweenPosition(position, left, right)

Bases: int, Bio.SeqFeature.Position

Specify the position of a boundary between two coordinates (OBSOLETE?).

Arguments:
  • position - The default integer position

  • left - The start (left) position of the boundary

  • right - The end (right) position of the boundary

This allows dealing with a position like 123^456. This indicates that the start of the sequence is somewhere between 123 and 456. It is up to the parser to set the position argument to either boundary point (depending on if this is being used as a start or end of the feature). For example as a feature end:

>>> p = BetweenPosition(456, 123, 456)
>>> p
BetweenPosition(456, left=123, right=456)
>>> print(p)
(123^456)
>>> int(p)
456

Integer equality and comparison use the given position,

>>> p == 456
True
>>> p in [455, 456, 457]
True
>>> p > 300
True

The old legacy properties of position and extension give the starting/lower/left position as an integer, and the distance to the ending/higher/right position as an integer. Note that the position object will act like either the left or the right end-point depending on how it was created:

>>> p2 = BetweenPosition(123, left=123, right=456)
>>> int(p) == int(p2)
False
>>> p == 456
True
>>> p2 == 123
True

Note this potentially surprising behavior:

>>> BetweenPosition(123, left=123, right=456) == ExactPosition(123)
True
>>> BetweenPosition(123, left=123, right=456) == BeforePosition(123)
True
>>> BetweenPosition(123, left=123, right=456) == AfterPosition(123)
True

i.e. For equality (and sorting) the position objects behave like integers.

static __new__(cls, position, left, right)

Create a new instance in BetweenPosition object.

__getnewargs__()

Return the arguments accepted by __new__.

Necessary to allow pickling and unpickling of class instances.

__repr__()

Represent the BetweenPosition object as a string for debugging.

__str__()

Return a representation of the BetweenPosition object (with python counting).

property position

Legacy attribute to get (left) position as integer (DEPRECATED).

property extension

Legacy attribute to get the between-position’s ‘width’ as an integer (DEPRECATED).

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.BeforePosition(position, extension=0)

Bases: int, Bio.SeqFeature.Position

Specify a position where the actual location occurs before it.

Arguments:
  • position - The upper boundary of where the location can occur.

  • extension - An optional argument which must be zero since we don’t have an extension. The argument is provided so that the same number of arguments can be passed to all position types.

This is used to specify positions like (<10..100) where the location occurs somewhere before position 10.

>>> p = BeforePosition(5)
>>> p
BeforePosition(5)
>>> print(p)
<5
>>> int(p)
5
>>> p + 10
BeforePosition(15)

Note this potentially surprising behavior:

>>> p == ExactPosition(5)
True
>>> p == AfterPosition(5)
True

Just remember that for equality and sorting the position objects act like integers.

static __new__(cls, position, extension=0)

Create a new instance in BeforePosition object.

__repr__()

Represent the location as a string for debugging.

__str__()

Return a representation of the BeforePosition object (with python counting).

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.AfterPosition(position, extension=0)

Bases: int, Bio.SeqFeature.Position

Specify a position where the actual location is found after it.

Arguments:
  • position - The lower boundary of where the location can occur.

  • extension - An optional argument which must be zero since we don’t have an extension. The argument is provided so that the same number of arguments can be passed to all position types.

This is used to specify positions like (>10..100) where the location occurs somewhere after position 10.

>>> p = AfterPosition(7)
>>> p
AfterPosition(7)
>>> print(p)
>7
>>> int(p)
7
>>> p + 10
AfterPosition(17)
>>> isinstance(p, AfterPosition)
True
>>> isinstance(p, Position)
True
>>> isinstance(p, int)
True

Note this potentially surprising behavior:

>>> p == ExactPosition(7)
True
>>> p == BeforePosition(7)
True

Just remember that for equality and sorting the position objects act like integers.

static __new__(cls, position, extension=0)

Create a new instance of the AfterPosition object.

__repr__()

Represent the location as a string for debugging.

__str__()

Return a representation of the AfterPosition object (with python counting).

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

__abstractmethods__ = frozenset({})
class Bio.SeqFeature.OneOfPosition(position, choices)

Bases: int, Bio.SeqFeature.Position

Specify a position where the location can be multiple positions.

This models the GenBank ‘one-of(1888,1901)’ function, and tries to make this fit within the Biopython Position models. If this was a start position it should act like 1888, but as an end position 1901.

>>> p = OneOfPosition(1888, [ExactPosition(1888), ExactPosition(1901)])
>>> p
OneOfPosition(1888, choices=[ExactPosition(1888), ExactPosition(1901)])
>>> int(p)
1888

Integer comparisons and operators act like using int(p),

>>> p == 1888
True
>>> p <= 1888
True
>>> p > 1888
False
>>> p + 100
OneOfPosition(1988, choices=[ExactPosition(1988), ExactPosition(2001)])
>>> isinstance(p, OneOfPosition)
True
>>> isinstance(p, Position)
True
>>> isinstance(p, int)
True
static __new__(cls, position, choices)

Initialize with a set of possible positions.

choices is a list of Position derived objects, specifying possible locations.

position is an integer specifying the default behavior.

__getnewargs__()

Return the arguments accepted by __new__.

Necessary to allow pickling and unpickling of class instances.

property position

Legacy attribute to get (left) position as integer (DEPRECATED).

property extension

Legacy attribute to get the one-of-position’s ‘width’ as an integer (DEPRECATED).

__abstractmethods__ = frozenset({})
__repr__()

Represent the OneOfPosition object as a string for debugging.

__str__()

Return a representation of the OneOfPosition object (with python counting).

__add__(offset)

Return a copy of the position object with its location shifted (PRIVATE).

class Bio.SeqFeature.PositionGap(gap_size)

Bases: object

Simple class to hold information about a gap between positions (DEPRECATED).

__init__(gap_size)

Initialize with a position object containing the gap information.

__repr__()

Represent the position gap as a string for debugging.

__str__()

Return a representation of the PositionGap object (with python counting).