1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
|
# Copyright 2009 by Peter Cock. All rights reserved.
# This code is part of the Biopython distribution and governed by its
# license. Please see the LICENSE file that should have been included
# as part of this package.
"""SeqFeature related tests for SeqRecord objects from Bio.SeqIO.
Initially this takes matched tests of GenBank and FASTA files from the NCBI
and confirms they are consistent using our different parsers.
"""
import unittest
from Bio.Alphabet import generic_dna, generic_rna, generic_protein
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord
from Bio.SeqFeature import SeqFeature, FeatureLocation, ExactPosition
from Bio.SeqFeature import WithinPosition, BeforePosition, AfterPosition, OneOfPosition
class SeqRecordCreation(unittest.TestCase):
"""Test basic creation of SeqRecords."""
def test_annotations(self):
"""Pass in annotations to SeqRecords"""
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test")
self.assertEqual(rec.annotations, {})
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test",
annotations={"test" : ["a test"]})
self.assertEqual(rec.annotations["test"], ["a test"])
def test_letter_annotations(self):
"""Pass in letter annotations to SeqRecords"""
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test")
self.assertEqual(rec.annotations, {})
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test",
letter_annotations={"test" : [1, 2, 3, 4]})
self.assertEqual(rec.letter_annotations["test"], [1, 2, 3, 4])
#Now try modifying it to a bad value...
try:
rec.letter_annotations["bad"] = "abc"
self.assert_(False, "Adding a bad letter_annotation should fail!")
except (TypeError, ValueError), e:
pass
#Now try setting it afterwards to a bad value...
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test")
try:
rec.letter_annotations={"test" : [1, 2, 3]}
self.assert_(False, "Changing to bad letter_annotations should fail!")
except (TypeError, ValueError), e:
pass
#Now try setting it at creation time to a bad value...
try:
rec = SeqRecord(Seq("ACGT", generic_dna),
id="Test", name="Test", description="Test",
letter_annotations={"test" : [1, 2, 3]})
self.assert_(False, "Wrong length letter_annotations should fail!")
except (TypeError, ValueError), e:
pass
class SeqRecordMethods(unittest.TestCase):
"""Test SeqRecord methods."""
def setUp(self) :
f0 = SeqFeature(FeatureLocation(0,26), type="source",
qualifiers={"mol_type":["fake protein"]})
f1 = SeqFeature(FeatureLocation(0,ExactPosition(10)))
f2 = SeqFeature(FeatureLocation(WithinPosition(12,3),BeforePosition(22)))
f3 = SeqFeature(FeatureLocation(AfterPosition(16),
OneOfPosition([ExactPosition(25),AfterPosition(26)])))
self.record = SeqRecord(Seq("ABCDEFGHIJKLMNOPQRSTUVWZYX", generic_protein),
id="TestID", name="TestName", description="TestDescr",
dbxrefs=["TestXRef"], annotations={"k":"v"},
letter_annotations = {"fake":"X"*26},
features = [f0,f1,f2,f3])
def test_slice_variantes(self):
"""Simple slices using different start/end values"""
for start in range(-30,30)+[None] :
for end in range(-30,30)+[None] :
if start is None and end is None : continue
rec = self.record[start:end]
seq = self.record.seq[start:end]
seq_str = str(self.record.seq)[start:end]
self.assertEqual(seq_str, str(seq))
self.assertEqual(seq_str, str(rec.seq))
self.assertEqual("X"*len(seq_str), rec.letter_annotations["fake"])
def test_slice_simple(self):
"""Simple slice"""
rec = self.record
self.assertEqual(len(rec), 26)
left = rec[:10]
self.assertEqual(str(left.seq), str(rec.seq[:10]))
right = rec[-10:]
self.assertEqual(str(right.seq), str(rec.seq[-10:]))
mid = rec[12:22]
self.assertEqual(str(mid.seq), str(rec.seq[12:22]))
for sub in [left, right, mid] :
self.assertEqual(len(sub), 10)
self.assertEqual(sub.id, "TestID")
self.assertEqual(sub.name, "TestName")
self.assertEqual(sub.description, "TestDescr")
self.assertEqual(sub.letter_annotations, {"fake":"X"*10})
self.assertEqual(sub.dbxrefs, []) # May change this...
self.assertEqual(sub.annotations, {}) # May change this...
self.assertEqual(len(sub.features), 1)
#By construction, each feature matches the full sliced region:
self.assertEqual(str(sub.features[0].extract(sub.seq)), str(sub.seq))
self.assertEqual(sub.features[0].extract(str(sub.seq)), str(sub.seq))
def test_add_simple(self):
"""Simple addition"""
rec = self.record + self.record
self.assertEqual(len(rec), 52)
self.assertEqual(rec.id, "TestID")
self.assertEqual(rec.name, "TestName")
self.assertEqual(rec.description, "TestDescr")
self.assertEqual(rec.dbxrefs, ["TestXRef"])
self.assertEqual(rec.annotations, {"k":"v"})
self.assertEqual(rec.letter_annotations, {"fake":"X"*52})
self.assertEqual(len(rec.features), 2*len(self.record.features))
def test_add_seq(self):
"""Simple addition of Seq or string"""
for other in [Seq("BIO"), "BIO"] :
rec = self.record + other # will use SeqRecord's __add__ method
self.assertEqual(len(rec), 26+3)
self.assertEqual(str(rec.seq), str(self.record.seq)+"BIO")
self.assertEqual(rec.id, "TestID")
self.assertEqual(rec.name, "TestName")
self.assertEqual(rec.description, "TestDescr")
self.assertEqual(rec.dbxrefs, ["TestXRef"])
self.assertEqual(rec.annotations, {"k":"v"})
self.assertEqual(rec.letter_annotations, {})
self.assertEqual(len(rec.features), len(self.record.features))
self.assertEqual(rec.features[0].type, "source")
self.assertEqual(rec.features[0].location.nofuzzy_start, 0)
self.assertEqual(rec.features[0].location.nofuzzy_end, 26) #not +3
def test_add_seq_left(self):
"""Simple left addition of Seq or string"""
for other in [Seq("BIO"), "BIO"] :
rec = other + self.record # will use SeqRecord's __radd__ method
self.assertEqual(len(rec), 26+3)
self.assertEqual(str(rec.seq), "BIO"+str(self.record.seq))
self.assertEqual(rec.id, "TestID")
self.assertEqual(rec.name, "TestName")
self.assertEqual(rec.description, "TestDescr")
self.assertEqual(rec.dbxrefs, ["TestXRef"])
self.assertEqual(rec.annotations, {"k":"v"})
self.assertEqual(rec.letter_annotations, {})
self.assertEqual(len(rec.features), len(self.record.features))
self.assertEqual(rec.features[0].type, "source")
self.assertEqual(rec.features[0].location.nofuzzy_start, 3)
self.assertEqual(rec.features[0].location.nofuzzy_end, 26+3)
def test_slice_add_simple(self):
"""Simple slice and add"""
for cut in range(27) :
rec = self.record[:cut] + self.record[cut:]
self.assertEqual(str(rec.seq), str(self.record.seq))
self.assertEqual(len(rec), 26)
self.assertEqual(rec.id, "TestID")
self.assertEqual(rec.name, "TestName")
self.assertEqual(rec.description, "TestDescr")
self.assertEqual(rec.dbxrefs, []) # May change this...
self.assertEqual(rec.annotations, {}) # May change this...
self.assertEqual(rec.letter_annotations, {"fake":"X"*26})
self.assert_(len(rec.features) <= len(self.record.features))
def test_slice_add_shift(self):
"""Simple slice and add to shift"""
for cut in range(27) :
rec = self.record[cut:] + self.record[:cut]
self.assertEqual(str(rec.seq), str(self.record.seq[cut:] + self.record.seq[:cut]))
self.assertEqual(len(rec), 26)
self.assertEqual(rec.id, "TestID")
self.assertEqual(rec.name, "TestName")
self.assertEqual(rec.description, "TestDescr")
self.assertEqual(rec.dbxrefs, []) # May change this...
self.assertEqual(rec.annotations, {}) # May change this...
self.assertEqual(rec.letter_annotations, {"fake":"X"*26})
self.assert_(len(rec.features) <= len(self.record.features))
if __name__ == "__main__":
runner = unittest.TextTestRunner(verbosity = 2)
unittest.main(testRunner=runner)
|