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

Source Code for Package Bio.Align

  1  # Copyright 2008-2011 by Peter Cock. 
  2  # All rights reserved. 
  3  # This code is part of the Biopython distribution and governed by its 
  4  # license.  Please see the LICENSE file that should have been included 
  5  # as part of this package. 
  6  """Code for dealing with sequence alignments. 
  7   
  8  One of the most important things in this module is the MultipleSeqAlignment 
  9  class, used in the Bio.AlignIO module. 
 10   
 11  """ 
 12  from __future__ import print_function 
 13   
 14  from Bio.Seq import Seq 
 15  from Bio.SeqRecord import SeqRecord 
 16  from Bio import Alphabet 
 17   
 18  # We only import this and subclass it for some limited backward compatibility. 
 19  from Bio.Align.Generic import Alignment as _Alignment 
 20   
 21   
22 -class MultipleSeqAlignment(_Alignment):
23 """Represents a classical multiple sequence alignment (MSA). 24 25 By this we mean a collection of sequences (usually shown as rows) which 26 are all the same length (usually with gap characters for insertions or 27 padding). The data can then be regarded as a matrix of letters, with well 28 defined columns. 29 30 You would typically create an MSA by loading an alignment file with the 31 AlignIO module: 32 33 >>> from Bio import AlignIO 34 >>> align = AlignIO.read("Clustalw/opuntia.aln", "clustal") 35 >>> print(align) 36 SingleLetterAlphabet() alignment with 7 rows and 156 columns 37 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191 38 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191 39 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191 40 TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191 41 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191 42 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191 43 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191 44 45 In some respects you can treat these objects as lists of SeqRecord objects, 46 each representing a row of the alignment. Iterating over an alignment gives 47 the SeqRecord object for each row: 48 49 >>> len(align) 50 7 51 >>> for record in align: 52 ... print("%s %i" % (record.id, len(record))) 53 gi|6273285|gb|AF191659.1|AF191 156 54 gi|6273284|gb|AF191658.1|AF191 156 55 gi|6273287|gb|AF191661.1|AF191 156 56 gi|6273286|gb|AF191660.1|AF191 156 57 gi|6273290|gb|AF191664.1|AF191 156 58 gi|6273289|gb|AF191663.1|AF191 156 59 gi|6273291|gb|AF191665.1|AF191 156 60 61 You can also access individual rows as SeqRecord objects via their index: 62 63 >>> print(align[0].id) 64 gi|6273285|gb|AF191659.1|AF191 65 >>> print(align[-1].id) 66 gi|6273291|gb|AF191665.1|AF191 67 68 And extract columns as strings: 69 70 >>> print(align[:, 1]) 71 AAAAAAA 72 73 Or, take just the first ten columns as a sub-alignment: 74 75 >>> print(align[:, :10]) 76 SingleLetterAlphabet() alignment with 7 rows and 10 columns 77 TATACATTAA gi|6273285|gb|AF191659.1|AF191 78 TATACATTAA gi|6273284|gb|AF191658.1|AF191 79 TATACATTAA gi|6273287|gb|AF191661.1|AF191 80 TATACATAAA gi|6273286|gb|AF191660.1|AF191 81 TATACATTAA gi|6273290|gb|AF191664.1|AF191 82 TATACATTAA gi|6273289|gb|AF191663.1|AF191 83 TATACATTAA gi|6273291|gb|AF191665.1|AF191 84 85 Combining this alignment slicing with alignment addition allows you to 86 remove a section of the alignment. For example, taking just the first 87 and last ten columns: 88 89 >>> print(align[:, :10] + align[:, -10:]) 90 SingleLetterAlphabet() alignment with 7 rows and 20 columns 91 TATACATTAAGTGTACCAGA gi|6273285|gb|AF191659.1|AF191 92 TATACATTAAGTGTACCAGA gi|6273284|gb|AF191658.1|AF191 93 TATACATTAAGTGTACCAGA gi|6273287|gb|AF191661.1|AF191 94 TATACATAAAGTGTACCAGA gi|6273286|gb|AF191660.1|AF191 95 TATACATTAAGTGTACCAGA gi|6273290|gb|AF191664.1|AF191 96 TATACATTAAGTATACCAGA gi|6273289|gb|AF191663.1|AF191 97 TATACATTAAGTGTACCAGA gi|6273291|gb|AF191665.1|AF191 98 99 Note - This object is intended to replace the existing Alignment object 100 defined in module Bio.Align.Generic but is not fully backwards compatible 101 with it. 102 103 Note - This object does NOT attempt to model the kind of alignments used 104 in next generation sequencing with multiple sequencing reads which are 105 much shorter than the alignment, and where there is usually a consensus or 106 reference sequence with special status. 107 """ 108
109 - def __init__(self, records, alphabet=None, 110 annotations=None):
111 """Initialize a new MultipleSeqAlignment object. 112 113 Arguments: 114 - records - A list (or iterator) of SeqRecord objects, whose 115 sequences are all the same length. This may be an be an 116 empty list. 117 - alphabet - The alphabet for the whole alignment, typically a gapped 118 alphabet, which should be a super-set of the individual 119 record alphabets. If omitted, a consensus alphabet is 120 used. 121 - annotations - Information about the whole alignment (dictionary). 122 123 You would normally load a MSA from a file using Bio.AlignIO, but you 124 can do this from a list of SeqRecord objects too: 125 126 >>> from Bio.Alphabet import generic_dna 127 >>> from Bio.Seq import Seq 128 >>> from Bio.SeqRecord import SeqRecord 129 >>> from Bio.Align import MultipleSeqAlignment 130 >>> a = SeqRecord(Seq("AAAACGT", generic_dna), id="Alpha") 131 >>> b = SeqRecord(Seq("AAA-CGT", generic_dna), id="Beta") 132 >>> c = SeqRecord(Seq("AAAAGGT", generic_dna), id="Gamma") 133 >>> align = MultipleSeqAlignment([a, b, c], annotations={"tool": "demo"}) 134 >>> print(align) 135 DNAAlphabet() alignment with 3 rows and 7 columns 136 AAAACGT Alpha 137 AAA-CGT Beta 138 AAAAGGT Gamma 139 >>> align.annotations 140 {'tool': 'demo'} 141 142 NOTE - The older Bio.Align.Generic.Alignment class only accepted a 143 single argument, an alphabet. This is still supported via a backwards 144 compatible "hack" so as not to disrupt existing scripts and users, but 145 is deprecated and will be removed in a future release. 146 """ 147 if isinstance(records, (Alphabet.Alphabet, Alphabet.AlphabetEncoder)): 148 if alphabet is None: 149 # TODO - Remove this backwards compatible mode! 150 alphabet = records 151 records = [] 152 import warnings 153 from Bio import BiopythonDeprecationWarning 154 warnings.warn("Invalid records argument: While the old " 155 "Bio.Align.Generic.Alignment class only " 156 "accepted a single argument (the alphabet), the " 157 "newer Bio.Align.MultipleSeqAlignment class " 158 "expects a list/iterator of SeqRecord objects " 159 "(which can be an empty list) and an optional " 160 "alphabet argument", BiopythonDeprecationWarning) 161 else: 162 raise ValueError("Invalid records argument") 163 if alphabet is not None: 164 if not isinstance(alphabet, (Alphabet.Alphabet, Alphabet.AlphabetEncoder)): 165 raise ValueError("Invalid alphabet argument") 166 self._alphabet = alphabet 167 else: 168 # Default while we add sequences, will take a consensus later 169 self._alphabet = Alphabet.single_letter_alphabet 170 171 self._records = [] 172 if records: 173 self.extend(records) 174 if alphabet is None: 175 # No alphabet was given, take a consensus alphabet 176 self._alphabet = Alphabet._consensus_alphabet(rec.seq.alphabet for 177 rec in self._records 178 if rec.seq is not None) 179 180 # Annotations about the whole alignment 181 if annotations is None: 182 annotations = {} 183 elif not isinstance(annotations, dict): 184 raise TypeError("annotations argument should be a dict") 185 self.annotations = annotations
186
187 - def extend(self, records):
188 """Add more SeqRecord objects to the alignment as rows. 189 190 They must all have the same length as the original alignment, and have 191 alphabets compatible with the alignment's alphabet. For example, 192 193 >>> from Bio.Alphabet import generic_dna 194 >>> from Bio.Seq import Seq 195 >>> from Bio.SeqRecord import SeqRecord 196 >>> from Bio.Align import MultipleSeqAlignment 197 >>> a = SeqRecord(Seq("AAAACGT", generic_dna), id="Alpha") 198 >>> b = SeqRecord(Seq("AAA-CGT", generic_dna), id="Beta") 199 >>> c = SeqRecord(Seq("AAAAGGT", generic_dna), id="Gamma") 200 >>> d = SeqRecord(Seq("AAAACGT", generic_dna), id="Delta") 201 >>> e = SeqRecord(Seq("AAA-GGT", generic_dna), id="Epsilon") 202 203 First we create a small alignment (three rows): 204 205 >>> align = MultipleSeqAlignment([a, b, c]) 206 >>> print(align) 207 DNAAlphabet() alignment with 3 rows and 7 columns 208 AAAACGT Alpha 209 AAA-CGT Beta 210 AAAAGGT Gamma 211 212 Now we can extend this alignment with another two rows: 213 214 >>> align.extend([d, e]) 215 >>> print(align) 216 DNAAlphabet() alignment with 5 rows and 7 columns 217 AAAACGT Alpha 218 AAA-CGT Beta 219 AAAAGGT Gamma 220 AAAACGT Delta 221 AAA-GGT Epsilon 222 223 Because the alignment object allows iteration over the rows as 224 SeqRecords, you can use the extend method with a second alignment 225 (provided its sequences have the same length as the original alignment). 226 """ 227 if len(self): 228 # Use the standard method to get the length 229 expected_length = self.get_alignment_length() 230 else: 231 # Take the first record's length 232 records = iter(records) # records arg could be list or iterator 233 try: 234 rec = next(records) 235 except StopIteration: 236 # Special case, no records 237 return 238 expected_length = len(rec) 239 self._append(rec, expected_length) 240 # Now continue to the rest of the records as usual 241 242 for rec in records: 243 self._append(rec, expected_length)
244
245 - def append(self, record):
246 """Add one more SeqRecord object to the alignment as a new row. 247 248 This must have the same length as the original alignment (unless this is 249 the first record), and have an alphabet compatible with the alignment's 250 alphabet. 251 252 >>> from Bio import AlignIO 253 >>> align = AlignIO.read("Clustalw/opuntia.aln", "clustal") 254 >>> print(align) 255 SingleLetterAlphabet() alignment with 7 rows and 156 columns 256 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191 257 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191 258 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191 259 TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191 260 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191 261 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191 262 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191 263 >>> len(align) 264 7 265 266 We'll now construct a dummy record to append as an example: 267 268 >>> from Bio.Seq import Seq 269 >>> from Bio.SeqRecord import SeqRecord 270 >>> dummy = SeqRecord(Seq("N"*156), id="dummy") 271 272 Now append this to the alignment, 273 274 >>> align.append(dummy) 275 >>> print(align) 276 SingleLetterAlphabet() alignment with 8 rows and 156 columns 277 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273285|gb|AF191659.1|AF191 278 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273284|gb|AF191658.1|AF191 279 TATACATTAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273287|gb|AF191661.1|AF191 280 TATACATAAAAGAAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273286|gb|AF191660.1|AF191 281 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273290|gb|AF191664.1|AF191 282 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273289|gb|AF191663.1|AF191 283 TATACATTAAAGGAGGGGGATGCGGATAAATGGAAAGGCGAAAG...AGA gi|6273291|gb|AF191665.1|AF191 284 NNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNN...NNN dummy 285 >>> len(align) 286 8 287 288 """ 289 if self._records: 290 self._append(record, self.get_alignment_length()) 291 else: 292 self._append(record)
293
294 - def _append(self, record, expected_length=None):
295 """Helper function (PRIVATE).""" 296 if not isinstance(record, SeqRecord): 297 raise TypeError("New sequence is not a SeqRecord object") 298 299 # Currently the get_alignment_length() call is expensive, so we need 300 # to avoid calling it repeatedly for __init__ and extend, hence this 301 # private _append method 302 if expected_length is not None and len(record) != expected_length: 303 # TODO - Use the following more helpful error, but update unit tests 304 # raise ValueError("New sequence is not of length %i" \ 305 # % self.get_alignment_length()) 306 raise ValueError("Sequences must all be the same length") 307 308 # Using not self.alphabet.contains(record.seq.alphabet) needs fixing 309 # for AlphabetEncoders (e.g. gapped versus ungapped). 310 if not Alphabet._check_type_compatible([self._alphabet, record.seq.alphabet]): 311 raise ValueError("New sequence's alphabet is incompatible") 312 self._records.append(record)
313
314 - def __add__(self, other):
315 """Combines two alignments with the same number of rows by adding them. 316 317 If you have two multiple sequence alignments (MSAs), there are two ways to think 318 about adding them - by row or by column. Using the extend method adds by row. 319 Using the addition operator adds by column. For example, 320 321 >>> from Bio.Alphabet import generic_dna 322 >>> from Bio.Seq import Seq 323 >>> from Bio.SeqRecord import SeqRecord 324 >>> from Bio.Align import MultipleSeqAlignment 325 >>> a1 = SeqRecord(Seq("AAAAC", generic_dna), id="Alpha") 326 >>> b1 = SeqRecord(Seq("AAA-C", generic_dna), id="Beta") 327 >>> c1 = SeqRecord(Seq("AAAAG", generic_dna), id="Gamma") 328 >>> a2 = SeqRecord(Seq("GT", generic_dna), id="Alpha") 329 >>> b2 = SeqRecord(Seq("GT", generic_dna), id="Beta") 330 >>> c2 = SeqRecord(Seq("GT", generic_dna), id="Gamma") 331 >>> left = MultipleSeqAlignment([a1, b1, c1], 332 ... annotations={"tool": "demo", "name": "start"}) 333 >>> right = MultipleSeqAlignment([a2, b2, c2], 334 ... annotations={"tool": "demo", "name": "end"}) 335 336 Now, let's look at these two alignments: 337 338 >>> print(left) 339 DNAAlphabet() alignment with 3 rows and 5 columns 340 AAAAC Alpha 341 AAA-C Beta 342 AAAAG Gamma 343 >>> print(right) 344 DNAAlphabet() alignment with 3 rows and 2 columns 345 GT Alpha 346 GT Beta 347 GT Gamma 348 349 And add them: 350 351 >>> combined = left + right 352 >>> print(combined) 353 DNAAlphabet() alignment with 3 rows and 7 columns 354 AAAACGT Alpha 355 AAA-CGT Beta 356 AAAAGGT Gamma 357 358 For this to work, both alignments must have the same number of records (here 359 they both have 3 rows): 360 361 >>> len(left) 362 3 363 >>> len(right) 364 3 365 >>> len(combined) 366 3 367 368 The individual rows are SeqRecord objects, and these can be added together. Refer 369 to the SeqRecord documentation for details of how the annotation is handled. This 370 example is a special case in that both original alignments shared the same names, 371 meaning when the rows are added they also get the same name. 372 373 Any common annotations are preserved, but differing annotation is lost. This is 374 the same behaviour used in the SeqRecord annotations and is designed to prevent 375 accidental propagation of inappropriate values: 376 377 >>> combined.annotations 378 {'tool': 'demo'} 379 380 """ 381 if not isinstance(other, MultipleSeqAlignment): 382 raise NotImplementedError 383 if len(self) != len(other): 384 raise ValueError("When adding two alignments they must have the same length" 385 " (i.e. same number or rows)") 386 alpha = Alphabet._consensus_alphabet([self._alphabet, other._alphabet]) 387 merged = (left + right for left, right in zip(self, other)) 388 # Take any common annotation: 389 annotations = dict() 390 for k, v in self.annotations.items(): 391 if k in other.annotations and other.annotations[k] == v: 392 annotations[k] = v 393 return MultipleSeqAlignment(merged, alpha, annotations)
394
395 - def __getitem__(self, index):
396 """Access part of the alignment. 397 398 Depending on the indices, you can get a SeqRecord object 399 (representing a single row), a Seq object (for a single columns), 400 a string (for a single characters) or another alignment 401 (representing some part or all of the alignment). 402 403 align[r,c] gives a single character as a string 404 align[r] gives a row as a SeqRecord 405 align[r,:] gives a row as a SeqRecord 406 align[:,c] gives a column as a Seq (using the alignment's alphabet) 407 408 align[:] and align[:,:] give a copy of the alignment 409 410 Anything else gives a sub alignment, e.g. 411 align[0:2] or align[0:2,:] uses only row 0 and 1 412 align[:,1:3] uses only columns 1 and 2 413 align[0:2,1:3] uses only rows 0 & 1 and only cols 1 & 2 414 415 We'll use the following example alignment here for illustration: 416 417 >>> from Bio.Alphabet import generic_dna 418 >>> from Bio.Seq import Seq 419 >>> from Bio.SeqRecord import SeqRecord 420 >>> from Bio.Align import MultipleSeqAlignment 421 >>> a = SeqRecord(Seq("AAAACGT", generic_dna), id="Alpha") 422 >>> b = SeqRecord(Seq("AAA-CGT", generic_dna), id="Beta") 423 >>> c = SeqRecord(Seq("AAAAGGT", generic_dna), id="Gamma") 424 >>> d = SeqRecord(Seq("AAAACGT", generic_dna), id="Delta") 425 >>> e = SeqRecord(Seq("AAA-GGT", generic_dna), id="Epsilon") 426 >>> align = MultipleSeqAlignment([a, b, c, d, e], generic_dna) 427 428 You can access a row of the alignment as a SeqRecord using an integer 429 index (think of the alignment as a list of SeqRecord objects here): 430 431 >>> first_record = align[0] 432 >>> print("%s %s" % (first_record.id, first_record.seq)) 433 Alpha AAAACGT 434 >>> last_record = align[-1] 435 >>> print("%s %s" % (last_record.id, last_record.seq)) 436 Epsilon AAA-GGT 437 438 You can also access use python's slice notation to create a sub-alignment 439 containing only some of the SeqRecord objects: 440 441 >>> sub_alignment = align[2:5] 442 >>> print(sub_alignment) 443 DNAAlphabet() alignment with 3 rows and 7 columns 444 AAAAGGT Gamma 445 AAAACGT Delta 446 AAA-GGT Epsilon 447 448 This includes support for a step, i.e. align[start:end:step], which 449 can be used to select every second sequence: 450 451 >>> sub_alignment = align[::2] 452 >>> print(sub_alignment) 453 DNAAlphabet() alignment with 3 rows and 7 columns 454 AAAACGT Alpha 455 AAAAGGT Gamma 456 AAA-GGT Epsilon 457 458 Or to get a copy of the alignment with the rows in reverse order: 459 460 >>> rev_alignment = align[::-1] 461 >>> print(rev_alignment) 462 DNAAlphabet() alignment with 5 rows and 7 columns 463 AAA-GGT Epsilon 464 AAAACGT Delta 465 AAAAGGT Gamma 466 AAA-CGT Beta 467 AAAACGT Alpha 468 469 You can also use two indices to specify both rows and columns. Using simple 470 integers gives you the entry as a single character string. e.g. 471 472 >>> align[3, 4] 473 'C' 474 475 This is equivalent to: 476 477 >>> align[3][4] 478 'C' 479 480 or: 481 482 >>> align[3].seq[4] 483 'C' 484 485 To get a single column (as a string) use this syntax: 486 487 >>> align[:, 4] 488 'CCGCG' 489 490 Or, to get part of a column, 491 492 >>> align[1:3, 4] 493 'CG' 494 495 However, in general you get a sub-alignment, 496 497 >>> print(align[1:5, 3:6]) 498 DNAAlphabet() alignment with 4 rows and 3 columns 499 -CG Beta 500 AGG Gamma 501 ACG Delta 502 -GG Epsilon 503 504 This should all seem familiar to anyone who has used the NumPy 505 array or matrix objects. 506 """ 507 if isinstance(index, int): 508 # e.g. result = align[x] 509 # Return a SeqRecord 510 return self._records[index] 511 elif isinstance(index, slice): 512 # e.g. sub_align = align[i:j:k] 513 return MultipleSeqAlignment(self._records[index], self._alphabet) 514 elif len(index) != 2: 515 raise TypeError("Invalid index type.") 516 517 # Handle double indexing 518 row_index, col_index = index 519 if isinstance(row_index, int): 520 # e.g. row_or_part_row = align[6, 1:4], gives a SeqRecord 521 return self._records[row_index][col_index] 522 elif isinstance(col_index, int): 523 # e.g. col_or_part_col = align[1:5, 6], gives a string 524 return "".join(rec[col_index] for rec in self._records[row_index]) 525 else: 526 # e.g. sub_align = align[1:4, 5:7], gives another alignment 527 return MultipleSeqAlignment((rec[col_index] for rec in self._records[row_index]), 528 self._alphabet)
529
530 - def sort(self, key=None, reverse=False):
531 """Sort the rows (SeqRecord objects) of the alignment in place. 532 533 This sorts the rows alphabetically using the SeqRecord object id by 534 default. The sorting can be controlled by supplying a key function 535 which must map each SeqRecord to a sort value. 536 537 This is useful if you want to add two alignments which use the same 538 record identifiers, but in a different order. For example, 539 540 >>> from Bio.Alphabet import generic_dna 541 >>> from Bio.Seq import Seq 542 >>> from Bio.SeqRecord import SeqRecord 543 >>> from Bio.Align import MultipleSeqAlignment 544 >>> align1 = MultipleSeqAlignment([ 545 ... SeqRecord(Seq("ACGT", generic_dna), id="Human"), 546 ... SeqRecord(Seq("ACGG", generic_dna), id="Mouse"), 547 ... SeqRecord(Seq("ACGC", generic_dna), id="Chicken"), 548 ... ]) 549 >>> align2 = MultipleSeqAlignment([ 550 ... SeqRecord(Seq("CGGT", generic_dna), id="Mouse"), 551 ... SeqRecord(Seq("CGTT", generic_dna), id="Human"), 552 ... SeqRecord(Seq("CGCT", generic_dna), id="Chicken"), 553 ... ]) 554 555 If you simple try and add these without sorting, you get this: 556 557 >>> print(align1 + align2) 558 DNAAlphabet() alignment with 3 rows and 8 columns 559 ACGTCGGT <unknown id> 560 ACGGCGTT <unknown id> 561 ACGCCGCT Chicken 562 563 Consult the SeqRecord documentation which explains why you get a 564 default value when annotation like the identifier doesn't match up. 565 However, if we sort the alignments first, then add them we get the 566 desired result: 567 568 >>> align1.sort() 569 >>> align2.sort() 570 >>> print(align1 + align2) 571 DNAAlphabet() alignment with 3 rows and 8 columns 572 ACGCCGCT Chicken 573 ACGTCGTT Human 574 ACGGCGGT Mouse 575 576 As an example using a different sort order, you could sort on the 577 GC content of each sequence. 578 579 >>> from Bio.SeqUtils import GC 580 >>> print(align1) 581 DNAAlphabet() alignment with 3 rows and 4 columns 582 ACGC Chicken 583 ACGT Human 584 ACGG Mouse 585 >>> align1.sort(key = lambda record: GC(record.seq)) 586 >>> print(align1) 587 DNAAlphabet() alignment with 3 rows and 4 columns 588 ACGT Human 589 ACGC Chicken 590 ACGG Mouse 591 592 There is also a reverse argument, so if you wanted to sort by ID 593 but backwards: 594 595 >>> align1.sort(reverse=True) 596 >>> print(align1) 597 DNAAlphabet() alignment with 3 rows and 4 columns 598 ACGG Mouse 599 ACGT Human 600 ACGC Chicken 601 602 """ 603 if key is None: 604 self._records.sort(key=lambda r: r.id, reverse=reverse) 605 else: 606 self._records.sort(key=key, reverse=reverse)
607
608 - def get_column(self, col):
609 """Returns a string containing a given column (DEPRECATED). 610 611 This is a method provided for backwards compatibility with the old 612 Bio.Align.Generic.Alignment object. Please use the slice notation 613 instead, since get_column is likely to be removed in a future release 614 of Biopython.. 615 """ 616 import warnings 617 from Bio import BiopythonDeprecationWarning 618 warnings.warn("This method is deprecated and is provided for backwards " 619 "compatibility with the old Bio.Align.Generic.Alignment " 620 "object. Please use the slice notation instead, as " 621 "get_column is likely to be removed in a future release " 622 "of Biopython.", BiopythonDeprecationWarning) 623 return _Alignment.get_column(self, col)
624
625 - def add_sequence(self, descriptor, sequence, start=None, end=None, 626 weight=1.0):
627 """Add a sequence to the alignment (DEPRECATED). 628 629 The start, end, and weight arguments are not supported! This method 630 only provides limited backwards compatibility with the old 631 Bio.Align.Generic.Alignment object. Please use the append method with 632 a SeqRecord instead, since add_sequence is likely to be removed in a 633 future release of Biopython. 634 """ 635 import warnings 636 from Bio import BiopythonDeprecationWarning 637 warnings.warn("The start, end, and weight arguments are not supported! " 638 "This method only provides limited backwards " 639 "compatibility with the old Bio.Align.Generic.Alignment " 640 "object. Please use the append method with a SeqRecord " 641 "instead, as the add_sequence method is likely to be " 642 "removed in a future release of Biopython.", 643 BiopythonDeprecationWarning) 644 # Should we handle start/end/strand information somehow? What for? 645 # TODO - Should we handle weights somehow? See also AlignInfo code... 646 if start is not None or end is not None or weight != 1.0: 647 raise ValueError("The add_Sequence method is obsolete, and only " 648 "provides limited backwards compatibily. The" 649 "start, end and weight arguments are not " 650 "supported.") 651 self.append(SeqRecord(Seq(sequence, self._alphabet), 652 id=descriptor, description=descriptor))
653 654 655 if __name__ == "__main__": 656 from Bio._utils import run_doctest 657 run_doctest() 658