Package Bio :: Package AlignIO :: Module PhylipIO
[hide private]
[frames] | no frames]

Source Code for Module Bio.AlignIO.PhylipIO

  1  # Copyright 2006-2008 by Peter Cock.  All rights reserved. 
  2  # This code is part of the Biopython distribution and governed by its 
  3  # license.  Please see the LICENSE file that should have been included 
  4  # as part of this package. 
  5   
  6  """ 
  7  Note 
  8  ==== 
  9  In TREE_PUZZLE (Schmidt et al. 2003) and PHYML (Guindon and Gascuel 2003) 
 10  a dot/period (".") in a sequence is interpreted as meaning the same 
 11  character as in the first sequence.  The PHYLIP 3.6 documentation says: 
 12   
 13     "a period was also previously allowed but it is no longer allowed, 
 14     because it sometimes is used in different senses in other programs" 
 15   
 16  At the time of writing, we do nothing special with a dot/period. 
 17  """ 
 18   
 19  from Bio.Alphabet import single_letter_alphabet 
 20  from sets import Set 
 21   
 22  from Bio.Align.Generic import Alignment 
 23  from Interfaces import AlignmentIterator, SequentialAlignmentWriter 
 24   
25 -class PhylipWriter(SequentialAlignmentWriter) :
26 """Phylip alignment writer."""
27 - def write_alignment(self, alignment) :
28 """Use this to write (another) single alignment to an open file. 29 30 This code will write interlaced alignments (when the sequences are 31 longer than 50 characters). 32 33 Note that record identifiers are strictly truncated at 10 characters. 34 35 For more information on the file format, please see: 36 http://evolution.genetics.washington.edu/phylip/doc/sequence.html 37 http://evolution.genetics.washington.edu/phylip/doc/main.html#inputfiles 38 """ 39 truncate=10 40 records = alignment.get_all_seqs() 41 handle = self.handle 42 43 if len(records)==0 : 44 raise ValueError("Must have at least one sequence") 45 length_of_seqs = alignment.get_alignment_length() 46 for record in records : 47 if length_of_seqs <> len(record.seq) : 48 raise ValueError("Sequences must all be the same length") 49 if length_of_seqs <= 0 : 50 raise ValueError("Non-empty sequences are required") 51 52 if len(records) > len(Set([r.id[:truncate] for r in records])) : 53 raise ValueError("Repeated identifier, possibly due to truncation") 54 55 56 # From experimentation, the use of tabs is not understood by the 57 # EMBOSS suite. The nature of the expected white space is not 58 # defined in the PHYLIP documentation, simply "These are in free 59 # format, separated by blanks". We'll use spaces to keep EMBOSS 60 # happy. 61 handle.write(" %i %s\n" % (len(records), length_of_seqs)) 62 block=0 63 while True : 64 for record in records : 65 if block==0 : 66 #Write name (truncated/padded to 10 characters) 67 """ 68 Quoting the PHYLIP version 3.6 documentation: 69 70 The name should be ten characters in length, filled out to 71 the full ten characters by blanks if shorter. Any printable 72 ASCII/ISO character is allowed in the name, except for 73 parentheses ("(" and ")"), square brackets ("[" and "]"), 74 colon (":"), semicolon (";") and comma (","). If you forget 75 to extend the names to ten characters in length by blanks, 76 the program [i.e. PHYLIP] will get out of synchronization 77 with the contents of the data file, and an error message will 78 result. 79 80 Note that Tab characters count as only one character in the 81 species names. Their inclusion can cause trouble. 82 """ 83 name = record.id.strip() 84 #Either remove the banned characters, or map them to something 85 #else like an underscore "_" or pipe "|" character... 86 for char in "[]()," : 87 name = name.replace(char,"") 88 for char in ":;" : 89 name = name.replace(char,"|") 90 91 #Now truncate and right pad to expected length. 92 handle.write(name[:truncate].ljust(truncate)) 93 else : 94 #write 10 space indent 95 handle.write(" "*truncate) 96 #Write five chunks of ten letters per line... 97 for chunk in range(0,5) : 98 i = block*50 + chunk*10 99 seq_segment = record.seq.tostring()[i:i+10] 100 #TODO - Force any gaps to be '-' character? Look at the alphabet... 101 #TODO - How to cope with '?' or '.' in the sequence? 102 handle.write(" %s" % seq_segment) 103 if i+10 > length_of_seqs : break 104 handle.write("\n") 105 block=block+1 106 if block*50 > length_of_seqs : break 107 handle.write("\n")
108
109 -class PhylipIterator(AlignmentIterator) :
110 """Reads a Phylip alignment file returning an Alignment object iterator. 111 112 Record identifiers are limited to at most 10 characters. 113 114 It only copes with interlaced phylip files! Sequential files won't work 115 where the sequences are split over multiple lines. 116 117 For more information on the file format, please see: 118 http://evolution.genetics.washington.edu/phylip/doc/sequence.html 119 http://evolution.genetics.washington.edu/phylip/doc/main.html#inputfiles 120 """ 121
122 - def _is_header(self, line) :
123 line = line.strip() 124 parts = filter(None, line.split()) 125 if len(parts)<>2 : 126 return False # First line should have two integers 127 try : 128 number_of_seqs = int(parts[0]) 129 length_of_seqs = int(parts[1]) 130 return True 131 except ValueError: 132 return False # First line should have two integers
133
134 - def next(self) :
135 handle = self.handle 136 137 try : 138 #Header we saved from when we were parsing 139 #the previous alignment. 140 line = self._header 141 del self._header 142 except AttributeError : 143 line = handle.readline() 144 145 if not line: return 146 line = line.strip() 147 parts = filter(None, line.split()) 148 if len(parts)<>2 : 149 raise ValueError("First line should have two integers") 150 try : 151 number_of_seqs = int(parts[0]) 152 length_of_seqs = int(parts[1]) 153 except ValueError: 154 raise ValueError("First line should have two integers") 155 156 assert self._is_header(line) 157 158 if self.records_per_alignment is not None \ 159 and self.records_per_alignment <> number_of_seqs : 160 raise ValueError("Found %i records in this alignment, told to expect %i" \ 161 % (number_of_seqs, self.records_per_alignment)) 162 163 ids = [] 164 seqs = [] 165 166 #Expects STRICT truncation/padding to 10 characters 167 #Does not require any white space between name and seq. 168 for i in range(0,number_of_seqs) : 169 line = handle.readline().rstrip() 170 ids.append(line[:10].strip()) #first ten characters 171 seqs.append([line[10:].strip().replace(" ","")]) 172 173 #Look for further blocks 174 line="" 175 while True : 176 #Skip any blank lines between blocks... 177 while ""==line.strip(): 178 line = handle.readline() 179 if not line : break #end of file 180 if not line : break #end of file 181 182 if self._is_header(line) : 183 #Looks like the start of a concatenated alignment 184 self._header = line 185 break 186 187 #print "New block..." 188 for i in range(0,number_of_seqs) : 189 seqs[i].append(line.strip().replace(" ","")) 190 line = handle.readline() 191 if (not line) and i+1 < number_of_seqs : 192 raise ValueError("End of file mid-block") 193 if not line : break #end of file 194 195 alignment = Alignment(self.alphabet) 196 for i in range(0,number_of_seqs) : 197 seq = "".join(seqs[i]) 198 if len(seq)<>length_of_seqs : 199 raise ValueError("Sequence %i length %i, expected length %i" \ 200 % (i+1, len(seq), length_of_seqs)) 201 alignment.add_sequence(ids[i], seq) 202 203 record = alignment.get_all_seqs()[-1] 204 assert ids[i] == record.id or ids[i] == record.description 205 record.id = ids[i] 206 record.name = ids[i] 207 record.description = ids[i] 208 return alignment
209 210 if __name__=="__main__" : 211 print "Running short mini-test" 212 213 phylip_text=""" 8 286 214 V_Harveyi_ --MKNWIKVA VAAIA--LSA A--------- ---------T VQAATEVKVG 215 B_subtilis MKMKKWTVLV VAALLAVLSA CG-------- ----NGNSSS KEDDNVLHVG 216 B_subtilis MKKALLALFM VVSIAALAAC GAGNDNQSKD NAKDGDLWAS IKKKGVLTVG 217 YA80_HAEIN MKKLLFTTAL LTGAIAFSTF ---------- -SHAGEIADR VEKTKTLLVG 218 FLIY_ECOLI MKLAHLGRQA LMGVMAVALV AG---MSVKS FADEG-LLNK VKERGTLLVG 219 E_coli_Gln --MKSVLKVS LAALTLAFAV S--------- ---------S HAADKKLVVA 220 Deinococcu -MKKSLLSLK LSGLLVPSVL ALS------- -LSACSSPSS TLNQGTLKIA 221 HISJ_E_COL MKKLVLSLSL VLAFSSATAA F--------- ---------- AAIPQNIRIG 222 223 MSGRYFPFTF VKQ--DKLQG FEVDMWDEIG KRNDYKIEYV TANFSGLFGL 224 ATGQSYPFAY KEN--GKLTG FDVEVMEAVA KKIDMKLDWK LLEFSGLMGE 225 TEGTYEPFTY HDKDTDKLTG YDVEVITEVA KRLGLKVDFK ETQWGSMFAG 226 TEGTYAPFTF HDK-SGKLTG FDVEVIRKVA EKLGLKVEFK ETQWDAMYAG 227 LEGTYPPFSF QGD-DGKLTG FEVEFAQQLA KHLGVEASLK PTKWDGMLAS 228 TDTAFVPFEF KQG--DKYVG FDVDLWAAIA KELKLDYELK PMDFSGIIPA 229 MEGTYPPFTS KNE-QGELVG FDVDIAKAVA QKLNLKPEFV LTEWSGILAG 230 TDPTYAPFES KNS-QGELVG FDIDLAKELC KRINTQCTFV ENPLDALIPS 231 232 LETGRIDTIS NQITMTDARK AKYLFADPYV VDG-AQITVR KGNDSIQGVE 233 LQTGKLDTIS NQVAVTDERK ETYNFTKPYA YAG-TQIVVK KDNTDIKSVD 234 LNSKRFDVVA NQVG-KTDRE DKYDFSDKYT TSR-AVVVTK KDNNDIKSEA 235 LNAKRFDVIA NQTNPSPERL KKYSFTTPYN YSG-GVIVTK SSDNSIKSFE 236 LDSKRIDVVI NQVTISDERK KKYDFSTPYT ISGIQALVKK GNEGTIKTAD 237 LQTKNVDLAL AGITITDERK KAIDFSDGYY KSG-LLVMVK ANNNDVKSVK 238 LQANKYDVIV NQVGITPERQ NSIGFSQPYA YSRPEIIVAK NNTFNPQSLA 239 LKAKKIDAIM SSLSITEKRQ QEIAFTDKLY AADSRLVVAK NSDIQP-TVE 240 241 DLAGKTVAVN LGSNFEQLLR DYDKDGKINI KTYDT--GIE HDVALGRADA 242 DLKGKTVAAV LGSNHAKNLE SKDPDKKINI KTYETQEGTL KDVAYGRVDA 243 DVKGKTSAQS LTSNYNKLAT N----AGAKV EGVEGMAQAL QMIQQARVDM 244 DLKGRKSAQS ATSNWGKDAK A----AGAQI LVVDGLAQSL ELIKQGRAEA 245 DLKGKKVGVG LGTNYEEWLR QNV--QGVDV RTYDDDPTKY QDLRVGRIDA 246 DLDGKVVAVK SGTGSVDYAK AN--IKTKDL RQFPNIDNAY MELGTNRADA 247 DLKGKRVGST LGSNYEKQLI DTG---DIKI VTYPGAPEIL ADLVAGRIDA 248 SLKGKRVGVL QGTTQETFGN EHWAPKGIEI VSYQGQDNIY SDLTAGRIDA 249 250 FIMDRLSALE -LIKKT-GLP LQLAGEPFET I-----QNAW PFVDNEKGRK 251 YVNSRTVLIA -QIKKT-GLP LKLAGDPIVY E-----QVAF PFAKDDAHDK 252 TYNDKLAVLN -YLKTSGNKN VKIAFETGEP Q-----STYF TFRKGS--GE 253 TINDKLAVLD -YFKQHPNSG LKIAYDRGDK T-----PTAF AFLQGE--DA 254 ILVDRLAALD -LVKKT-NDT LAVTGEAFSR Q-----ESGV ALRKGN--ED 255 VLHDTPNILY -FIKTAGNGQ FKAVGDSLEA Q-----QYGI AFPKGS--DE 256 AYNDRLVVNY -IINDQ-KLP VRGAGQIGDA A-----PVGI ALKKGN--SA 257 AFQDEVAASE GFLKQPVGKD YKFGGPSVKD EKLFGVGTGM GLRKED--NE 258 259 LQAEVNKALA EMRADGTVEK ISVKWFGADI TK---- 260 LRKKVNKALD ELRKDGTLKK LSEKYFNEDI TVEQKH 261 VVDQVNKALK EMKEDGTLSK ISKKWFGEDV SK---- 262 LITKFNQVLE ALRQDGTLKQ ISIEWFGYDI TQ---- 263 LLKAVNDAIA EMQKDGTLQA LSEKWFGADV TK---- 264 LRDKVNGALK TLRENGTYNE IYKKWFGTEP K----- 265 LKDQIDKALT EMRSDGTFEK ISQKWFGQDV GQP--- 266 LREALNKAFA EMRADGTYEK LAKKYFDFDV YGG--- 267 """ 268 269 from cStringIO import StringIO 270 handle = StringIO(phylip_text) 271 count=0 272 for alignment in PhylipIterator(handle) : 273 for record in alignment.get_all_seqs() : 274 count=count+1 275 print record.id 276 #print record.seq.tostring() 277 assert count == 8 278 279 expected="""mkklvlslsl vlafssataa faaipqniri gtdptyapfe sknsqgelvg 280 fdidlakelc krintqctfv enpldalips lkakkidaim sslsitekrq qeiaftdkly 281 aadsrlvvak nsdiqptves lkgkrvgvlq gttqetfgne hwapkgieiv syqgqdniys 282 dltagridaafqdevaaseg flkqpvgkdy kfggpsvkde klfgvgtgmg lrkednelre 283 alnkafaemradgtyeklak kyfdfdvygg""".replace(" ","").replace("\n","").upper() 284 assert record.seq.tostring().replace("-","") == expected 285 286 #From here: 287 #http://atgc.lirmm.fr/phyml/usersguide.html 288 phylip_text2="""5 60 289 Tax1 CCATCTCACGGTCGGTACGATACACCTGCTTTTGGCAG 290 Tax2 CCATCTCACGGTCAGTAAGATACACCTGCTTTTGGCGG 291 Tax3 CCATCTCCCGCTCAGTAAGATACCCCTGCTGTTGGCGG 292 Tax4 TCATCTCATGGTCAATAAGATACTCCTGCTTTTGGCGG 293 Tax5 CCATCTCACGGTCGGTAAGATACACCTGCTTTTGGCGG 294 295 GAAATGGTCAATATTACAAGGT 296 GAAATGGTCAACATTAAAAGAT 297 GAAATCGTCAATATTAAAAGGT 298 GAAATGGTCAATCTTAAAAGGT 299 GAAATGGTCAATATTAAAAGGT""" 300 301 phylip_text3="""5 60 302 Tax1 CCATCTCACGGTCGGTACGATACACCTGCTTTTGGCAGGAAATGGTCAATATTACAAGGT 303 Tax2 CCATCTCACGGTCAGTAAGATACACCTGCTTTTGGCGGGAAATGGTCAACATTAAAAGAT 304 Tax3 CCATCTCCCGCTCAGTAAGATACCCCTGCTGTTGGCGGGAAATCGTCAATATTAAAAGGT 305 Tax4 TCATCTCATGGTCAATAAGATACTCCTGCTTTTGGCGGGAAATGGTCAATCTTAAAAGGT 306 Tax5 CCATCTCACGGTCGGTAAGATACACCTGCTTTTGGCGGGAAATGGTCAATATTAAAAGGT""" 307 308 handle = StringIO(phylip_text2) 309 list2 = list(PhylipIterator(handle)) 310 handle.close() 311 assert len(list2)==1 312 assert len(list2[0].get_all_seqs())==5 313 314 handle = StringIO(phylip_text3) 315 list3 = list(PhylipIterator(handle)) 316 handle.close() 317 assert len(list3)==1 318 assert len(list3[0].get_all_seqs())==5 319 320 for i in range(0,5) : 321 list2[0].get_all_seqs()[i].id == list3[0].get_all_seqs()[i].id 322 list2[0].get_all_seqs()[i].seq.tostring() == list3[0].get_all_seqs()[i].seq.tostring() 323 324 #From here: 325 #http://evolution.genetics.washington.edu/phylip/doc/sequence.html 326 #Note the lack of any white space between names 2 and 3 and their seqs. 327 phylip_text4=""" 5 42 328 Turkey AAGCTNGGGC ATTTCAGGGT 329 Salmo gairAAGCCTTGGC AGTGCAGGGT 330 H. SapiensACCGGTTGGC CGTTCAGGGT 331 Chimp AAACCCTTGC CGTTACGCTT 332 Gorilla AAACCCTTGC CGGTACGCTT 333 334 GAGCCCGGGC AATACAGGGT AT 335 GAGCCGTGGC CGGGCACGGT AT 336 ACAGGTTGGC CGTTCAGGGT AA 337 AAACCGAGGC CGGGACACTC AT 338 AAACCATTGC CGGTACGCTT AA""" 339 340 #From here: 341 #http://evolution.genetics.washington.edu/phylip/doc/sequence.html 342 phylip_text5=""" 5 42 343 Turkey AAGCTNGGGC ATTTCAGGGT 344 GAGCCCGGGC AATACAGGGT AT 345 Salmo gairAAGCCTTGGC AGTGCAGGGT 346 GAGCCGTGGC CGGGCACGGT AT 347 H. SapiensACCGGTTGGC CGTTCAGGGT 348 ACAGGTTGGC CGTTCAGGGT AA 349 Chimp AAACCCTTGC CGTTACGCTT 350 AAACCGAGGC CGGGACACTC AT 351 Gorilla AAACCCTTGC CGGTACGCTT 352 AAACCATTGC CGGTACGCTT AA""" 353 354 phylip_text5a=""" 5 42 355 Turkey AAGCTNGGGC ATTTCAGGGT GAGCCCGGGC AATACAGGGT AT 356 Salmo gairAAGCCTTGGC AGTGCAGGGT GAGCCGTGGC CGGGCACGGT AT 357 H. SapiensACCGGTTGGC CGTTCAGGGT ACAGGTTGGC CGTTCAGGGT AA 358 Chimp AAACCCTTGC CGTTACGCTT AAACCGAGGC CGGGACACTC AT 359 Gorilla AAACCCTTGC CGGTACGCTT AAACCATTGC CGGTACGCTT AA""" 360 361 handle = StringIO(phylip_text4) 362 list4 = list(PhylipIterator(handle)) 363 handle.close() 364 assert len(list4)==1 365 assert len(list4[0].get_all_seqs())==5 366 367 handle = StringIO(phylip_text5) 368 try : 369 list5 = list(PhylipIterator(handle)) 370 assert len(list5)==1 371 assert len(list5[0].get_all_seqs())==5 372 print "That should have failed..." 373 except ValueError : 374 print "Evil multiline non-interlaced example failed as expected" 375 handle.close() 376 377 handle = StringIO(phylip_text5a) 378 list5 = list(PhylipIterator(handle)) 379 handle.close() 380 assert len(list5)==1 381 assert len(list4[0].get_all_seqs())==5 382 383 print "Concatenation" 384 handle = StringIO(phylip_text4 + "\n" + phylip_text4) 385 assert len(list(PhylipIterator(handle))) == 2 386 387 handle = StringIO(phylip_text3 + "\n" + phylip_text4 + "\n\n\n" + phylip_text) 388 assert len(list(PhylipIterator(handle))) == 3 389 390 print "OK" 391 392 print "Checking write/read" 393 handle = StringIO() 394 PhylipWriter(handle).write_file(list5) 395 handle.seek(0) 396 list6 = list(PhylipIterator(handle)) 397 assert len(list5) == len(list6) 398 for a1,a2 in zip(list5, list6) : 399 assert len(a1.get_all_seqs()) == len(a2.get_all_seqs()) 400 for r1, r2 in zip(a1.get_all_seqs(), a2.get_all_seqs()) : 401 assert r1.id == r2.id 402 assert r1.seq.tostring() == r2.seq.tostring() 403 print "Done" 404