This module provides functions for PCR. Primers with 5’ tails as well as inverse PCR on circular templates are handled correctly.
Bases: pydna.dsdna.Dseqrecord
The Amplicon class holds information about a PCR reaction involving two primers and one template. This class is used by the Anneal class and is not meant to be instantiated directly.
Parameters: | forward_primer : SeqRecord(Biopython)
reverse_primer : SeqRecord(Biopython)
template : Dseqrecord
saltc : float, optional
forward_primer_concentration : float, optional
rc : float, optional
|
---|
This method returns a simple figure of the two primers binding to a part of the template.
5gctactacacacgtactgactg3
|||||||||||||||||||||| tm 52.6 (dbd) 58.3
5gctactacacacgtactgactg...caagatagagtcagtaaccaca3
3cgatgatgtgtgcatgactgac...gttctatctcagtcattggtgt5
|||||||||||||||||||||| tm 49.1 (dbd) 57.7
3gttctatctcagtcattggtgt5
Returns: | figure:string :
|
---|
Notes
tm in the figure above is the melting temperature (tm) for each primer calculated according to SantaLucia 1998 [1].
dbd is the tm calculation for enzymes with dsDNA binding domains like Pfu-Sso7d [2]. See [3] for more information.
References
[1] |
|
[2] |
|
[3] | http://www.thermoscientificbio.com/webtools/tmc/ |
Returns a Dseqrecord object containing flankdnlength bases downstream of the reverse primer footprint. Truncated if the template is not long enough.
<---- flankdn ------>
3actactgactatctTAATAA5
||||||||||||||
acgcattcagctactgtactactgactatctatcgtacatgtactatcgtat
Returns a Dseqrecord object containing flankuplength bases upstream of the forward primer footprint, Truncated if the template is not long enough.
<--- flankup --->
5TAATAAactactgactatct3
||||||||||||||
acgcattcagctactgtactactgactatctatcg
Returns a string containing a text representation of two proposed PCR programs. The first program is adapted for Taq poymerase, while the second is adapted for Pfu-Sso7d or other polymerases with DNA binding domain.
Taq (rate 30 nt/s)
Three-step| 30 cycles | |SantaLucia 1998
94.0°C |94.0°C | |SaltC 50mM
__________|_____ 72.0°C |72.0°C|
04min00s |30s \ ________|______|
| \ 46.0°C/ 0min 1s|10min |
| \_____/ | |
| 30s | |4-8°C
Pfu-Sso7d (rate 15s/kb)
Three-step| 30 cycles | |Breslauer1986,SantaLucia1998
98.0°C |98.0°C | |SaltC 50mM
__________|_____ 72.0°C |72.0°C|Primer1C 1µM
00min30s |10s \ 61.0°C ________|______|Primer2C 1µM
| \______/ 0min 0s|10min |
| 10s | |4-8°C
Bases: object
Parameters: | primers : iterable containing SeqRecord objects
template : Dseqrecord object
limit : int, optional
|
---|
Examples
>>> import pydna
>>> template = pydna.Dseqrecord("tacactcaccgtctatcattatctactatcgactgtatcatctgatagcac")
>>> from Bio.SeqRecord import SeqRecord
>>> p1 = pydna.read(">p1\ntacactcaccgtctatcattatc", ds = False)
>>> p2 = pydna.read(">p2\ngtgctatcagatgatacagtcg", ds = False)
>>> ann = pydna.Anneal((p1,p2), template)
>>> print ann.report()
Template <unknown name> 51 nt linear:
Primer p1 anneals forward at position 23
Primer p2 anneals reverse at position 29
>>> ann.products
[Amplicon(51)]
>>> amplicon_list = ann.products
>>> amplicon = amplicon_list.pop()
>>> amplicon
Amplicon(51)
>>> print amplicon.figure()
5tacactcaccgtctatcattatc...cgactgtatcatctgatagcac3
|||||||||||||||||||||| tm 50.6 (dbd) 60.5
3gctgacatagtagactatcgtg5
5tacactcaccgtctatcattatc3
||||||||||||||||||||||| tm 49.4 (dbd) 58.8
3atgtgagtggcagatagtaatag...gctgacatagtagactatcgtg5
>>> amplicon.annotations['date'] = '02-FEB-2013' # Set the date for this example to pass the doctest
>>> print amplicon
Dseqrecord
circular: False
size: 51
ID: 51bp U96+TO06Y6pFs74SQx8M1IVTBiY
Name: 51bp_PCR_prod
Description: Product_p1_p2
Number of features: 2
/date=02-FEB-2013
Dseq(-51)
taca..gcac
atgt..cgtg
>>> print amplicon.program()
Taq (rate 30 nt/s)
Three-step| 30 cycles | |SantaLucia 1998
94.0°C |94.0°C | |SaltC 50mM
__________|_____ 72.0°C |72.0°C|
04min00s |30s \ ________|______|
| \ 44.0°C/ 0min 1s|10min |
| \_____/ | |
| 30s | |4-8°C
Pfu-Sso7d (rate 15s/kb)
Three-step| 30 cycles | |Breslauer1986,SantaLucia1998
98.0°C |98.0°C | |SaltC 50mM
__________|_____ 72.0°C |72.0°C|Primer1C 1µM
00min30s |10s \ 62.0°C ________|______|Primer2C 1µM
| \______/ 0min 0s|10min |
| 10s | |4-8°C
>>>
Attributes
products: list | A list of Amplicon objects, one for each primer pair that may form a PCR product. |
Bases: Bio.SeqRecord.SeqRecord
This class holds information about a primer on a template
Returns the melting temperature (Tm) of the primer using the basic formula.
Parameters: | primer : string
|
---|---|
Returns: | tm : int |
Examples
>>> from pydna.amplify import basictm
>>> basictm("ggatcc")
20
>>>
pcr is a convenience function for Anneal to simplify its usage, especially from the command line. If more than one PCR product is formed, an exception is raised.
args is any iterable of sequences or an iterable of iterables of sequences. args will be greedily flattened.
Parameters: | args : iterable containing sequence objects
limit : int = 13, optional
|
---|---|
Returns: | product : Dseqrecord
|
Notes
sequences in args could be of type:
The last sequence will be interpreted as the template all preceeding sequences as primers.
This is a powerful function, use with care!
Examples
>>> import pydna
>>> template = pydna.Dseqrecord("tacactcaccgtctatcattatctactatcgactgtatcatctgatagcac")
>>> from Bio.SeqRecord import SeqRecord
>>> p1 = pydna.read(">p1\ntacactcaccgtctatcattatc", ds = False)
>>> p2 = pydna.read(">p2\ncgactgtatcatctgatagcac", ds = False).reverse_complement()
>>> pydna.pcr(p1, p2, template)
Amplicon(51)
>>> pydna.pcr([p1, p2], template)
Amplicon(51)
>>> pydna.pcr((p1,p2,), template)
Amplicon(51)
>>>
Returns the melting temperature (Tm) of the primer using the nearest neighbour algorithm. Formula and thermodynamic data is taken from Breslauer 1986. These data are no longer widely used.
Breslauer 1986, table 2 [4]
pair | dH | dS | dG |
---|---|---|---|
AA/TT | 9.1 | 24.0 | 1.9 |
AT/TA | 8.6 | 23.9 | 1.5 |
TA/AT | 6.0 | 16.9 | 0.9 |
CA/GT | 5.8 | 12.9 | 1.9 |
GT/CA | 6.5 | 17.3 | 1.3 |
CT/GA | 7.8 | 20.8 | 1.6 |
GA/CT | 5.6 | 13.5 | 1.6 |
CG/GC | 11.9 | 27.8 | 3.6 |
GC/CG | 11.1 | 26.7 | 3.1 |
GG/CC | 11.0 | 26.6 | 3.1 |
Parameters: | primer : string
|
---|---|
Returns: | tm : float |
References
[4] | K.J. Breslauer et al., “Predicting DNA Duplex Stability from the Base Sequence,” Proceedings of the National Academy of Sciences 83, no. 11 (1986): 3746. |
Examples
>>> from pydna.amplify import tmbreslauer86
>>> tmbreslauer86("ACGTCATCGACACTATCATCGAC")
64.28863985851899
Returns the tm for a primer using a formula adapted to polymerases with a DNA binding domain.
Parameters: | primer : string
primerc : float
saltc : float, optional
thermodynamics : bool, optional
|
---|---|
Returns: | tm : float
tm,dH,dS : tuple
|
Returns the melting temperature (Tm) of the primer using the nearest neighbour algorithm. Formula and thermodynamic data is taken from SantaLucia 1998 [1]. This implementation gives the same answer as the one provided by Biopython (See Examples).
Thermodynamic data used:
pair | dH | dS |
---|---|---|
AA/TT | 7.9 | 22.2 |
AT/TA | 7.2 | 20.4 |
TA/AT | 7.2 | 21.3 |
CA/GT | 8.5 | 22.7 |
GT/CA | 8.4 | 22.4 |
CT/GA | 7.8 | 21.0 |
GA/CT | 8.2 | 22.2 |
CG/GC | 10.6 | 27.2 |
GC/CG | 9.8 | 24.4 |
GG/CC | 8.0 | 19.9 |
Parameters: | primer : string
|
---|---|
Returns: | tm : float
|
Examples
>>> from pydna.amplify import tmstaluc98
>>> from Bio.SeqUtils.MeltingTemp import Tm_staluc
>>> tmstaluc98("ACGTCATCGACACTATCATCGAC")
54.55597724052518
>>> Tm_staluc("ACGTCATCGACACTATCATCGAC")
54.555977240525124
>>>
This module provides functions for assembly of sequences by homologous recombination and other related techniques. Given a list of sequences (Dseqrecords), all sequences will be analyzed for overlapping regions of DNA (common substrings).
The assembly algorithm is based on graph theory where each overlapping region forms a node and sequences separating the overlapping regions form edges.
Bases: object
Assembly of a list of linear DNA fragments into linear or circular constructs. Accepts a list of Dseqrecords (source fragments) to initiate an Assembly object. Several methods are available for analysis of overlapping sequences, graph construction and assembly.
Parameters: | dsrecs : list
|
---|
Examples
>>> from pydna import Assembly, Dseqrecord
>>> a = Dseqrecord("acgatgctatactgCCCCCtgtgctgtgctcta")
>>> b = Dseqrecord("tgtgctgtgctctaTTTTTtattctggctgtatc")
>>> c = Dseqrecord("tattctggctgtatcGGGGGtacgatgctatactg")
>>> x = Assembly((a,b,c), limit=14)
>>> x.analyze_overlaps()
'6 sequences analyzed of which 3 have shared homologies with totally 3 overlaps'
>>> x.create_graph()
"A graph with 5', 3' and 3 internal nodes was created"
>>> x.assemble_circular_from_graph()
'1 circular products were formed'
>>> x
Assembly object:
Sequences........................: 33, 34, 35
Sequences with shared homologies.: 33, 34, 35
Homology limit (bp)..............: 14
Number of overlaps...............: 3
Nodes in graph...................: 5
Assembly protocol................: overlaps
Circular products................: 59
Linear products..................: No linear products
>>> x.circular_products
[Dseqrecord(o59)]
>>> x.circular_products[0].seq.watson
'CCCCCtgtgctgtgctctaTTTTTtattctggctgtatcGGGGGtacgatgctatactg'
This method carries out a pairwise analysis for common sub strings (shared homologous sequences) among DNA fragments stored in the dsrecs property.
The analyzed_dsrec property will contain copies of the sequences in the dsrecs property that have common substrings. Common sub strings are added as features to the sequences in the analyzed_dsrec property.
Parameters: | limit : int
|
---|
This method carries out a pairwise analysis for shared sub sequences equal to or longer than limit. The difference between this method and analyze_overlaps() is that shared sequences must start or end at the extremities of the source fragments.
Parameters: | limit : int
|
---|
Assembles all circular products from the graph created by create_graph()
Convenience method to assembles all possible linear fusion PCR products. This method is a synonym of assemble_gibson_linear()
Convenience method to assembles all possible circular Gibson assembly products expected from in-vivo homologous recombination by the following steps:
The property “circular_products” is initiated with a list of Fragment objects.
Convenience method to assembles all possible linear Gibson assembly products expected from in-vivo homologous recombination by the following steps:
The property “linear_products” is initiated with a list of Fragment objects.
Convenience method to assembles all possible circular homologous recombination (hr) products expected from in-vivo homologous recombination by the following steps:
The property “circular_products” is initiated with a list of Fragment objects.
Convenience method to assembles all possible linear homologous recombination (hr) products expected from in-vivo homologous recombination by the following steps:
The property “linear_products” is initiated with a list of Fragment objects.
Assembles all linear products from the graph created by create_graph()
Creates a graph from analyzed source fragments. Nodes in the graph are shared sub sequences found as features in the source sequences. Edges are the sequences that connect the shared sub sequences in each source fragment. Two additional nodes called 5’ and 3’ are also added to represent DNA ends.
Sequences fed to this class is stored in this property
The shortest common sub strings to be considered
The max number of nodes allowed. This can be reset to some other value
Bases: pydna.dsdna.Dseqrecord
This class holds information about a DNA assembly. This class is instantiated by the Assembly class and is not meant to be instantiated directly.
Synonym of detailed_fig()
Returns a small ascii representation of the assembled fragments. Each fragment is represented by:
Size of common 5' substring|Name and size of DNA fragment| Size of common 5' substring
Linear:
frag20| 6
\/
/\
6|frag23| 6
\/
/\
6|frag14
Circular:
-|2577|61
| \/
| /\
| 61|5681|98
| \/
| /\
| 98|2389|557
| \/
| /\
| 557-
| |
--------------------------
Synonym of small_fig()
Bases: pydna.dsdna.Dseqrecord
This class holds information about a DNA fragment in an assembly. This class is instantiated by the Assembly class and is not meant to be instantiated directly.
Provides a class for downloading sequences from genbank.
Class to facilitate download from genbank.
Parameters: | users_email : string
proxy : string, optional
tool : string, optional
|
---|
Examples
>>> import pydna
>>> gb=pydna.Genbank("me@mail.se", proxy = "http://proxy.com:3128")
>>> rec = gb.nucleotide("L09137") # <- pUC19 from genbank
>>> print len(rec)
2686
Download a genbank nuclotide record.
Item is a string containing one genbank acession number [5] for a nucleotide file. Start and stop are intervals to be downloaded. This is useful as some genbank records are large. If strand is “c”, “C”, “crick”, “Crick”, “antisense”,”Antisense”, “2” or 2, the watson(antisense) strand is returned, otherwise the sense strand is returned.
Alternatively, item can be a string containing an url that returns a sequence in genbank or FASTA format.
The genbank E-utilities can provide such urls [6].
Result is returned as a Dseqrecord object.
Genbank nucleotide accession numbers have this format:
References
[5] | http://www.dsimb.inserm.fr/~fuchs/M2BI/AnalSeq/Annexes/Sequences/Accession_Numbers.htm |
[6] | http://www.ncbi.nlm.nih.gov/books/NBK25499/#chapter4.EFetch |
Provides two classes, Dseq and Dseqrecord, for handling double stranded DNA sequences. Dseq and Dseqrecord are subclasses of Biopythons Seq and SeqRecord classes, respectively. These classes support the notion of circular and linear DNA.
Bases: Bio.Seq.Seq
Dseq is a class designed to hold information for a double stranded DNA fragment. Dseq also holds information describing the topology of the DNA fragment (linear or circular).
Parameters: | watson : str
crick : str, optional
ovhg : int, optional
linear : bool, optional
circular : bool, optional
alphabet : Bio.Alphabet, optional
|
---|
See also
Examples
Dseq is a subclass of the Biopython Seq object. It stores two strings representing the watson (sense) and crick(antisense) strands. two properties called linear and circular, and a numeric value ovhg (overhang) describing the stagger for the watson and crick strand in the 5’ end of the fragment.
The most common usage is probably to create a Dseq object as a part of a Dseqrecord object (see pydna.dsdna.Dseqrecord).
There are three ways of creating a Dseq object directly:
Only one argument (string):
>>> import pydna
>>> pydna.Dseq("aaa")
Dseq(-3)
aaa
ttt
The given string will be interpreted as the watson strand of a blunt, linear double stranded sequence object. The crick strand is created automatically from the watson strand.
Two arguments (string, string):
>>> import pydna
>>> pydna.Dseq("gggaaat","ttt")
Dseq(-7)
gggaaat
ttt
If both watson and crick are given, but not ovhg an attempt will be made to find the best annealing between the strands. There are limitations to this! For long fragments it is quite slow. The length of the annealing sequences have to be at least half the length of the shortest of the strands.
Three arguments (string, string, ovhg=int):
The ovhg parameter controls the stagger at the five prime end:
ovhg=-2
XXXXX
XXXXX
ovhg=-1
XXXXX
XXXXX
ovhg=0
XXXXX
XXXXX
ovhg=1
XXXXX
XXXXX
ovhg=2
XXXXX
XXXXX
Example of creating Dseq objects with different amounts of stagger:
>>> pydna.Dseq(watson="agt",crick="actta",ovhg=-2)
Dseq(-7)
agt
attca
>>> pydna.Dseq(watson="agt",crick="actta",ovhg=-1)
Dseq(-6)
agt
attca
>>> pydna.Dseq(watson="agt",crick="actta",ovhg=0)
Dseq(-5)
agt
attca
>>> pydna.Dseq(watson="agt",crick="actta",ovhg=1)
Dseq(-5)
agt
attca
>>> pydna.Dseq(watson="agt",crick="actta",ovhg=2)
Dseq(-5)
agt
attca
If the ovhg parameter is psecified a crick strand also needs to be supplied, otherwise an exception is raised.
>>> pydna.Dseq(watson="agt",ovhg=2)
Traceback (most recent call last):
File "<stdin>", line 1, in <module>
File "/usr/local/lib/python2.7/dist-packages/pydna_/dsdna.py", line 169, in __init__
else:
Exception: ovhg defined without crick strand!
The default alphabet is set to Biopython IUPACAmbiguousDNA
The shape of the fragment is set by either:
Note that both ends of the DNA fragment has to be compatible to set circular = True (or linear = False).
>>> pydna.Dseq("aaa","ttt")
Dseq(-3)
aaa
ttt
>>> pydna.Dseq("aaa","ttt",ovhg=0)
Dseq(-3)
aaa
ttt
>>> pydna.Dseq("aaa", "ttt", linear = False ,ovhg=0)
Dseq(o3)
aaa
ttt
>>> pydna.Dseq("aaa", "ttt", circular = True , ovhg=0)
Dseq(o3)
aaa
ttt
Coercing to string
>>> a=pydna.Dseq("tttcccc","aaacccc")
>>> a
Dseq(-11)
tttcccc
ccccaaa
>>> str(a)
'ggggtttcccc'
The double stranded part is accessible through the dsdata property:
>>> a.dsdata
'ttt'
A Dseq object can be longer that either the watson or crick strands.
<-- length -->
GATCCTTT
AAAGCCTAG
The slicing of a linear Dseq object works mostly as it does for a string.
>>> s="ggatcc"
>>> s[2:3]
'a'
>>> s[2:4]
'at'
>>> s[2:4:-1]
''
>>> s[::2]
'gac'
>>> import pydna
>>> d=pydna.Dseq(s, linear=True)
>>> d[2:3]
Dseq(-1)
a
t
>>> d[2:4]
Dseq(-2)
at
ta
>>> d[2:4:-1]
Dseq(-0)
>>> d[::2]
Dseq(-3)
gac
ctg
The slicing of a circular Dseq object has a slightly different meaning.
>>> s="ggAtCc"
>>> d=pydna.Dseq(s, circular=True)
>>> d
Dseq(o6)
ggAtCc
ccTaGg
>>> d[4:3]
Dseq(-5)
CcggA
GgccT
The slice [X:X] produces an empty slice for a string, while this will return the linearized sequence starting at X:
>>> s="ggatcc"
>>> d=pydna.Dseq(s, circular=True)
>>> d
Dseq(o6)
ggatcc
cctagg
>>> d[3:3]
Dseq(-6)
tccgga
aggcct
>>>
Fill in of five prime protruding ends and chewing back of three prime protruding ends by a DNA polymerase providing both 5’-3’ DNA polymerase activity and 3’-5’ nuclease acitivty (such as T4 DNA polymerase). This can be done in presence of any combination of the four A, G, C or T. Default are all four nucleotides together.
Alias for the t4() method
Parameters: | nucleotides : str |
---|
Examples
>>> import pydna
>>> a=pydna.Dseq("gatcgatc")
>>> a
Dseq(-8)
gatcgatc
ctagctag
>>> a.T4()
Dseq(-8)
gatcgatc
ctagctag
>>> a.T4("t")
Dseq(-8)
gatcgat
tagctag
>>> a.T4("a")
Dseq(-8)
gatcga
agctag
>>> a.T4("g")
Dseq(-8)
gatcg
gctag
>>>
Returns a list of linear Dseq fragments produced in the digestion. If there are no cuts, an empty list is returned.
Parameters: | enzymes : enzyme object or iterable of such objects
|
---|---|
Returns: | frags : list
|
Examples
>>> from pydna import Dseq
>>> seq=Dseq("ggatccnnngaattc")
>>> seq
Dseq(-15)
ggatccnnngaattc
cctaggnnncttaag
>>> from Bio.Restriction import BamHI,EcoRI
>>> type(seq.cut(BamHI))
<type 'list'>
>>> for frag in seq.cut(BamHI):
... print frag.fig()
Dseq(-5)
g
cctag
Dseq(-14)
gatccnnngaattc
gnnncttaag
>>> seq.cut(EcoRI, BamHI) == seq.cut(BamHI, EcoRI)
True
>>> a,b,c = seq.cut(EcoRI, BamHI)
>>> a+b+c
Dseq(-15)
ggatccnnngaattc
cctaggnnncttaag
>>>
Returns a representation of the sequence, truncated if longer than 30 bp:
Examples
>>> import pydna
>>> a=pydna.Dseq("atcgcttactagcgtactgatcatctgact")
>>> a
Dseq(-30)
atcgcttactagcgtactgatcatctgact
tagcgaatgatcgcatgactagtagactga
>>> a+=Dseq("A")
>>> a
Dseq(-31)
atcg..actA
tagc..tgaT
Fill in of five prime protruding end with a DNA polymerase that has only DNA polymerase activity (such as exo-klenow [7]) and any combination of A, G, C or T. Default are all four nucleotides together.
Parameters: | nucleotides : str |
---|
References
[7] | http://en.wikipedia.org/wiki/Klenow_fragment#The_exo-_Klenow_fragment |
Examples
>>> import pydna
>>> a=pydna.Dseq("aaa", "ttt")
>>> a
Dseq(-3)
aaa
ttt
>>> a.fill_in()
Dseq(-3)
aaa
ttt
>>> b=pydna.Dseq("caaa", "cttt")
>>> b
Dseq(-5)
caaa
tttc
>>> b.fill_in()
Dseq(-5)
caaag
gtttc
>>> b.fill_in("g")
Dseq(-5)
caaag
gtttc
>>> b.fill_in("tac")
Dseq(-5)
caaa
tttc
>>> c=pydna.Dseq("aaac", "tttg")
>>> c
Dseq(-5)
aaac
gttt
>>> c.fill_in()
Dseq(-5)
aaac
gttt
>>>
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].
Returns -1 if the subsequence is NOT found.
Parameters: | sub : string or Seq object
start : int, optional
end : int, optional
|
---|
Examples
>>> import pydna
>>> seq = Dseq("atcgactgacgtgtt")
>>> seq
Dseq(-15)
atcgactgacgtgtt
tagctgactgcacaa
>>> seq.find("gac")
3
>>> seq = pydna.Dseq(watson="agt",crick="actta",ovhg=-2)
>>> seq
Dseq(-7)
agt
attca
>>> seq.find("taa")
2
Returns a tuple describing the structure of the 5’ end of the DNA fragment
See also
Examples
>>> import pydna
>>> a=pydna.Dseq("aaa", "ttt")
>>> a
Dseq(-3)
aaa
ttt
>>> a.five_prime_end()
('blunt', '')
>>> a=pydna.Dseq("aaa", "ttt", ovhg=1)
>>> a
Dseq(-4)
aaa
ttt
>>> a.five_prime_end()
("3'", 't')
>>> a=pydna.Dseq("aaa", "ttt", ovhg=-1)
>>> a
Dseq(-4)
aaa
ttt
>>> a.five_prime_end()
("5'", 'a')
>>>
Returns a circularized Dseq object. This can only be done if the two ends are compatible, otherwise a TypeError is raised.
Examples
>>> import pydna
>>> a=pydna.Dseq("catcgatc")
>>> a
Dseq(-8)
catcgatc
gtagctag
>>> a.looped()
Dseq(o8)
catcgatc
gtagctag
>>> a.T4("t")
Dseq(-8)
catcgat
tagctag
>>> a.T4("t").looped()
Dseq(o7)
catcgat
gtagcta
>>> a.T4("a")
Dseq(-8)
catcga
agctag
>>> a.T4("a").looped()
Traceback (most recent call last):
File "<stdin>", line 1, in <module>
File "/usr/local/lib/python2.7/dist-packages/pydna/dsdna.py", line 357, in looped
if type5 == type3 and str(sticky5) == str(rc(sticky3)):
TypeError: DNA cannot be circularized.
5' and 3' sticky ends not compatible!
>>>
Simulates treatment a nuclease with 5’-3’ and 3’-5’ single strand specific exonuclease activity (such as mung bean nuclease [8])
ggatcc -> gatcc
ctaggg ctagg
ggatcc -> ggatc
tcctag cctag
>>> import pydna
>>> b=pydna.Dseq("caaa", "cttt")
>>> b
Dseq(-5)
caaa
tttc
>>> b.mung()
Dseq(-3)
aaa
ttt
>>> c=pydna.Dseq("aaac", "tttg")
>>> c
Dseq(-5)
aaac
gttt
>>> c.mung()
Dseq(-3)
aaa
ttt
References
[8] | http://en.wikipedia.org/wiki/Mung_bean_nuclease |
Returns a Dseq object where watson and crick have switched places.
Examples
>>> import pydna
>>> a=pydna.Dseq("catcgatc")
>>> a
Dseq(-8)
catcgatc
gtagctag
>>> b=a.reverse_complement()
>>> b
Dseq(-8)
gatcgatg
ctagctac
>>>
Returns a tuple describing the structure of the 5’ end of the DNA fragment
>>> import pydna
>>> a=pydna.Dseq("aaa", "ttt")
>>> a
Dseq(-3)
aaa
ttt
>>> a.three_prime_end()
('blunt', '')
>>> a=pydna.Dseq("aaa", "ttt", ovhg=1)
>>> a
Dseq(-4)
aaa
ttt
>>> a.three_prime_end()
("3'", 'a')
>>> a=pydna.Dseq("aaa", "ttt", ovhg=-1)
>>> a
Dseq(-4)
aaa
ttt
>>> a.three_prime_end()
("5'", 't')
>>>
See also
Returns a blunt, linear copy of a circular Dseq object. This can only be done if the Dseq object is circular, otherwise a TypeError is raised.
Examples
>>> import pydna
>>> a=pydna.Dseq("catcgatc", circular=True)
>>> a
Dseq(o8)
catcgatc
gtagctag
>>> a.tolinear()
Dseq(-8)
catcgatc
gtagctag
>>>
Bases: Bio.SeqRecord.SeqRecord
Dseqrecord is a double stranded version of the Biopython SeqRecord [9] class. The Dseqrecord object holds a Dseq object describing the sequence. Additionally, Dseqrecord hold meta information about the sequence in the from of a list of SeqFeatures, in the same way as the SeqRecord does. The Dseqrecord can be initialized with a string, Seq, Dseq, SeqRecord or another Dseqrecord. The sequence information will be stored in a Dseq object in all cases. Dseqrecord objects can be read or parsed from sequences in Fasta, Embl or Genbank format.
There is a short representation associated with the Dseqrecord. Dseqrecord(-3) represents a linear sequence of length 2 while Dseqrecord(o7) represents a circular sequence of length 7.
Dseqrecord and Dseq share the same concept of length
<-- length -->
GATCCTTT
AAAGCCTAG
Parameters: | record : string, Seq, SeqRecord, Dseq or other Dseqrecord object
circular : bool, optional
linear : bool, optional
|
---|
References
[9] | http://biopython.org/wiki/SeqRecord |
Examples
>>> from pydna import Dseqrecord
>>> a=Dseqrecord("aaa")
>>> a
Dseqrecord(-3)
>>> a.seq
Dseq(-3)
aaa
ttt
>>> from Bio.Seq import Seq
>>> b=Dseqrecord(Seq("aaa"))
>>> b
Dseqrecord(-3)
>>> b.seq
Dseq(-3)
aaa
ttt
>>> from Bio.SeqRecord import SeqRecord
>>> c=Dseqrecord(SeqRecord(Seq("aaa")))
>>> c
Dseqrecord(-3)
>>> c.seq
Dseq(-3)
aaa
ttt
>>> a.seq.alphabet
IUPACAmbiguousDNA()
>>> b.seq.alphabet
IUPACAmbiguousDNA()
>>> c.seq.alphabet
IUPACAmbiguousDNA()
>>>
Returns the cSEGUID for the sequence. The cSEGUID is the SEGUID checksum calculated for the lexicographically minimal string rotation of the sequence or its reverse complement.
Only defined for circular sequences.
The cSEGUID checksum uniqely identifies a circular sequence.
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaat", circular=True)
>>> a.cseguid()
'oopV+6158nHJqedi8lsshIfcqYA'
>>> a=pydna.Dseqrecord("ataa", circular=True)
>>> a.cseguid()
'oopV+6158nHJqedi8lsshIfcqYA'
Digest the Dseqrecord object with one or more restriction enzymes. returns a list of linear Dseqrecords. If there are no cuts, an empty list is returned.
See also Dseq.cut()
Parameters: | enzymes : enzyme object or iterable of such objects
|
---|---|
Returns: | Dseqrecord_frags : list
|
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("ggatcc")
>>> from Bio.Restriction import BamHI
>>> a.cut(BamHI)
[Dseqrecord(-5), Dseqrecord(-5)]
>>> frag1, frag2 = a.cut(BamHI)
>>> frag1.seq
Dseq(-5)
g
cctag
>>> frag2.seq
Dseq(-5)
gatcc
g
Returns the sequence as a string using a format supported by Biopython SeqIO [10]. Default is “gb” which is short for Genbank.
References
[10] | http://biopython.org/wiki/SeqIO |
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaa")
>>> a.annotations['date'] = '02-FEB-2013'
>>> a
Dseqrecord(-3)
>>> print a.format()
LOCUS . 3 bp DNA linear UNK 02-FEB-2013
DEFINITION .
ACCESSION <unknown id>
VERSION <unknown id>
KEYWORDS .
SOURCE .
ORGANISM .
.
FEATURES Location/Qualifiers
ORIGIN
1 aaa
//
This method is similar to cut() but throws an exception if there is not excactly on cut i.e. none or more than one digestion products.
Returns a circular version of the Dseqrecord object. The underlying Dseq object has to have compatible ends.
See also
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaa")
>>> a
Dseqrecord(-3)
>>> b=a.looped()
>>> b
Dseqrecord(o3)
>>>
Returns the reverse complement.
See also
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("ggaatt")
>>> a
Dseqrecord(-6)
>>> a.seq
Dseq(-6)
ggaatt
ccttaa
>>> a.reverse_complement().seq
Dseq(-6)
aattcc
ttaagg
>>>
Returns the SEGUID [11] for the sequence.
References
[11] | http://wiki.christophchamp.com/index.php/SEGUID |
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaa")
>>> a.seguid()
'YG7G6b2Kj/KtFOX63j8mRHHoIlE'
Returns a circular Dseqrecord with a new origin <shift>. This only works on circular Dseqrecords. If we consider the following circular sequence:
The T and the G on the watson strand are linked together as well as the A and the C of the of the crick strand.
if shift is 1, this indicates a new origin at position 1:
new sequence:
Shift is always positive and 0<shift<length, so in the example below, permissible values of shift are 1,2 and 3
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaat",circular=True)
>>> a
Dseqrecord(o4)
>>> a.seq
Dseq(o4)
aaat
ttta
>>> b=a.shifted(1)
>>> b
Dseqrecord(o4)
>>> b.seq
Dseq(o4)
aata
ttat
Adds a seguid stamp to the description property. This will show in the genbank format. The following string:
SEGUID <seguid>
will be appended to the description property of the Dseqrecord object (string).
The stamp can be verified with verify_stamp()
See also
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaa")
>>> a.stamp()
>>> a.description
'<unknown description> SEGUID YG7G6b2Kj/KtFOX63j8mRHHoIlE'
>>> a.verify_stamp()
True
This function returns a new circular sequence (Dseqrecord object), which has ben rotated in such a way that there is maximum overlap between the sequence and ref, which may be a string, Biopython Seq, SeqRecord object or another Dseqrecord object.
The reason for using this could be to rotate a recombinant plasmid so that it starts at the same position as the vector backbone after cloning.
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("gaat",circular=True)
>>> a.seq
Dseq(o4)
gaat
ctta
>>> d = a[2:] + a[:2]
>>> d.seq
Dseq(-4)
atga
tact
>>> insert=pydna.Dseqrecord("CCC")
>>> recombinant = (d+insert).looped()
>>> recombinant.seq
Dseq(o7)
atgaCCC
tactGGG
>>> recombinant.synced(a).seq
Dseq(o7)
gaCCCat
ctGGGta
Returns a linear, blunt copy of a circular Dseqrecord object. The underlying Dseq object has to be circular.
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("aaa", circular = True)
>>> a
Dseqrecord(o3)
>>> b=a.tolinear()
>>> b
Dseqrecord(-3)
>>>
Verifies the SEGUID stamp in the description property is valid. True if stamp match the sequid calculated from the sequence. Exception raised if no stamp can be found.
>>> import pydna
>>> a=pydna.Dseqrecord("aaa")
>>> a.annotations['date'] = '02-FEB-2013'
>>> a.seguid()
'YG7G6b2Kj/KtFOX63j8mRHHoIlE'
>>> print a.format("gb")
LOCUS . 3 bp DNA linear UNK 02-FEB-2013
DEFINITION .
ACCESSION <unknown id>
VERSION <unknown id>
KEYWORDS .
SOURCE .
ORGANISM .
.
FEATURES Location/Qualifiers
ORIGIN
1 aaa
//
>>> a.stamp()
>>> a
Dseqrecord(-3)
>>> print a.format("gb")
LOCUS . 3 bp DNA linear UNK 02-FEB-2013
DEFINITION <unknown description> SEGUID YG7G6b2Kj/KtFOX63j8mRHHoIlE
ACCESSION <unknown id>
VERSION <unknown id>
KEYWORDS .
SOURCE .
ORGANISM .
.
FEATURES Location/Qualifiers
ORIGIN
1 aaa
//
>>> a.verify_stamp()
True
>>>
See also
Writes the Dseqrecord to a file using the format f, which must be a format supported by Biopython SeqIO for writing [12]. Default is “gb” which is short for Genbank. Note that Biopython SeqIO reads more formats than it writes.
Filename is the path to the file where the sequece is to be written. The filename is optional, if it is not given, the description property (string) is used together with the format.
If obj is the Dseqrecord object, the default file name will be:
<obj.description>.<f>
Where <f> is “gb” by default. If the filename already exists and AND the sequence it contains is different, a new file name will be used so that the old file is not lost:
<obj.description>_NEW.<f>
References
[12] | http://biopython.org/wiki/SeqIO |
This function returns all DNA sequences found in data. If no sequences are found, an empty list is returned. This is a greedy function, use carefully.
Parameters: | data : string or iterable
ds : bool
|
---|---|
Returns: | list :
|
See also
References
[13] | http://biopython.org/wiki/SeqRecord |
returns the reverse complement of sequence (string) accepts mixed DNA/RNA
This function is similar the parse() function but expects one and only one sequence or and exception is thrown.
Parameters: | data : string
ds : bool
|
---|---|
Returns: | Dseqrecord :
|
See also
Notes
The data parameter is similar to the data parameter for parse().
This module provides a class for opening a sequence using an editor that accepts a file as a command line argument.
ApE - A plasmid Editor [14] is and excellent editor for this purpose.
[14] | http://biologylabs.utah.edu/jorgensen/wayned/ape/ |
The Editor class needs to be instantiated before use.
Parameters: | shell_command_for_editor : str
tmpdir : str, optional
|
---|
Examples
>>> import pydna
>>> ape = pydna.Editor("tclsh8.6 /home/bjorn/.ApE/apeextractor/ApE.vfs/lib/app-AppMain/AppMain.tcl")
>>> ape.open("aaa") # This command opens the sequence in the ApE editor
Finds all the common substrings between stringx and stringy longer than limit. This function is case sensitive.
returns a list of tuples describing the substrings The list is sorted longest -> shortest.
Examples
[(startx1,starty1,length1),(startx2,starty2,length2), ...]
startx1 = position in x where substring 1 starts starty1 = position in y where substring 1 starts length = lenght of substring
This module contain experimental functionality for smulation of agarose gel electrophoresis of DNA.
This code is in a very early stage of development and documentation. Read the source code to get an idea for how it works.
A Gel contains lanes into which DNA can be loaded and run. The Gel is exposured to see the location and intensity of DNA.
This module contain functions for primer design.
This function return primer pairs that are useful for fusion of DNA sequences given in template. Given two sequences that we wish to fuse (a and b) to form fragment c.
_________ a _________ __________ b ________
/ \ / \
agcctatcatcttggtctctgca <--> TTTATATCGCATGACTCTTCTTT
||||||||||||||||||||||| |||||||||||||||||||||||
tcggatagtagaaccagagacgt <--> AAATATAGCGTACTGAGAAGAAA
agcctatcatcttggtctctgcaTTTATATCGCATGACTCTTCTTT
||||||||||||||||||||||||||||||||||||||||||||||
tcggatagtagaaccagagacgtAAATATAGCGTACTGAGAAGAAA
\___________________ c ______________________/
We can design tailed primers to fuse a and b by fusion PCR, Gibson assembly or in-vivo homologous recombination. The basic requirements for the primers for the three techniques are the same.
Design tailed primers incorporating a part of the next or previous fragment to be assembled.
agcctatcatcttggtctctgca
|||||||||||||||||||||||
gagacgtAAATATA
|||||||||||||||||||||||
tcggatagtagaaccagagacgt
TTTATATCGCATGACTCTTCTTT
|||||||||||||||||||||||
ctctgcaTTTATAT
|||||||||||||||||||||||
AAATATAGCGTACTGAGAAGAAA
PCR products with flanking sequences are formed in the PCR process.
agcctatcatcttggtctctgcaTTTATAT
||||||||||||||||||||||||||||||
tcggatagtagaaccagagacgtAAATATA
\____________/
identical
sequences
____________
/ \
ctctgcaTTTATATCGCATGACTCTTCTTT
||||||||||||||||||||||||||||||
gagacgtAAATATAGCGTACTGAGAAGAAA
The fragments can be fused by any of the techniques mentioned earlier.
agcctatcatcttggtctctgcaTTTATATCGCATGACTCTTCTTT
||||||||||||||||||||||||||||||||||||||||||||||
tcggatagtagaaccagagacgtAAATATAGCGTACTGAGAAGAAA
Parameters: | templates : list of Dseqrecord
minlength : int, optional
maxlength : int, optional
tot_length : int, optional
target_tm : float, optional
primerc : float, optional
saltc : float, optional
formula : function
|
---|---|
Returns: | primer_pairs : list of tuples of Bio.Seqrecord objects
|
Examples
>>> import pydna
>>> a=pydna.Dseqrecord("atgactgctaacccttccttggtgttgaacaagatcgacgacatttcgttcgaaacttacgatg")
>>> b=pydna.Dseqrecord("ccaaacccaccaggtaccttatgtaagtacttcaagtcgccagaagacttcttggtcaagttgcc")
>>> c=pydna.Dseqrecord("tgtactggtgctgaaccttgtatcaagttgggtgttgacgccattgccccaggtggtcgtttcgtt")
>>> primer_pairs = pydna.assembly_primers([a,b,c])
>>> p=[]
>>> for t, (f,r) in zip([a,b,c], primer_pairs): p.append(pydna.pcr(f,r,t))
>>> p
[Amplicon(87), Amplicon(110), Amplicon(90)]
>>> assemblyobj = pydna.Assembly(p)
>>> assemblyobj
Assembly object:
Sequences........................: 87, 110, 90
Sequences with shared homologies.: No analyzed sequences
Homology limit (bp)..............: 25
Number of overlaps...............: No overlaps
Nodes in graph...................: No graph
Assembly protocol................: No protocol
Circular products................: No circular products
Linear products..................: No linear products
>>>
>>> assemblyobj.assemble_gibson_linear()
'5 linear products were formed with the following sizes (bp): 195, 155, 150, 47, 45'
>>> assemblyobj.linear_products
[Dseqrecord(-195), Dseqrecord(-155), Dseqrecord(-150), Dseqrecord(-47), Dseqrecord(-45)]
>>> assemblyobj.linear_products[0].seguid()
'1eNv3d/1PqDPP8qJZIVoA45Www8'
>>> (a+b+c).seguid()
'1eNv3d/1PqDPP8qJZIVoA45Www8'
>>> pydna.eq(a+b+c, assemblyobj.linear_products[0])
True
>>>
>>> print assemblyobj.linear_products[0].small_fig()
87bp_PCR_prod|47
\/
/\
47|110bp_PCR_prod|45
\/
/\
45|90bp_PCR_prod
>>>
This function can design primers for PCR amplification of a given sequence. This function accepts a Dseqrecord object containing the template sequence and returns a tuple cntaining two ::mod`Bio.SeqRecord.SeqRecord` objects describing the primers.
Primer tails can optionally be given in the form of strings.
An predesigned primer can be given, either the forward or reverse primers. In this case this function tries to design a primer with a Tm to match the given primer.
Parameters: | template : Dseqrecord
minlength : int, optional
maxlength : int, optional
fp, rp : SeqRecord, optional
fp_tail, rp_tail : string, optional
target_tm : float, optional
primerc : float, optional
saltc : float, optional
formula : function
|
---|---|
Returns: | fp, rp : tuple
|
Examples
>>> import pydna
>>> t=pydna.Dseqrecord("atgactgctaacccttccttggtgttgaacaagatcgacgacatttcgttcgaaacttacgatg")
>>> t
Dseqrecord(-64)
>>> pf,pr = pydna.cloning_primers(t)
>>> pf
SeqRecord(seq=Seq('atgactgctaacccttc', Alphabet()), id='pfw64', name='pfw64', description='pfw64', dbxrefs=[])
>>> pr
SeqRecord(seq=Seq('catcgtaagtttcgaac', Alphabet()), id='prv64', name='prv64', description='prv64', dbxrefs=[])
>>> pcr_prod = pydna.pcr(pf, pr, t)
>>> pcr_prod
Amplicon(64)
>>>
>>> print pcr_prod.figure()
5atgactgctaacccttc...gttcgaaacttacgatg3
||||||||||||||||| tm 42.4 (dbd) 52.9
3caagctttgaatgctac5
5atgactgctaacccttc3
||||||||||||||||| tm 44.5 (dbd) 54.0
3tactgacgattgggaag...caagctttgaatgctac5
>>> pf,pr = pydna.cloning_primers(t, fp_tail="GGATCC", rp_tail="GAATTC")
>>> pf
SeqRecord(seq=Seq('GGATCCatgactgctaacccttc', Alphabet()), id='pfw64', name='pfw64', description='pfw64', dbxrefs=[])
>>> pr
SeqRecord(seq=Seq('GAATTCcatcgtaagtttcgaac', Alphabet()), id='prv64', name='prv64', description='prv64', dbxrefs=[])
>>> pcr_prod = pydna.pcr(pf, pr, t)
>>> print pcr_prod.figure()
5atgactgctaacccttc...gttcgaaacttacgatg3
||||||||||||||||| tm 42.4 (dbd) 52.9
3caagctttgaatgctacCTTAAG5
5GGATCCatgactgctaacccttc3
||||||||||||||||| tm 44.5 (dbd) 54.0
3tactgacgattgggaag...caagctttgaatgctac5
>>> print pcr_prod.seq
GGATCCatgactgctaacccttccttggtgttgaacaagatcgacgacatttcgttcgaaacttacgatgGAATTC
>>>
>>> from Bio.Seq import Seq
>>> from Bio.SeqRecord import SeqRecord
>>> pf = SeqRecord(Seq("atgactgctaacccttccttggtgttg"))
>>> pf,pr = pydna.cloning_primers(t, fp = pf, fp_tail="GGATCC", rp_tail="GAATTC")
>>> pf
SeqRecord(seq=Seq('GGATCCatgactgctaacccttccttggtgttg', Alphabet()), id='pfw64', name='pfw64', description='pfw64', dbxrefs=[])
>>> pr
SeqRecord(seq=Seq('GAATTCcatcgtaagtttcgaacgaaatgtcgtc', Alphabet()), id='prv64', name='prv64', description='prv64', dbxrefs=[])
>>> ampl = pydna.pcr(pf,pr,t)
>>> print ampl.figure()
5atgactgctaacccttccttggtgttg...gacgacatttcgttcgaaacttacgatg3
|||||||||||||||||||||||||||| tm 57.5 (dbd) 72.2
3ctgctgtaaagcaagctttgaatgctacCTTAAG5
5GGATCCatgactgctaacccttccttggtgttg3
||||||||||||||||||||||||||| tm 59.0 (dbd) 72.3
3tactgacgattgggaaggaaccacaac...ctgctgtaaagcaagctttgaatgctac5
>>>
This module provides miscellaneous functions.
Decompose s into Lyndon words according to the Chen-Fox-Lyndon theorem. The arguments are the same as for ChenFoxLyndonBreakpoints but the return values are subsequences of s rather than indices of breakpoints.
Algorithms on strings and sequences based on Lyndon words. David Eppstein, October 2011.
Find starting positions of Chen-Fox-Lyndon decomposition of s. The decomposition is a set of Lyndon words that start at 0 and continue until the next position. 0 itself is not output, but the final breakpoint at the end of s is. The argument s must be of a type that can be indexed (e.g. a list, tuple, or string). The algorithm follows Duval, J. Algorithms 1983, but uses 0-based indexing rather than Duval’s choice of 1-based indexing.
Algorithms on strings and sequences based on Lyndon words. David Eppstein, October 2011.
Find the rotation of s that is smallest in lexicographic order. Duval 1983 describes how to modify his algorithm to do so but I think it’s cleaner and more general to work from the ChenFoxLyndon output.
Algorithms on strings and sequences based on Lyndon words. David Eppstein, October 2011.
This function tries to copy all features in source_seq and copy them to target_seq. Source_sr and target_sr are objects with a features property, such as Dseqrecord or Biopython SeqRecord.
Parameters: | source_seq : SeqRecord or Dseqrecord
target_seq : SeqRecord or Dseqrecord
|
---|---|
Returns: | bool : True
|
Returns the cSEGUID for the sequence. The cSEGUID is the SEGUID checksum calculated for the lexicographically minimal string rotation of a DNA sequence. Only defined for circular sequences.
Compares two or more DNA sequences for equality i.e. they represent the same DNA molecule. Comparisons are case insensitive.
Parameters: | args : iterable
circular : bool, optional
linear : bool, optional
|
---|---|
Returns: | eq : bool
|
Notes
Compares two or more DNA sequences for equality i.e. if they represent the same DNA molecule.
Two linear sequences are considiered equal if either:
Two circular sequences are considered equal if they are circular permutations:
The topology for the comparison can be set using one of the keywords linear or circular to True or False.
If circular or linear is not set, it will be deduced from the topology of each sequence for sequences that have a linear or circular attribute (like Dseq and Dseqrecord).
Examples
>>> from pydna import eq, Dseqrecord
>>> eq("aaa","AAA")
True
>>> eq("aaa","AAA","TTT")
True
>>> eq("aaa","AAA","TTT","tTt")
True
>>> eq("aaa","AAA","TTT","tTt", linear=True)
True
>>> eq("Taaa","aTaa", linear = True)
False
>>> eq("Taaa","aTaa", circular = True)
True
>>> a=Dseqrecord("Taaa")
>>> b=Dseqrecord("aTaa")
>>> eq(a,b)
False
>>> eq(a,b,circular=True)
True
>>> a=a.looped()
>>> b=b.looped()
>>> eq(a,b)
True
>>> eq(a,b,circular=False)
False
>>> eq(a,b,linear=True)
False
>>> eq(a,b,linear=False)
True
>>> eq("ggatcc","GGATCC")
True
>>> eq("ggatcca","GGATCCa")
True
>>> eq("ggatcca","tGGATCC")
True
Shift the origin of seq which is assumed to be a circular sequence.
Parameters: | seq : string, Biopython Seq, Biopython SeqRecord, Dseq or Dseqrecord
|
---|---|
Returns: | new_seq : string, Biopython Seq, Biopython SeqRecord, Dseq or Dseqrecord
|
See also
Examples
>>> import pydna
>>> pydna.shift_origin("taaa",1)
'aaat'
>>> pydna.shift_origin("taaa",0)
'taaa'
>>> pydna.shift_origin("taaa",2)
'aata'
>>> pydna.shift_origin("gatc",2)
'tcga'
Synchronize two circular sequences. This function tries to rotate the circular sequence seq so that it has a maximum overlap with ref.
Parameters: | seq : string, Biopython Seq, Biopython SeqRecord, Dseq or Dseqrecord
ref : string, Biopython Seq, Biopython SeqRecord, Dseq or Dseqrecord
|
---|---|
Returns: | sequence : string, Biopython Seq, Biopython SeqRecord, Dseq or Dseqrecord
|
See also
Examples
>>> import pydna
>>> pydna.sync("taaatc","aaataa")
'aaatct'
>>> pydna.sync("taaatc","aaataa")
'aaatct'
>>> pydna.sync("taaat","aaataa")
'aaatt'