Package Bio :: Package SeqIO :: Module AceIO
[hide private]
[frames] | no frames]

Source Code for Module Bio.SeqIO.AceIO

  1  # Copyright 2008-2010 by Peter Cock.  All rights reserved. 
  2  # 
  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   
  7  """Bio.SeqIO support for the "ace" file format. 
  8   
  9  You are expected to use this module via the Bio.SeqIO functions. 
 10  See also the Bio.Sequencing.Ace module which offers more than just accessing 
 11  the contig consensus sequences in an ACE file as SeqRecord objects. 
 12  """ 
 13   
 14  from __future__ import print_function 
 15   
 16  from Bio.Seq import Seq 
 17  from Bio.SeqRecord import SeqRecord 
 18  from Bio.Alphabet import generic_nucleotide, generic_dna, generic_rna, Gapped 
 19  from Bio.Sequencing import Ace 
 20   
 21   
22 -def AceIterator(handle):
23 """Returns SeqRecord objects from an ACE file. 24 25 This uses the Bio.Sequencing.Ace module to do the hard work. Note that 26 by iterating over the file in a single pass, we are forced to ignore any 27 WA, CT, RT or WR footer tags. 28 29 Ace files include the base quality for each position, which are taken 30 to be PHRED style scores. Just as if you had read in a FASTQ or QUAL file 31 using PHRED scores using Bio.SeqIO, these are stored in the SeqRecord's 32 letter_annotations dictionary under the "phred_quality" key. 33 34 >>> from Bio import SeqIO 35 >>> with open("Ace/consed_sample.ace", "rU") as handle: 36 ... for record in SeqIO.parse(handle, "ace"): 37 ... print("%s %s... %i" % (record.id, record.seq[:10], len(record))) 38 ... print(max(record.letter_annotations["phred_quality"])) 39 Contig1 agccccgggc... 1475 40 90 41 42 However, ACE files do not include a base quality for any gaps in the 43 consensus sequence, and these are represented in Biopython with a quality 44 of zero. Using zero is perhaps misleading as there may be very strong 45 evidence to support the gap in the consensus. Previous versions of 46 Biopython therefore used None instead, but this complicated usage, and 47 prevented output of the gapped sequence as FASTQ format. 48 49 >>> from Bio import SeqIO 50 >>> with open("Ace/contig1.ace", "rU") as handle: 51 ... for record in SeqIO.parse(handle, "ace"): 52 ... print("%s ...%s..." % (record.id, record.seq[85:95])) 53 ... print(record.letter_annotations["phred_quality"][85:95]) 54 ... print(max(record.letter_annotations["phred_quality"])) 55 Contig1 ...AGAGG-ATGC... 56 [57, 57, 54, 57, 57, 0, 57, 72, 72, 72] 57 90 58 Contig2 ...GAATTACTAT... 59 [68, 68, 68, 68, 68, 68, 68, 68, 68, 68] 60 90 61 62 """ 63 for ace_contig in Ace.parse(handle): 64 #Convert the ACE contig record into a SeqRecord... 65 consensus_seq_str = ace_contig.sequence 66 #Assume its DNA unless there is a U in it, 67 if "U" in consensus_seq_str: 68 if "T" in consensus_seq_str: 69 #Very odd! Error? 70 alpha = generic_nucleotide 71 else: 72 alpha = generic_rna 73 else: 74 alpha = generic_dna 75 76 if "*" in consensus_seq_str: 77 #For consistency with most other file formats, map 78 #any * gaps into - gaps. 79 assert "-" not in consensus_seq_str 80 consensus_seq = Seq(consensus_seq_str.replace("*", "-"), 81 Gapped(alpha, gap_char="-")) 82 else: 83 consensus_seq = Seq(consensus_seq_str, alpha) 84 85 #TODO? - Base segments (BS lines) which indicates which read 86 #phrap has chosen to be the consensus at a particular position. 87 #Perhaps as SeqFeature objects? 88 89 #TODO - Supporting reads (RD lines, plus perhaps QA and DS lines) 90 #Perhaps as SeqFeature objects? 91 92 seq_record = SeqRecord(consensus_seq, 93 id=ace_contig.name, 94 name=ace_contig.name) 95 96 #Consensus base quality (BQ lines). Note that any gaps (originally 97 #as * characters) in the consensus do not get a quality entry, so 98 #we assign a quality of None (zero would be missleading as there may 99 #be excelent support for having a gap here). 100 quals = [] 101 i = 0 102 for base in consensus_seq: 103 if base == "-": 104 quals.append(0) 105 else: 106 quals.append(ace_contig.quality[i]) 107 i += 1 108 assert i == len(ace_contig.quality) 109 seq_record.letter_annotations["phred_quality"] = quals 110 111 yield seq_record
112 #All done 113 114 115 if __name__ == "__main__": 116 from Bio._utils import run_doctest 117 run_doctest() 118