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

Source Code for Package Bio.AlignIO

  1  # Copyright 2008-2010 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  """Multiple sequence alignment input/output as alignment objects. 
  7   
  8  The Bio.AlignIO interface is deliberately very similar to Bio.SeqIO, and in 
  9  fact the two are connected internally.  Both modules use the same set of file 
 10  format names (lower case strings).  From the user's perspective, you can read 
 11  in a PHYLIP file containing one or more alignments using Bio.AlignIO, or you 
 12  can read in the sequences within these alignmenta using Bio.SeqIO. 
 13   
 14  Bio.AlignIO is also documented at http://biopython.org/wiki/AlignIO and by 
 15  a whole chapter in our tutorial: 
 16   
 17  * `HTML Tutorial`_ 
 18  * `PDF Tutorial`_ 
 19   
 20  .. _`HTML Tutorial`: http://biopython.org/DIST/docs/tutorial/Tutorial.html 
 21  .. _`PDF Tutorial`: http://biopython.org/DIST/docs/tutorial/Tutorial.pdf 
 22   
 23  Input 
 24  ----- 
 25  For the typical special case when your file or handle contains one and only 
 26  one alignment, use the function Bio.AlignIO.read().  This takes an input file 
 27  handle (or in recent versions of Biopython a filename as a string), format 
 28  string and optional number of sequences per alignment.  It will return a single 
 29  MultipleSeqAlignment object (or raise an exception if there isn't just one 
 30  alignment): 
 31   
 32  >>> from Bio import AlignIO 
 33  >>> align = AlignIO.read("Phylip/interlaced.phy", "phylip") 
 34  >>> print(align) 
 35  SingleLetterAlphabet() alignment with 3 rows and 384 columns 
 36  -----MKVILLFVLAVFTVFVSS---------------RGIPPE...I-- CYS1_DICDI 
 37  MAHARVLLLALAVLATAAVAVASSSSFADSNPIRPVTDRAASTL...VAA ALEU_HORVU 
 38  ------MWATLPLLCAGAWLLGV--------PVCGAAELSVNSL...PLV CATH_HUMAN 
 39   
 40  For the general case, when the handle could contain any number of alignments, 
 41  use the function Bio.AlignIO.parse(...) which takes the same arguments, but 
 42  returns an iterator giving MultipleSeqAlignment objects (typically used in a 
 43  for loop). If you want random access to the alignments by number, turn this 
 44  into a list: 
 45   
 46  >>> from Bio import AlignIO 
 47  >>> alignments = list(AlignIO.parse("Emboss/needle.txt", "emboss")) 
 48  >>> print(alignments[2]) 
 49  SingleLetterAlphabet() alignment with 2 rows and 120 columns 
 50  -KILIVDDQYGIRILLNEVFNKEGYQTFQAANGLQALDIVTKER...--- ref_rec 
 51  LHIVVVDDDPGTCVYIESVFAELGHTCKSFVRPEAAEEYILTHP...HKE gi|94967506|receiver 
 52   
 53  Most alignment file formats can be concatenated so as to hold as many 
 54  different multiple sequence alignments as possible.  One common example 
 55  is the output of the tool seqboot in the PHLYIP suite.  Sometimes there 
 56  can be a file header and footer, as seen in the EMBOSS alignment output. 
 57   
 58  Output 
 59  ------ 
 60  Use the function Bio.AlignIO.write(...), which takes a complete set of 
 61  Alignment objects (either as a list, or an iterator), an output file handle 
 62  (or filename in recent versions of Biopython) and of course the file format:: 
 63   
 64    from Bio import AlignIO 
 65    alignments = ... 
 66    count = SeqIO.write(alignments, "example.faa", "fasta") 
 67   
 68  If using a handle make sure to close it to flush the data to the disk:: 
 69   
 70    from Bio import AlignIO 
 71    alignments = ... 
 72    with open("example.faa", "w") as handle: 
 73      count = SeqIO.write(alignments, handle, "fasta") 
 74   
 75  In general, you are expected to call this function once (with all your 
 76  alignments) and then close the file handle.  However, for file formats 
 77  like PHYLIP where multiple alignments are stored sequentially (with no file 
 78  header and footer), then multiple calls to the write function should work as 
 79  expected when using handles. 
 80   
 81  If you are using a filename, the repeated calls to the write functions will 
 82  overwrite the existing file each time. 
 83   
 84  Conversion 
 85  ---------- 
 86  The Bio.AlignIO.convert(...) function allows an easy interface for simple 
 87  alignnment file format conversions. Additionally, it may use file format 
 88  specific optimisations so this should be the fastest way too. 
 89   
 90  In general however, you can combine the Bio.AlignIO.parse(...) function with 
 91  the Bio.AlignIO.write(...) function for sequence file conversion. Using 
 92  generator expressions provides a memory efficient way to perform filtering or 
 93  other extra operations as part of the process. 
 94   
 95  File Formats 
 96  ------------ 
 97  When specifying the file format, use lowercase strings.  The same format 
 98  names are also used in Bio.SeqIO and include the following: 
 99   
