Sophie

Sophie

distrib > * > cooker > x86_64 > by-pkgid > 635dc0b7819f4e396a16d64269572c71 > files > 1775

biopython-doc-1.58-1.x86_64.rpm

<?xml version="1.0" encoding="ascii"?>
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN"
          "DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml" xml:lang="en" lang="en">
<head>
  <title>Bio.SeqFeature.SeqFeature</title>
  <link rel="stylesheet" href="epydoc.css" type="text/css" />
  <script type="text/javascript" src="epydoc.js"></script>
</head>

<body bgcolor="white" text="black" link="blue" vlink="#204080"
      alink="#204080">
<!-- ==================== NAVIGATION BAR ==================== -->
<table class="navbar" border="0" width="100%" cellpadding="0"
       bgcolor="#a0c0ff" cellspacing="0">
  <tr valign="middle">

  <!-- Tree link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="module-tree.html">Trees</a>&nbsp;&nbsp;&nbsp;</th>

  <!-- Index link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="identifier-index.html">Indices</a>&nbsp;&nbsp;&nbsp;</th>

  <!-- Help link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="help.html">Help</a>&nbsp;&nbsp;&nbsp;</th>

      <th class="navbar" width="100%"></th>
  </tr>
</table>
<table width="100%" cellpadding="0" cellspacing="0">
  <tr valign="top">
    <td width="100%">
      <span class="breadcrumbs">
        <a href="Bio-module.html">Package&nbsp;Bio</a> ::
        <a href="Bio.SeqFeature-module.html">Module&nbsp;SeqFeature</a> ::
        Class&nbsp;SeqFeature
      </span>
    </td>
    <td>
      <table cellpadding="0" cellspacing="0">
        <!-- hide/show private -->
        <tr><td align="right"><span class="options">[<a href="javascript:void(0);" class="privatelink"
    onclick="toggle_private();">hide&nbsp;private</a>]</span></td></tr>
        <tr><td align="right"><span class="options"
            >[<a href="frames.html" target="_top">frames</a
            >]&nbsp;|&nbsp;<a href="Bio.SeqFeature.SeqFeature-class.html"
            target="_top">no&nbsp;frames</a>]</span></td></tr>
      </table>
    </td>
  </tr>
</table>
<!-- ==================== CLASS DESCRIPTION ==================== -->
<h1 class="epydoc">Class SeqFeature</h1><p class="nomargin-top"><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature">source&nbsp;code</a></span></p>
<pre class="base-tree">
object --+
         |
        <strong class="uidshort">SeqFeature</strong>
</pre>

<hr />
<p>Represent a Sequence Feature on an object.</p>
  <p>Attributes: o location - the location of the feature on the sequence 
  (FeatureLocation) o type - the specified type of the feature (ie. CDS, 
  exon, repeat...) o 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 &quot;join&quot; o 
  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) o id - A string identifier for the 
  feature. o ref - A reference to another sequence. This could be an 
  accession number for some different sequence. o ref_db - A different 
  database for the reference accession number. o qualifiers - A dictionary 
  of qualifiers on the feature. These are analagous to the qualifiers from 
  a GenBank feature table. The keys of the dictionary are qualifier names, 
  the values are the qualifier values. o sub_features - Additional 
  SeqFeatures which fall under this 'parent' feature. For instance, if we 
  having something like:</p>
  <p>CDS    join(1..10,30..40,50..60)</p>
  <p>Then the top level feature would be of type 'CDS' from 1 to 60 
  (actually 0 to 60 in Python counting) with location_operator='join', and 
  the three sub- features would also be of type 'CDS', and would be from 1 
  to 10, 30 to 40 and 50 to 60, respectively (although actually using 
  Python counting).</p>
  <p>To get the nucleotide sequence for this CDS, you would need to take 
  the parent sequence and do seq[0:10]+seq[29:40]+seq[49:60] (Python 
  counting). Things are more complicated with strands and fuzzy positions. 
  To save you dealing with all these special cases, the SeqFeature provides
  an extract method to do this for you.</p>

<!-- ==================== INSTANCE METHODS ==================== -->
<a name="section-InstanceMethods"></a>
<table class="summary" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr bgcolor="#70b0f0" class="table-header">
  <td colspan="2" class="table-header">
    <table border="0" cellpadding="0" cellspacing="0" width="100%">
      <tr valign="top">
        <td align="left"><span class="table-header">Instance Methods</span></td>
        <td align="right" valign="top"
         ><span class="options">[<a href="#section-InstanceMethods"
         class="privatelink" onclick="toggle_private();"
         >hide private</a>]</span></td>
      </tr>
    </table>
  </td>
</tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__init__" class="summary-sig-name">__init__</a>(<span class="summary-sig-arg">self</span>,
        <span class="summary-sig-arg">location</span>=<span class="summary-sig-default">None</span>,
        <span class="summary-sig-arg">type</span>=<span class="summary-sig-default"><code class="variable-quote">'</code><code class="variable-string"></code><code class="variable-quote">'</code></span>,
        <span class="summary-sig-arg">location_operator</span>=<span class="summary-sig-default"><code class="variable-quote">'</code><code class="variable-string"></code><code class="variable-quote">'</code></span>,
        <span class="summary-sig-arg">strand</span>=<span class="summary-sig-default">None</span>,
        <span class="summary-sig-arg">id</span>=<span class="summary-sig-default"><code class="variable-quote">'</code><code class="variable-string">&lt;unknown id&gt;</code><code class="variable-quote">'</code></span>,
        <span class="summary-sig-arg">qualifiers</span>=<span class="summary-sig-default">None</span>,
        <span class="summary-sig-arg">sub_features</span>=<span class="summary-sig-default">None</span>,
        <span class="summary-sig-arg">ref</span>=<span class="summary-sig-default">None</span>,
        <span class="summary-sig-arg">ref_db</span>=<span class="summary-sig-default">None</span>)</span><br />
      Initialize a SeqFeature on a Sequence.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__init__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__repr__" class="summary-sig-name">__repr__</a>(<span class="summary-sig-arg">self</span>)</span><br />
      A string representation of the record for debugging.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__repr__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__str__" class="summary-sig-name">__str__</a>(<span class="summary-sig-arg">self</span>)</span><br />
      A readable summary of the feature intended to be printed to screen.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__str__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr class="private">
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#_shift" class="summary-sig-name" onclick="show_private();">_shift</a>(<span class="summary-sig-arg">self</span>,
        <span class="summary-sig-arg">offset</span>)</span><br />
      Returns a copy of the feature with its location shifted (PRIVATE).</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature._shift">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr class="private">
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#_flip" class="summary-sig-name" onclick="show_private();">_flip</a>(<span class="summary-sig-arg">self</span>,
        <span class="summary-sig-arg">length</span>)</span><br />
      Returns a copy of the feature with its location flipped (PRIVATE).</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature._flip">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#extract" class="summary-sig-name">extract</a>(<span class="summary-sig-arg">self</span>,
        <span class="summary-sig-arg">parent_sequence</span>)</span><br />
      Extract feature sequence from the supplied parent sequence.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.extract">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__nonzero__" class="summary-sig-name">__nonzero__</a>(<span class="summary-sig-arg">self</span>)</span><br />
      Returns True regardless of the length of the feature.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__nonzero__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__len__" class="summary-sig-name">__len__</a>(<span class="summary-sig-arg">self</span>)</span><br />
      Returns the length of the region described by a feature.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__len__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__iter__" class="summary-sig-name">__iter__</a>(<span class="summary-sig-arg">self</span>)</span><br />
      Iterate over the parent positions within the feature.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__iter__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
<tr>
    <td width="15%" align="right" valign="top" class="summary">
      <span class="summary-type">&nbsp;</span>
    </td><td class="summary">
      <table width="100%" cellpadding="0" cellspacing="0" border="0">
        <tr>
          <td><span class="summary-sig"><a href="Bio.SeqFeature.SeqFeature-class.html#__contains__" class="summary-sig-name">__contains__</a>(<span class="summary-sig-arg">self</span>,
        <span class="summary-sig-arg">value</span>)</span><br />
      Check if an integer position is within the feature.</td>
          <td align="right" valign="top">
            <span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__contains__">source&nbsp;code</a></span>
            
          </td>
        </tr>
      </table>
      
    </td>
  </tr>
  <tr>
    <td colspan="2" class="summary">
    <p class="indent-wrapped-lines"><b>Inherited from <code>object</code></b>:
      <code>__delattr__</code>,
      <code>__format__</code>,
      <code>__getattribute__</code>,
      <code>__hash__</code>,
      <code>__new__</code>,
      <code>__reduce__</code>,
      <code>__reduce_ex__</code>,
      <code>__setattr__</code>,
      <code>__sizeof__</code>,
      <code>__subclasshook__</code>
      </p>
    </td>
  </tr>
