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

Source Code for Module Bio.Align.AlignInfo

  1  # This code is part of the Biopython distribution and governed by its 
  2  # license.  Please see the LICENSE file that should have been included 
  3  # as part of this package. 
  4  # 
  5   
  6  """Extract information from alignment objects. 
  7   
  8  In order to try and avoid huge alignment objects with tons of functions, 
  9  functions which return summary type information about alignments should 
 10  be put into classes in this module. 
 11  """ 
 12   
 13  from __future__ import print_function 
 14   
 15  import math 
 16  import sys 
 17   
 18  from Bio import Alphabet 
 19  from Bio.Alphabet import IUPAC 
 20  from Bio.Seq import Seq 
 21  from Bio.SubsMat import FreqTable 
 22   
 23  __docformat__ = "restructuredtext en" 
 24   
 25  # Expected random distributions for 20-letter protein, and 
 26  # for 4-letter nucleotide alphabets 
 27  Protein20Random = 0.05 
 28  Nucleotide4Random = 0.25 
 29   
 30   
31 -class SummaryInfo(object):
32 """Calculate summary info about the alignment. 33 34 This class should be used to caclculate information summarizing the 35 results of an alignment. This may either be straight consensus info 36 or more complicated things. 37 """
38 - def __init__(self, alignment):
39 """Initialize with the alignment to calculate information on. 40 41 ic_vector attribute. A dictionary. Keys: column numbers. Values: 42 """ 43 self.alignment = alignment 44 self.ic_vector = {}
45
46 - def dumb_consensus(self, threshold=.7, ambiguous="X", 47 consensus_alpha=None, require_multiple=0):
48 """Output a fast consensus sequence of the alignment. 49 50 This doesn't do anything fancy at all. It will just go through the 51 sequence residue by residue and count up the number of each type 52 of residue (ie. A or G or T or C for DNA) in all sequences in the 53 alignment. If the percentage of the most common residue type is 54 greater then the passed threshold, then we will add that residue type, 55 otherwise an ambiguous character will be added. 56 57 This could be made a lot fancier (ie. to take a substitution matrix 58 into account), but it just meant for a quick and dirty consensus. 59 60 Arguments: 61 - threshold - The threshold value that is required to add a particular 62 atom. 63 - ambiguous - The ambiguous character to be added when the threshold is 64 not reached. 65 - consensus_alpha - The alphabet to return for the consensus sequence. 66 If this is None, then we will try to guess the alphabet. 67 - require_multiple - If set as 1, this will require that more than 68 1 sequence be part of an alignment to put it in the consensus (ie. 69 not just 1 sequence and gaps). 70 """ 71 # Iddo Friedberg, 1-JUL-2004: changed ambiguous default to "X" 72 consensus = '' 73 74 # find the length of the consensus we are creating 75 con_len = self.alignment.get_alignment_length() 76 77 # go through each seq item 78 for n in range(con_len): 79 # keep track of the counts of the different atoms we get 80 atom_dict = {} 81 num_atoms = 0 82 83 for record in self.alignment._records: 84 # make sure we haven't run past the end of any sequences 85 # if they are of different lengths 86 if n < len(record.seq): 87 if record.seq[n] != '-' and record.seq[n] != '.': 88 if record.seq[n] not in atom_dict: 89 atom_dict[record.seq[n]] = 1 90 else: 91 atom_dict[record.seq[n]] += 1 92 93 num_atoms = num_atoms + 1 94 95 max_atoms = [] 96 max_size = 0 97 98 for atom in atom_dict: 99 if atom_dict[atom] > max_size: 100 max_atoms = [atom] 101 max_size = atom_dict[atom] 102 elif atom_dict[atom] == max_size: 103 max_atoms.append(atom) 104 105 if require_multiple and num_atoms == 1: 106 consensus += ambiguous 107 elif (len(max_atoms) == 1) and ((float(max_size) / float(num_atoms)) 108 >= threshold): 109 consensus += max_atoms[0] 110 else: 111 consensus += ambiguous 112 113 # we need to guess a consensus alphabet if one isn't specified 114 if consensus_alpha is None: 115 consensus_alpha = self._guess_consensus_alphabet(ambiguous) 116 117 return Seq(consensus, consensus_alpha)
118
119 - def gap_consensus(self, threshold=.7, ambiguous="X", 120 consensus_alpha=None, require_multiple=0):
121 """Same as dumb_consensus(), but allows gap on the output. 122 123 Things to do: 124 - Let the user define that with only one gap, the result 125 character in consensus is gap. 126 - Let the user select gap character, now 127 it takes the same as input. 128 """ 129 # Iddo Friedberg, 1-JUL-2004: changed ambiguous default to "X" 130 consensus = '' 131 132 # find the length of the consensus we are creating 133 con_len = self.alignment.get_alignment_length() 134 135 # go through each seq item 136 for n in range(con_len): 137 # keep track of the counts of the different atoms we get 138 atom_dict = {} 139 num_atoms = 0 140 141 for record in self.alignment._records: 142 # make sure we haven't run past the end of any sequences 143 # if they are of different lengths 144 if n < len(record.seq): 145 if record.seq[n] not in atom_dict: 146 atom_dict[record.seq[n]] = 1 147 else: 148 atom_dict[record.seq[n]] += 1 149 150 num_atoms += 1 151 152 max_atoms = [] 153 max_size = 0 154 155 for atom in atom_dict: 156 if atom_dict[atom] > max_size: 157 max_atoms = [atom] 158 max_size = atom_dict[atom] 159 elif atom_dict[atom] == max_size: 160 max_atoms.append(atom) 161 162 if require_multiple and num_atoms == 1: 163 consensus += ambiguous 164 elif (len(max_atoms) == 1) and ((float(max_size) / float(num_atoms)) 165 >= threshold): 166 consensus += max_atoms[0] 167 else: 168 consensus += ambiguous 169 170 # we need to guess a consensus alphabet if one isn't specified 171 if consensus_alpha is None: 172 # TODO - Should we make this into a Gapped alphabet? 173 consensus_alpha = self._guess_consensus_alphabet(ambiguous) 174 175 return Seq(consensus, consensus_alpha)
176
177 - def _guess_consensus_alphabet(self, ambiguous):
178 """Pick an (ungapped) alphabet for an alignment consesus sequence. 179 180 This just looks at the sequences we have, checks their type, and 181 returns as appropriate type which seems to make sense with the 182 sequences we've got. 183 """ 184 # Start with the (un-gapped version of) the alignment alphabet 185 a = Alphabet._get_base_alphabet(self.alignment._alphabet) 186 187 # Now check its compatible with all the rest of the sequences 188 for record in self.alignment: 189 # Get the (un-gapped version of) the sequence's alphabet 190 alt = Alphabet._get_base_alphabet(record.seq.alphabet) 191 if not isinstance(alt, a.__class__): 192 raise ValueError("Alignment contains a sequence with \ 193 an incompatible alphabet.") 194 195 # Check the ambiguous character we are going to use in the consensus 196 # is in the alphabet's list of valid letters (if defined). 197 if hasattr(a, "letters") and a.letters is not None \ 198 and ambiguous not in a.letters: 199 # We'll need to pick a more generic alphabet... 200 if isinstance(a, IUPAC.IUPACUnambiguousDNA): 201 if ambiguous in IUPAC.IUPACUnambiguousDNA().letters: 202 a = IUPAC.IUPACUnambiguousDNA() 203 else: 204 a = Alphabet.generic_dna 205 elif isinstance(a, IUPAC.IUPACUnambiguousRNA): 206 if ambiguous in IUPAC.IUPACUnambiguousRNA().letters: 207 a = IUPAC.IUPACUnambiguousRNA() 208 else: 209 a = Alphabet.generic_rna 210 elif isinstance(a, IUPAC.IUPACProtein): 211 if ambiguous in IUPAC.ExtendedIUPACProtein().letters: 212 a = IUPAC.ExtendedIUPACProtein() 213 else: 214 a = Alphabet.generic_protein 215 else: 216 a = Alphabet.single_letter_alphabet 217 return a
218
219 - def replacement_dictionary(self, skip_chars=[]):
220 """Generate a replacement dictionary to plug into a substitution matrix 221 222 This should look at an alignment, and be able to generate the number 223 of substitutions of different residues for each other in the 224 aligned object. 225 226 Will then return a dictionary with this information:: 227 228 {('A', 'C') : 10, ('C', 'A') : 12, ('G', 'C') : 15 ....} 229 230 This also treats weighted sequences. The following example shows how 231 we calculate the replacement dictionary. Given the following 232 multiple sequence alignment:: 233 234 GTATC 0.5 235 AT--C 0.8 236 CTGTC 1.0 237 238 For the first column we have:: 239 240 ('A', 'G') : 0.5 * 0.8 = 0.4 241 ('C', 'G') : 0.5 * 1.0 = 0.5 242 ('A', 'C') : 0.8 * 1.0 = 0.8 243 244 We then continue this for all of the columns in the alignment, summing 245 the information for each substitution in each column, until we end 246 up with the replacement dictionary. 247 248 Arguments: 249 - skip_chars - A list of characters to skip when creating the dictionary. 250 251 For instance, you might have Xs (screened stuff) or Ns, and not want 252 to include the ambiguity characters in the dictionary. 253 """ 254 # get a starting dictionary based on the alphabet of the alignment 255 rep_dict, skip_items = self._get_base_replacements(skip_chars) 256 257 # iterate through each record 258 for rec_num1 in range(len(self.alignment._records)): 259 # iterate through each record from one beyond the current record 260 # to the end of the list of records 261 for rec_num2 in range(rec_num1 + 1, len(self.alignment._records)): 262 # for each pair of records, compare the sequences and add 263 # the pertinent info to the dictionary 264 rep_dict = self._pair_replacement( 265 self.alignment._records[rec_num1].seq, 266 self.alignment._records[rec_num2].seq, 267 self.alignment._records[rec_num1].annotations.get('weight', 1.0), 268 self.alignment._records[rec_num2].annotations.get('weight', 1.0), 269 rep_dict, skip_items) 270 271 return rep_dict
272
273 - def _pair_replacement(self, seq1, seq2, weight1, weight2, 274 start_dict, ignore_chars):
275 """Compare two sequences and generate info on the replacements seen. 276 277 Arguments: 278 - seq1, seq2 - The two sequences to compare. 279 - weight1, weight2 - The relative weights of seq1 and seq2. 280 - start_dict - The dictionary containing the starting replacement 281 info that we will modify. 282 - ignore_chars - A list of characters to ignore when calculating 283 replacements (ie. '-'). 284 285 Returns: 286 - A replacment dictionary which is modified from initial_dict with 287 the information from the sequence comparison. 288 """ 289 # loop through each residue in the sequences 290 for residue_num in range(len(seq1)): 291 residue1 = seq1[residue_num] 292 try: 293 residue2 = seq2[residue_num] 294 # if seq2 is shorter, then we just stop looking at replacements 295 # and return the information 296 except IndexError: 297 return start_dict 298 299 # if the two residues are characters we want to count 300 if (residue1 not in ignore_chars) and (residue2 not in ignore_chars): 301 try: 302 # add info about the replacement to the dictionary, 303 # modified by the sequence weights 304 start_dict[(residue1, residue2)] += weight1 * weight2 305 306 # if we get a key error, then we've got a problem with alphabets 307 except KeyError: 308 raise ValueError("Residues %s, %s not found in alphabet %s" 309 % (residue1, residue2, 310 self.alignment._alphabet)) 311 312 return start_dict
313
314 - def _get_all_letters(self):
315 """Returns a string containing the expected letters in the alignment.""" 316 all_letters = self.alignment._alphabet.letters 317 if all_letters is None \ 318 or (isinstance(self.alignment._alphabet, Alphabet.Gapped) 319 and all_letters == self.alignment._alphabet.gap_char): 320 # We are dealing with a generic alphabet class where the 321 # letters are not defined! We must build a list of the 322 # letters used... 323 set_letters = set() 324 for record in self.alignment: 325 # Note the built in set does not have a union_update 326 # which was provided by the sets module's Set 327 set_letters = set_letters.union(record.seq) 328 list_letters = sorted(set_letters) 329 all_letters = "".join(list_letters) 330 return all_letters
331
332 - def _get_base_replacements(self, skip_items=[]):
333 """Get a zeroed dictionary of all possible letter combinations. 334 335 This looks at the type of alphabet and gets the letters for it. 336 It then creates a dictionary with all possible combinations of these 337 letters as keys (ie. ('A', 'G')) and sets the values as zero. 338 339 Returns: 340 - The base dictionary created 341 - A list of alphabet items to skip when filling the dictionary. 342 (Right now the only thing I can imagine in this list is gap 343 characters, but maybe X's or something else might be useful later. 344 This will also include any characters that are specified to be 345 skipped.) 346 """ 347 base_dictionary = {} 348 all_letters = self._get_all_letters() 349 350 # if we have a gapped alphabet we need to find the gap character 351 # and drop it out 352 if isinstance(self.alignment._alphabet, Alphabet.Gapped): 353 skip_items.append(self.alignment._alphabet.gap_char) 354 all_letters = all_letters.replace(self.alignment._alphabet.gap_char, '') 355 356 # now create the dictionary 357 for first_letter in all_letters: 358 for second_letter in all_letters: 359 if first_letter not in skip_items and \ 360 second_letter not in skip_items: 361 base_dictionary[(first_letter, second_letter)] = 0 362 363 return base_dictionary, skip_items
364
365 - def pos_specific_score_matrix(self, axis_seq=None, 366 chars_to_ignore=[]):
367 """Create a position specific score matrix object for the alignment. 368 369 This creates a position specific score matrix (pssm) which is an 370 alternative method to look at a consensus sequence. 371 372 Arguments: 373 - chars_to_ignore - A listing of all characters not to include in 374 the pssm. If the alignment alphabet declares a gap character, 375 then it will be excluded automatically. 376 - axis_seq - An optional argument specifying the sequence to 377 put on the axis of the PSSM. This should be a Seq object. If nothing 378 is specified, the consensus sequence, calculated with default 379 parameters, will be used. 380 381 Returns: 382 - A PSSM (position specific score matrix) object. 383 """ 384 # determine all of the letters we have to deal with 385 all_letters = self._get_all_letters() 386 assert all_letters 387 388 if not isinstance(chars_to_ignore, list): 389 raise TypeError("chars_to_ignore should be a list.") 390 391 # if we have a gap char, add it to stuff to ignore 392 if isinstance(self.alignment._alphabet, Alphabet.Gapped): 393 chars_to_ignore.append(self.alignment._alphabet.gap_char) 394 395 for char in chars_to_ignore: 396 all_letters = all_letters.replace(char, '') 397 398 if axis_seq: 399 left_seq = axis_seq 400 assert len(axis_seq) == self.alignment.get_alignment_length() 401 else: 402 left_seq = self.dumb_consensus() 403 404 pssm_info = [] 405 # now start looping through all of the sequences and getting info 406 for residue_num in range(len(left_seq)): 407 score_dict = self._get_base_letters(all_letters) 408 for record in self.alignment._records: 409 try: 410 this_residue = record.seq[residue_num] 411 # if we hit an index error we've run out of sequence and 412 # should not add new residues 413 except IndexError: 414 this_residue = None 415 416 if this_residue and this_residue not in chars_to_ignore: 417 weight = record.annotations.get('weight', 1.0) 418 try: 419 score_dict[this_residue] += weight 420 # if we get a KeyError then we have an alphabet problem 421 except KeyError: 422 raise ValueError("Residue %s not found in alphabet %s" 423 % (this_residue, 424 self.alignment._alphabet)) 425 426 pssm_info.append((left_seq[residue_num], 427 score_dict)) 428 429 return PSSM(pssm_info)
430
431 - def _get_base_letters(self, letters):
432 """Create a zeroed dictionary with all of the specified letters. 433 """ 434 base_info = {} 435 for letter in letters: 436 base_info[letter] = 0 437 438 return base_info
439
440 - def information_content(self, start=0, 441 end=None, 442 e_freq_table=None, log_base=2, 443 chars_to_ignore=[]):
444 """Calculate the information content for each residue along an alignment. 445 446 Arguments: 447 - start, end - The starting an ending points to calculate the 448 information content. These points should be relative to the first 449 sequence in the alignment, starting at zero (ie. even if the 'real' 450 first position in the seq is 203 in the initial sequence, for 451 the info content, we need to use zero). This defaults to the entire 452 length of the first sequence. 453 - e_freq_table - A FreqTable object specifying the expected frequencies 454 for each letter in the alphabet we are using (e.g. {'G' : 0.4, 455 'C' : 0.4, 'T' : 0.1, 'A' : 0.1}). Gap characters should not be 456 included, since these should not have expected frequencies. 457 - log_base - The base of the logathrim to use in calculating the 458 information content. This defaults to 2 so the info is in bits. 459 - chars_to_ignore - A listing of characterw which should be ignored 460 in calculating the info content. 461 462 Returns: 463 - A number representing the info content for the specified region. 464 465 Please see the Biopython manual for more information on how information 466 content is calculated. 467 """ 468 # if no end was specified, then we default to the end of the sequence 469 if end is None: 470 end = len(self.alignment._records[0].seq) 471 472 if start < 0 or end > len(self.alignment._records[0].seq): 473 raise ValueError("Start (%s) and end (%s) are not in the \ 474 range %s to %s" 475 % (start, end, 0, len(self.alignment._records[0].seq))) 476 # determine random expected frequencies, if necessary 477 random_expected = None 478 if not e_freq_table: 479 # TODO - What about ambiguous alphabets? 480 base_alpha = Alphabet._get_base_alphabet(self.alignment._alphabet) 481 if isinstance(base_alpha, Alphabet.ProteinAlphabet): 482 random_expected = Protein20Random 483 elif isinstance(base_alpha, Alphabet.NucleotideAlphabet): 484 random_expected = Nucleotide4Random 485 else: 486 errstr = "Error in alphabet: not Nucleotide or Protein, " 487 errstr += "supply expected frequencies" 488 raise ValueError(errstr) 489 del base_alpha 490 elif not isinstance(e_freq_table, FreqTable.FreqTable): 491 raise ValueError("e_freq_table should be a FreqTable object") 492 493 # determine all of the letters we have to deal with 494 all_letters = self._get_all_letters() 495 for char in chars_to_ignore: 496 all_letters = all_letters.replace(char, '') 497 498 info_content = {} 499 for residue_num in range(start, end): 500 freq_dict = self._get_letter_freqs(residue_num, 501 self.alignment._records, 502 all_letters, chars_to_ignore) 503 # print freq_dict, 504 column_score = self._get_column_info_content(freq_dict, 505 e_freq_table, 506 log_base, 507 random_expected) 508 509 info_content[residue_num] = column_score 510 # sum up the score 511 total_info = sum(info_content.values()) 512 # fill in the ic_vector member: holds IC for each column 513 for i in info_content: 514 self.ic_vector[i] = info_content[i] 515 return total_info
516
517 - def _get_letter_freqs(self, residue_num, all_records, letters, to_ignore):
518 """Determine the frequency of specific letters in the alignment. 519 520 Arguments: 521 - residue_num - The number of the column we are getting frequencies 522 from. 523 - all_records - All of the SeqRecords in the alignment. 524 - letters - The letters we are interested in getting the frequency 525 for. 526 - to_ignore - Letters we are specifically supposed to ignore. 527 528 This will calculate the frequencies of each of the specified letters 529 in the alignment at the given frequency, and return this as a 530 dictionary where the keys are the letters and the values are the 531 frequencies. 532 """ 533 freq_info = self._get_base_letters(letters) 534 535 total_count = 0 536 # collect the count info into the dictionary for all the records 537 for record in all_records: 538 try: 539 if record.seq[residue_num] not in to_ignore: 540 weight = record.annotations.get('weight', 1.0) 541 freq_info[record.seq[residue_num]] += weight 542 total_count += weight 543 # getting a key error means we've got a problem with the alphabet 544 except KeyError: 545 raise ValueError("Residue %s not found in alphabet %s" 546 % (record.seq[residue_num], 547 self.alignment._alphabet)) 548 549 if total_count == 0: 550 # This column must be entirely ignored characters 551 for letter in freq_info: 552 assert freq_info[letter] == 0 553 # TODO - Map this to NA or NaN? 554 else: 555 # now convert the counts into frequencies 556 for letter in freq_info: 557 freq_info[letter] = freq_info[letter] / total_count 558 559 return freq_info
560
561 - def _get_column_info_content(self, obs_freq, e_freq_table, log_base, 562 random_expected):
563 """Calculate the information content for a column. 564 565 Arguments: 566 - obs_freq - The frequencies observed for each letter in the column. 567 - e_freq_table - An optional argument specifying the expected 568 frequencies for each letter. This is a SubsMat.FreqTable instance. 569 - log_base - The base of the logathrim to use in calculating the 570 info content. 571 """ 572 try: 573 gap_char = self.alignment._alphabet.gap_char 574 except AttributeError: 575 # The alphabet doesn't declare a gap - there could be none 576 # in the sequence... or just a vague alphabet. 577 gap_char = "-" # Safe? 578 579 if e_freq_table: 580 if not isinstance(e_freq_table, FreqTable.FreqTable): 581 raise ValueError("e_freq_table should be a FreqTable object") 582 # check the expected freq information to make sure it is good 583 for key in obs_freq: 584 if (key != gap_char and key not in e_freq_table): 585 raise ValueError("Expected frequency letters %s " 586 "do not match observed %s" 587 % (list(e_freq_table), 588 list(obs_freq) - [gap_char])) 589 590 total_info = 0.0 591 592 for letter in obs_freq: 593 inner_log = 0.0 594 # if we have expected frequencies, modify the log value by them 595 # gap characters do not have expected frequencies, so they 596 # should just be the observed frequency. 597 if letter != gap_char: 598 if e_freq_table: 599 inner_log = obs_freq[letter] / e_freq_table[letter] 600 else: 601 inner_log = obs_freq[letter] / random_expected 602 # if the observed frequency is zero, we don't add any info to the 603 # total information content 604 if inner_log > 0: 605 letter_info = (obs_freq[letter] * 606 math.log(inner_log) / math.log(log_base)) 607 total_info += letter_info 608 return total_info
609
610 - def get_column(self, col):
611 # TODO - Deprecate this and implement slicing? 612 return self.alignment[:, col]
613 614
615 -class PSSM(object):
616 """Represent a position specific score matrix. 617 618 This class is meant to make it easy to access the info within a PSSM 619 and also make it easy to print out the information in a nice table. 620 621 Let's say you had an alignment like this:: 622 623 GTATC 624 AT--C 625 CTGTC 626 627 The position specific score matrix (when printed) looks like:: 628 629 G A T C 630 G 1 1 0 1 631 T 0 0 3 0 632 A 1 1 0 0 633 T 0 0 2 0 634 C 0 0 0 3 635 636 You can access a single element of the PSSM using the following:: 637 638 your_pssm[sequence_number][residue_count_name] 639 640 For instance, to get the 'T' residue for the second element in the 641 above alignment you would need to do: 642 643 your_pssm[1]['T'] 644 """
645 - def __init__(self, pssm):
646 """Initialize with pssm data to represent. 647 648 The pssm passed should be a list with the following structure: 649 650 list[0] - The letter of the residue being represented (for instance, 651 from the example above, the first few list[0]s would be GTAT... 652 list[1] - A dictionary with the letter substitutions and counts. 653 """ 654 self.pssm = pssm
655
656 - def __getitem__(self, pos):
657 return self.pssm[pos][1]
658
659 - def __str__(self):
660 out = " " 661 all_residues = sorted(self.pssm[0][1]) 662 663 # first print out the top header 664 for res in all_residues: 665 out += " %s" % res 666 out += "\n" 667 668 # for each item, write out the substitutions 669 for item in self.pssm: 670 out += "%s " % item[0] 671 for res in all_residues: 672 out += " %.1f" % item[1][res] 673 674 out += "\n" 675 return out
676
677 - def get_residue(self, pos):
678 """Return the residue letter at the specified position. 679 """ 680 return self.pssm[pos][0]
681 682 693 694 if __name__ == "__main__": 695 print("Quick test") 696 from Bio import AlignIO 697 from Bio.Align.Generic import Alignment 698 699 filename = "../../Tests/GFF/multi.fna" 700 format = "fasta" 701 expected = FreqTable.FreqTable({"A": 0.25, "G": 0.25, "T": 0.25, "C": 0.25}, 702 FreqTable.FREQ, 703 IUPAC.unambiguous_dna) 704 705 alignment = AlignIO.read(open(filename), format) 706 for record in alignment: 707 print(record.seq) 708 print("=" * alignment.get_alignment_length()) 709 710 summary = SummaryInfo(alignment) 711 consensus = summary.dumb_consensus(ambiguous="N") 712 print(consensus) 713 consensus = summary.gap_consensus(ambiguous="N") 714 print(consensus) 715 print("") 716 print(summary.pos_specific_score_matrix(chars_to_ignore=['-'], 717 axis_seq=consensus)) 718 print("") 719 # Have a generic alphabet, without a declared gap char, so must tell 720 # provide the frequencies and chars to ignore explicitly. 721 print(summary.information_content(e_freq_table=expected, 722 chars_to_ignore=['-'])) 723 print("") 724 print("Trying a protein sequence with gaps and stops") 725 726 alpha = Alphabet.HasStopCodon(Alphabet.Gapped(Alphabet.generic_protein, "-"), "*") 727 a = Alignment(alpha) 728 a.add_sequence("ID001", "MHQAIFIYQIGYP*LKSGYIQSIRSPEYDNW-") 729 a.add_sequence("ID002", "MH--IFIYQIGYAYLKSGYIQSIRSPEY-NW*") 730 a.add_sequence("ID003", "MHQAIFIYQIGYPYLKSGYIQSIRSPEYDNW*") 731 print(a) 732 print("=" * a.get_alignment_length()) 733 734 s = SummaryInfo(a) 735 c = s.dumb_consensus(ambiguous="X") 736 print(c) 737 c = s.gap_consensus(ambiguous="X") 738 print(c) 739 print("") 740 print(s.pos_specific_score_matrix(chars_to_ignore=['-', '*'], axis_seq=c)) 741 742 print(s.information_content(chars_to_ignore=['-', '*'])) 743 744 print("Done") 745