100    - clustal -   Output from Clustal W or X, see also the module Bio.Clustalw 
101      which can be used to run the command line tool from Biopython. 
102    - emboss    - EMBOSS tools' "pairs" and "simple" alignment formats. 
103    - fasta     - The generic sequence file format where each record starts with 
104      an identifer line starting with a ">" character, followed by 
105      lines of sequence. 
106    - fasta-m10 - For the pairswise alignments output by Bill Pearson's FASTA 
107      tools when used with the -m 10 command line option for machine 
108      readable output. 
109    - ig        - The IntelliGenetics file format, apparently the same as the 
110      MASE alignment format. 
111    - nexus     - Output from NEXUS, see also the module Bio.Nexus which can also 
112      read any phylogenetic trees in these files. 
113    - phylip    - Interlaced PHYLIP, as used by the PHLIP tools. 
114    - phylip-sequential - Sequential PHYLIP. 
115    - phylip-relaxed - PHYLIP like format allowing longer names. 
116    - stockholm - A richly annotated alignment file format used by PFAM. 
117   
118  Note that while Bio.AlignIO can read all the above file formats, it cannot 
119  write to all of them. 
120   
121  You can also use any file format supported by Bio.SeqIO, such as "fasta" or 
122  "ig" (which are listed above), PROVIDED the sequences in your file are all the 
123  same length. 
124  """ 
125   
126   
127  from __future__ import print_function 
128  from Bio._py3k import basestring 
129   
130  __docformat__ = "restructuredtext en"  # not just plaintext 
131   
132  # TODO 
133  # - define policy on reading aligned sequences with gaps in 
134  #   (e.g. - and . characters) including how the alphabet interacts 
135  # 
136  # - Can we build the to_alignment(...) functionality 
137  #   into the generic Alignment class instead? 
138  # 
139  # - How best to handle unique/non unique record.id when writing. 
140  #   For most file formats reading such files is fine; The stockholm 
141  #   parser would fail. 
142  # 
143  # - MSF multiple alignment format, aka GCG, aka PileUp format (*.msf) 
144  #   http://www.bioperl.org/wiki/MSF_multiple_alignment_format 
145   
146  from Bio.Align import MultipleSeqAlignment 
147  from Bio.Align.Generic import Alignment 
148  from Bio.Alphabet import Alphabet, AlphabetEncoder, _get_base_alphabet 
149  from Bio.File import as_handle 
150   
151  from . import StockholmIO 
152  from . import ClustalIO 
153  from . import NexusIO 
154  from . import PhylipIO 
155  from . import EmbossIO 
156  from . import FastaIO 
157   
158  # Convention for format names is "mainname-subtype" in lower case. 
159  # Please use the same names as BioPerl and EMBOSS where possible. 
160   
161  _FormatToIterator = {  # "fasta" is done via Bio.SeqIO 
162                       "clustal": ClustalIO.ClustalIterator, 
163                       "emboss": EmbossIO.EmbossIterator, 
164                       "fasta-m10": FastaIO.FastaM10Iterator, 
165                       "nexus": NexusIO.NexusIterator, 
166                       "phylip": PhylipIO.PhylipIterator, 
167                       "phylip-sequential": PhylipIO.SequentialPhylipIterator, 
168                       "phylip-relaxed": PhylipIO.RelaxedPhylipIterator, 
169                       "stockholm": StockholmIO.StockholmIterator, 
170                       } 
171   
172  _FormatToWriter = {  # "fasta" is done via Bio.SeqIO 
173                       # "emboss" : EmbossIO.EmbossWriter, (unfinished) 
174                     "nexus": NexusIO.NexusWriter, 
175                     "phylip": PhylipIO.PhylipWriter, 
176                     "phylip-sequential": PhylipIO.SequentialPhylipWriter, 
177                     "phylip-relaxed": PhylipIO.RelaxedPhylipWriter, 
178                     "stockholm": StockholmIO.StockholmWriter, 
179                     "clustal": ClustalIO.ClustalWriter, 
180                     } 
181   
182   
183 -def write(alignments, handle, format):
184 """Write complete set of alignments to a file. 185 186 Arguments: 187 - alignments - A list (or iterator) of Alignment objects (ideally the 188 new MultipleSeqAlignment objects), or (if using Biopython 189 1.54 or later) a single alignment object. 190 - handle - File handle object to write to, or filename as string 191 (note older versions of Biopython only took a handle). 192 - format - lower case string describing the file format to write. 193 194 You should close the handle after calling this function. 195 196 Returns the number of alignments written (as an integer). 197 """ 198 from Bio import SeqIO 199 200 # Try and give helpful error messages: 201 if not isinstance(format, basestring): 202 raise TypeError("Need a string for the file format (lower case)") 203 if not format: 204 raise ValueError("Format required (lower case string)") 205 if format != format.lower(): 206 raise ValueError("Format string '%s' should be lower case" % format) 207 208 if isinstance(alignments, Alignment): 209 # This raised an exception in older versions of Biopython 210 alignments = [alignments] 211 212 with as_handle(handle, 'w') as fp: 213 # Map the file format to a writer class 214 if format in _FormatToWriter: 215 writer_class = _FormatToWriter[format] 216 count = writer_class(fp).write_file(alignments) 217 elif format in SeqIO._FormatToWriter: 218 # Exploit the existing SeqIO parser to do the dirty work! 219 # TODO - Can we make one call to SeqIO.write() and count the alignments? 220 count = 0 221 for alignment in alignments: 222 if not isinstance(alignment, Alignment): 223 raise TypeError( 224 "Expect a list or iterator of Alignment objects.") 225 SeqIO.write(alignment, fp, format) 226 count += 1 227 elif format in _FormatToIterator or format in SeqIO._FormatToIterator: 228 raise ValueError("Reading format '%s' is supported, but not writing" 229 % format) 230 else: 231 raise ValueError("Unknown format '%s'" % format) 232 233 assert isinstance(count, int), "Internal error - the underlying %s " \ 234 "writer should have returned the alignment count, not %s" \ 235 % (format, repr(count)) 236 237 return count
238 239 240 # This is a generator function!
241 -def _SeqIO_to_alignment_iterator(handle, format, alphabet=None, seq_count=None):
242 """Uses Bio.SeqIO to create an MultipleSeqAlignment iterator (PRIVATE). 243 244 Arguments: 245 - handle - handle to the file. 246 - format - string describing the file format. 247 - alphabet - optional Alphabet object, useful when the sequence type 248 cannot be automatically inferred from the file itself 249 (e.g. fasta, phylip, clustal) 250 - seq_count - Optional integer, number of sequences expected in each 251 alignment. Recommended for fasta format files. 252 253 If count is omitted (default) then all the sequences in the file are 254 combined into a single MultipleSeqAlignment. 255 """ 256 from Bio import SeqIO 257 assert format in SeqIO._FormatToIterator 258 259 if seq_count: 260 # Use the count to split the records into batches. 261 seq_record_iterator = SeqIO.parse(handle, format, alphabet) 262 263 records = [] 264 for record in seq_record_iterator: 265 records.append(record) 266 if len(records) == seq_count: 267 yield MultipleSeqAlignment(records, alphabet) 268 records = [] 269 if len(records) > 0: 270 raise ValueError("Check seq_count argument, not enough sequences?") 271 else: 272 # Must assume that there is a single alignment using all 273 # the SeqRecord objects: 274 records = list(SeqIO.parse(handle, format, alphabet)) 275 if records: 276 yield MultipleSeqAlignment(records, alphabet) 277 raise StopIteration
278 279
280 -def _force_alphabet(alignment_iterator, alphabet):
281 """Iterate over alignments, over-riding the alphabet (PRIVATE).""" 282 # Assume the alphabet argument has been pre-validated 283 given_base_class = _get_base_alphabet(alphabet).__class__ 284 for align in alignment_iterator: 285 if not isinstance(_get_base_alphabet(align._alphabet), 286 given_base_class): 287 raise ValueError("Specified alphabet %s clashes with " 288 "that determined from the file, %s" 289 % (repr(alphabet), repr(align._alphabet))) 290 for record in align: 291 if not isinstance(_get_base_alphabet(record.seq.alphabet), 292 given_base_class): 293 raise ValueError("Specified alphabet %s clashes with " 294 "that determined from the file, %s" 295 % (repr(alphabet), repr(record.seq.alphabet))) 296 record.seq.alphabet = alphabet 297 align._alphabet = alphabet 298 yield align
299 300
301 -def parse(handle, format, seq_count=None, alphabet=None):
302 """Iterate over an alignment file as MultipleSeqAlignment objects. 303 304 Arguments: 305 - handle - handle to the file, or the filename as a string 306 (note older versions of Biopython only took a handle). 307 - format - string describing the file format. 308 - alphabet - optional Alphabet object, useful when the sequence type 309 cannot be automatically inferred from the file itself 310 (e.g. fasta, phylip, clustal) 311 - seq_count - Optional integer, number of sequences expected in each 312 alignment. Recommended for fasta format files. 313 314 If you have the file name in a string 'filename', use: 315 316 >>> from Bio import AlignIO 317 >>> filename = "Emboss/needle.txt" 318 >>> format = "emboss" 319 >>> for alignment in AlignIO.parse(filename, format): 320 ... print("Alignment of length %i" % alignment.get_alignment_length()) 321 Alignment of length 124 322 Alignment of length 119 323 Alignment of length 120 324 Alignment of length 118 325 Alignment of length 125 326 327 If you have a string 'data' containing the file contents, use:: 328 329 from Bio import AlignIO 330 from StringIO import StringIO 331 my_iterator = AlignIO.parse(StringIO(data), format) 332 333 Use the Bio.AlignIO.read() function when you expect a single record only. 334 """ 335 from Bio import SeqIO 336 337 # Try and give helpful error messages: 338 if not isinstance(format, basestring): 339 raise TypeError("Need a string for the file format (lower case)") 340 if not format: 341 raise ValueError("Format required (lower case string)") 342 if format != format.lower(): 343 raise ValueError("Format string '%s' should be lower case" % format) 344 if alphabet is not None and not (isinstance(alphabet, Alphabet) or 345 isinstance(alphabet, AlphabetEncoder)): 346 raise ValueError("Invalid alphabet, %s" % repr(alphabet)) 347 if seq_count is not None and not isinstance(seq_count, int): 348 raise TypeError("Need integer for seq_count (sequences per alignment)") 349 350 with as_handle(handle, 'rU') as fp: 351 # Map the file format to a sequence iterator: 352 if format in _FormatToIterator: 353 iterator_generator = _FormatToIterator[format] 354 if alphabet is None: 355 i = iterator_generator(fp, seq_count) 356 else: 357 try: 358 # Initially assume the optional alphabet argument is supported 359 i = iterator_generator(fp, seq_count, alphabet=alphabet) 360 except TypeError: 361 # It isn't supported. 362 i = _force_alphabet(iterator_generator(fp, seq_count), 363 alphabet) 364 365 elif format in SeqIO._FormatToIterator: 366 # Exploit the existing SeqIO parser to the dirty work! 367 i = _SeqIO_to_alignment_iterator(fp, format, 368 alphabet=alphabet, 369 seq_count=seq_count) 370 else: 371 raise ValueError("Unknown format '%s'" % format) 372 373 # This imposes some overhead... wait until we drop Python 2.4 to fix it 374 for a in i: 375 yield a
376 377
378 -def read(handle, format, seq_count=None, alphabet=None):
379 """Turns an alignment file into a single MultipleSeqAlignment object. 380 381 Arguments: 382 - handle - handle to the file, or the filename as a string 383 (note older versions of Biopython only took a handle). 384 - format - string describing the file format. 385 - alphabet - optional Alphabet object, useful when the sequence type 386 cannot be automatically inferred from the file itself 387 (e.g. fasta, phylip, clustal) 388 - seq_count - Optional integer, number of sequences expected in each 389 alignment. Recommended for fasta format files. 390 391 If the handle contains no alignments, or more than one alignment, 392 an exception is raised. For example, using a PFAM/Stockholm file 393 containing one alignment: 394 395 >>> from Bio import AlignIO 396 >>> filename = "Clustalw/protein.aln" 397 >>> format = "clustal" 398 >>> alignment = AlignIO.read(filename, format) 399 >>> print("Alignment of length %i" % alignment.get_alignment_length()) 400 Alignment of length 411 401 402 If however you want the first alignment from a file containing 403 multiple alignments this function would raise an exception. 404 405 >>> from Bio import AlignIO 406 >>> filename = "Emboss/needle.txt" 407 >>> format = "emboss" 408 >>> alignment = AlignIO.read(filename, format) 409 Traceback (most recent call last): 410 ... 411 ValueError: More than one record found in handle 412 413 Instead use: 414 415 >>> from Bio import AlignIO 416 >>> filename = "Emboss/needle.txt" 417 >>> format = "emboss" 418 >>> alignment = next(AlignIO.parse(filename, format)) 419 >>> print("First alignment has length %i" % alignment.get_alignment_length()) 420 First alignment has length 124 421 422 You must use the Bio.AlignIO.parse() function if you want to read multiple 423 records from the handle. 424 """ 425 iterator = parse(handle, format, seq_count, alphabet) 426 try: 427 first = next(iterator) 428 except StopIteration: 429 first = None 430 if first is None: 431 raise ValueError("No records found in handle") 432 try: 433 second = next(iterator) 434 except StopIteration: 435 second = None 436 if second is not None: 437 raise ValueError("More than one record found in handle") 438 if seq_count: 439 assert len(first) == seq_count 440 return first
441 442
443 -def convert(in_file, in_format, out_file, out_format, alphabet=None):
444 """Convert between two alignment files, returns number of alignments. 445 446 - in_file - an input handle or filename 447 - in_format - input file format, lower case string 448 - output - an output handle or filename 449 - out_file - output file format, lower case string 450 - alphabet - optional alphabet to assume 451 452 **NOTE** - If you provide an output filename, it will be opened which will 453 overwrite any existing file without warning. This may happen if even the 454 conversion is aborted (e.g. an invalid out_format name is given). 455 """ 456 # TODO - Add optimised versions of important conversions 457 # For now just off load the work to SeqIO parse/write 458 with as_handle(in_file, 'rU') as in_handle: 459 # Don't open the output file until we've checked the input is OK: 460 alignments = parse(in_handle, in_format, None, alphabet) 461 462 # This will check the arguments and issue error messages, 463 # after we have opened the file which is a shame. 464 with as_handle(out_file, 'w') as out_handle: 465 count = write(alignments, out_handle, out_format) 466 467 return count
468 469 470 if __name__ == "__main__": 471 from Bio._utils import run_doctest 472 run_doctest() 473