</table>
<!-- ==================== PROPERTIES ==================== -->
<a name="section-Properties"></a>
<table class="summary" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr bgcolor="#70b0f0" class="table-header">
  <td colspan="2" class="table-header">
    <table border="0" cellpadding="0" cellspacing="0" width="100%">
      <tr valign="top">
        <td align="left"><span class="table-header">Properties</span></td>
        <td align="right" valign="top"
         ><span class="options">[<a href="#section-Properties"
         class="privatelink" onclick="toggle_private();"
         >hide private</a>]</span></td>
      </tr>
    </table>
  </td>
</tr>
  <tr>
    <td colspan="2" class="summary">
    <p class="indent-wrapped-lines"><b>Inherited from <code>object</code></b>:
      <code>__class__</code>
      </p>
    </td>
  </tr>
</table>
<!-- ==================== METHOD DETAILS ==================== -->
<a name="section-MethodDetails"></a>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr bgcolor="#70b0f0" class="table-header">
  <td colspan="2" class="table-header">
    <table border="0" cellpadding="0" cellspacing="0" width="100%">
      <tr valign="top">
        <td align="left"><span class="table-header">Method Details</span></td>
        <td align="right" valign="top"
         ><span class="options">[<a href="#section-MethodDetails"
         class="privatelink" onclick="toggle_private();"
         >hide private</a>]</span></td>
      </tr>
    </table>
  </td>
</tr>
</table>
<a name="__init__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__init__</span>(<span class="sig-arg">self</span>,
        <span class="sig-arg">location</span>=<span class="sig-default">None</span>,
        <span class="sig-arg">type</span>=<span class="sig-default"><code class="variable-quote">'</code><code class="variable-string"></code><code class="variable-quote">'</code></span>,
        <span class="sig-arg">location_operator</span>=<span class="sig-default"><code class="variable-quote">'</code><code class="variable-string"></code><code class="variable-quote">'</code></span>,
        <span class="sig-arg">strand</span>=<span class="sig-default">None</span>,
        <span class="sig-arg">id</span>=<span class="sig-default"><code class="variable-quote">'</code><code class="variable-string">&lt;unknown id&gt;</code><code class="variable-quote">'</code></span>,
        <span class="sig-arg">qualifiers</span>=<span class="sig-default">None</span>,
        <span class="sig-arg">sub_features</span>=<span class="sig-default">None</span>,
        <span class="sig-arg">ref</span>=<span class="sig-default">None</span>,
        <span class="sig-arg">ref_db</span>=<span class="sig-default">None</span>)</span>
    <br /><em class="fname">(Constructor)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__init__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Initialize a SeqFeature on a Sequence.</p>
  <p>location can either be a FeatureLocation (with strand argument also 
  given if required), or None.</p>
  <p>e.g. With no strand, on the forward strand, and on the reverse 
  strand:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span>f1 = SeqFeature(FeatureLocation(5,10), type=<span class="py-string">&quot;domain&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span>f2 = SeqFeature(FeatureLocation(7,110), strand=1, type=<span class="py-string">&quot;CDS&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span>f3 = SeqFeature(FeatureLocation(9,108), strand=-1, type=<span class="py-string">&quot;CDS&quot;</span>)</pre>
  <p>An invalid strand will trigger an exception:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span>f4 = SeqFeature(FeatureLocation(50,60), strand=2)
<span class="py-except">Traceback (most recent call last):</span>
<span class="py-except">   ...</span>
<span class="py-except">ValueError: Strand should be +1, -1, 0 or None, not 2</span></pre>
  <p>For exact start/end positions, an integer can be used (as shown above)
  as shorthand for the ExactPosition object. For non-exact locations, the 
  FeatureLocation must be specified via the appropriate position 
  objects.</p>
  <dl class="fields">
    <dt>Overrides:
        object.__init__
    </dt>
  </dl>
</td></tr></table>
</div>
<a name="__repr__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__repr__</span>(<span class="sig-arg">self</span>)</span>
    <br /><em class="fname">(Representation operator)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__repr__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>A string representation of the record for debugging.</p>
  <dl class="fields">
    <dt>Overrides:
        object.__repr__
    </dt>
  </dl>
</td></tr></table>
</div>
<a name="__str__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__str__</span>(<span class="sig-arg">self</span>)</span>
    <br /><em class="fname">(Informal representation operator)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__str__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>A readable summary of the feature intended to be printed to 
  screen.</p>
  <dl class="fields">
    <dt>Overrides:
        object.__str__
    </dt>
  </dl>
</td></tr></table>
</div>
<a name="_shift"></a>
<div class="private">
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">_shift</span>(<span class="sig-arg">self</span>,
        <span class="sig-arg">offset</span>)</span>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature._shift">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Returns a copy of the feature with its location shifted (PRIVATE).</p>
  <p>The annotation qaulifiers are copied.</p>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="_flip"></a>
<div class="private">
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">_flip</span>(<span class="sig-arg">self</span>,
        <span class="sig-arg">length</span>)</span>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature._flip">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Returns a copy of the feature with its location flipped (PRIVATE).</p>
  <p>The argument length gives the length of the parent sequence. For 
  example a location 0..20 (+1 strand) with parent length 30 becomes after 
  flipping 10..30 (-1 strand). Strandless (None) or unknown strand (0) 
  remain like that - just their end points are changed.</p>
  <p>The annotation qaulifiers are copied.</p>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="extract"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">extract</span>(<span class="sig-arg">self</span>,
        <span class="sig-arg">parent_sequence</span>)</span>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.extract">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Extract feature sequence from the supplied parent sequence.</p>
  <p>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.</p>
  <p>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.</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.Seq <span class="py-keyword">import</span> Seq
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.Alphabet <span class="py-keyword">import</span> generic_protein
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span>seq = Seq(<span class="py-string">&quot;MKQHKAMIVALIVICITAVVAAL&quot;</span>, generic_protein)
<span class="py-prompt">&gt;&gt;&gt; </span>f = SeqFeature(FeatureLocation(8,15), type=<span class="py-string">&quot;domain&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span>f.extract(seq)
<span class="py-output">Seq('VALIVIC', ProteinAlphabet())</span></pre>
  <p>Note - currently only sub-features of type &quot;join&quot; are 
  supported.</p>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="__nonzero__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__nonzero__</span>(<span class="sig-arg">self</span>)</span>
    <br /><em class="fname">(Boolean test operator)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__nonzero__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Returns True regardless of the length of the feature.</p>
  <p>This behaviour is for backwards compatibility, since until the __len__
  method was added, a SeqFeature always evaluated as True.</p>
  <p>Note that in comparison, Seq objects, strings, lists, etc, will all 
  evaluate to False if they have length zero.</p>
  <p>WARNING: The SeqFeature may in future evaluate to False when its 
  length is zero (in order to better match normal python behaviour)!</p>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="__len__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__len__</span>(<span class="sig-arg">self</span>)</span>
    <br /><em class="fname">(Length operator)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__len__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Returns the length of the region described by a feature.</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.Seq <span class="py-keyword">import</span> Seq
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.Alphabet <span class="py-keyword">import</span> generic_protein
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span>seq = Seq(<span class="py-string">&quot;MKQHKAMIVALIVICITAVVAAL&quot;</span>, generic_protein)
<span class="py-prompt">&gt;&gt;&gt; </span>f = SeqFeature(FeatureLocation(8,15), type=<span class="py-string">&quot;domain&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span>len(f)
<span class="py-output">7</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span>f.extract(seq)
<span class="py-output">Seq('VALIVIC', ProteinAlphabet())</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span>len(f.extract(seq))
<span class="py-output">7</span></pre>
  <p>For simple features without subfeatures this is the same as the region
  spanned (end position minus start position). However, for a feature 
  defined by combining several subfeatures (e.g. a CDS as the join of 
  several exons) the gaps are not counted (e.g. introns). This ensures that
  len(f) == len(f.extract(parent_seq)), and also makes sure things work 
  properly with features wrapping the origin etc.</p>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="__iter__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__iter__</span>(<span class="sig-arg">self</span>)</span>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__iter__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Iterate over the parent positions within the feature.</p>
  <p>The iteration order is strand aware, and can be thought of as moving 
  along the feature using the parent sequence coordinates:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span>f = SeqFeature(FeatureLocation(5,10), type=<span class="py-string">&quot;domain&quot;</span>, strand=-1)
<span class="py-prompt">&gt;&gt;&gt; </span>len(f)
<span class="py-output">5</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">for</span> i <span class="py-keyword">in</span> f: <span class="py-keyword">print</span> i
<span class="py-output">9</span>
<span class="py-output">8</span>
<span class="py-output">7</span>
<span class="py-output">6</span>
<span class="py-output">5</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span>list(f)
<span class="py-output">[9, 8, 7, 6, 5]</span></pre>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<a name="__contains__"></a>
<div>
<table class="details" border="1" cellpadding="3"
       cellspacing="0" width="100%" bgcolor="white">
<tr><td>
  <table width="100%" cellpadding="0" cellspacing="0" border="0">
  <tr valign="top"><td>
  <h3 class="epydoc"><span class="sig"><span class="sig-name">__contains__</span>(<span class="sig-arg">self</span>,
        <span class="sig-arg">value</span>)</span>
    <br /><em class="fname">(In operator)</em>
  </h3>
  </td><td align="right" valign="top"
    ><span class="codelink"><a href="Bio.SeqFeature-pysrc.html#SeqFeature.__contains__">source&nbsp;code</a></span>&nbsp;
    </td>
  </tr></table>
  
  <p>Check if an integer position is within the feature.</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span>f = SeqFeature(FeatureLocation(5,10), type=<span class="py-string">&quot;domain&quot;</span>, strand=-1)
<span class="py-prompt">&gt;&gt;&gt; </span>len(f)
<span class="py-output">5</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span>[i <span class="py-keyword">for</span> i <span class="py-keyword">in</span> range(15) <span class="py-keyword">if</span> i <span class="py-keyword">in</span> f]
<span class="py-output">[5, 6, 7, 8, 9]</span></pre>
  <p>For example, to see which features include a SNP position, you could 
  use this:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio <span class="py-keyword">import</span> SeqIO
<span class="py-prompt">&gt;&gt;&gt; </span>record = SeqIO.read(<span class="py-string">&quot;GenBank/NC_000932.gb&quot;</span>, <span class="py-string">&quot;gb&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">for</span> f <span class="py-keyword">in</span> record.features:
<span class="py-more">... </span>    <span class="py-keyword">if</span> 1750 <span class="py-keyword">in</span> f:
<span class="py-more">... </span>        <span class="py-keyword">print</span> f.type, f.strand, f.location
<span class="py-output">source 1 [0:154478]</span>
<span class="py-output">gene -1 [1716:4347]</span>
<span class="py-output">tRNA -1 [1716:4347]</span></pre>
  <p>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:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">for</span> f <span class="py-keyword">in</span> record.features:
<span class="py-more">... </span>    <span class="py-keyword">if</span> 1760 <span class="py-keyword">in</span> f:
<span class="py-more">... </span>        <span class="py-keyword">print</span> f.type, f.strand, f.location
<span class="py-output">source 1 [0:154478]</span>
<span class="py-output">gene -1 [1716:4347]</span></pre>
  <p>Note that additional care may be required with fuzzy locations, for 
  example just before a BeforePosition:</p>
<pre class="py-doctest">
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> SeqFeature, FeatureLocation
<span class="py-prompt">&gt;&gt;&gt; </span><span class="py-keyword">from</span> Bio.SeqFeature <span class="py-keyword">import</span> BeforePosition
<span class="py-prompt">&gt;&gt;&gt; </span>f = SeqFeature(FeatureLocation(BeforePosition(3),8), type=<span class="py-string">&quot;domain&quot;</span>)
<span class="py-prompt">&gt;&gt;&gt; </span>len(f)
<span class="py-output">5</span>
<span class="py-output"></span><span class="py-prompt">&gt;&gt;&gt; </span>[i <span class="py-keyword">for</span> i <span class="py-keyword">in</span> range(10) <span class="py-keyword">if</span> i <span class="py-keyword">in</span> f]
<span class="py-output">[3, 4, 5, 6, 7]</span></pre>
  <dl class="fields">
  </dl>
</td></tr></table>
</div>
<br />
<!-- ==================== NAVIGATION BAR ==================== -->
<table class="navbar" border="0" width="100%" cellpadding="0"
       bgcolor="#a0c0ff" cellspacing="0">
  <tr valign="middle">

  <!-- Tree link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="module-tree.html">Trees</a>&nbsp;&nbsp;&nbsp;</th>

  <!-- Index link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="identifier-index.html">Indices</a>&nbsp;&nbsp;&nbsp;</th>

  <!-- Help link -->
      <th>&nbsp;&nbsp;&nbsp;<a
        href="help.html">Help</a>&nbsp;&nbsp;&nbsp;</th>

      <th class="navbar" width="100%"></th>
  </tr>
</table>
<table border="0" cellpadding="0" cellspacing="0" width="100%%">
  <tr>
    <td align="left" class="footer">
    Generated by Epydoc 3.0.1 on Thu Aug 18 18:22:16 2011
    </td>
    <td align="right" class="footer">
      <a target="mainFrame" href="http://epydoc.sourceforge.net"
        >http://epydoc.sourceforge.net</a>
    </td>
  </tr>
</table>

<script type="text/javascript">
  <!--
  // Private objects are initially displayed (because if
  // javascript is turned off then we want them to be
  // visible); but by default, we want to hide them.  So hide
  // them unless we have a cookie that says to show them.
  checkCookie();
  // -->
</script>
</body>
</html>