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

Source Code for Module Bio.SeqIO.PhdIO

  1  # Copyright 2008-2010 by Peter Cock.  All rights reserved. 
  2  # Revisions copyright 2009 by Cymon J. Cox.  All rights reserved. 
  3  # 
  4  # This code is part of the Biopython distribution and governed by its 
  5  # license.  Please see the LICENSE file that should have been included 
  6  # as part of this package. 
  7   
  8  """Bio.SeqIO support for the "phd" file format. 
  9   
 10  PHD files are output by PHRED and used by PHRAP and CONSED. 
 11   
 12  You are expected to use this module via the Bio.SeqIO functions, under the 
 13  format name "phd". See also the underlying Bio.Sequencing.Phd module. 
 14   
 15  For example, using Bio.SeqIO we can read in one of the example PHRED files 
 16  from the Biopython unit tests: 
 17   
 18      >>> from Bio import SeqIO 
 19      >>> for record in SeqIO.parse("Phd/phd1", "phd"): 
 20      ...     print(record.id) 
 21      ...     print("%s..." % record.seq[:10]) 
 22      ...     print("%s..." % record.letter_annotations["phred_quality"][:10]) 
 23      34_222_(80-A03-19).b.ab1 
 24      ctccgtcgga... 
 25      [9, 9, 10, 19, 22, 37, 28, 28, 24, 22]... 
 26      425_103_(81-A03-19).g.ab1 
 27      cgggatccca... 
 28      [14, 17, 22, 10, 10, 10, 15, 8, 8, 9]... 
 29      425_7_(71-A03-19).b.ab1 
 30      acataaatca... 
 31      [10, 10, 10, 10, 8, 8, 6, 6, 6, 6]... 
 32   
 33  Since PHRED files contain quality scores, you can save them as FASTQ or as 
 34  QUAL files, for example using Bio.SeqIO.write(...), or simply with the format 
 35  method of the SeqRecord object: 
 36   
 37      >>> print(record[:50].format("fastq")) 
 38      @425_7_(71-A03-19).b.ab1 
 39      acataaatcaaattactnaccaacacacaaaccngtctcgcgtagtggag 
 40      + 
 41      ++++))'''')(''')$!$''')''''(+.''$!$))))+)))''''''' 
 42      <BLANKLINE> 
 43   
 44  Or, 
 45   
 46      >>> print(record[:50].format("qual")) 
 47      >425_7_(71-A03-19).b.ab1 
 48      10 10 10 10 8 8 6 6 6 6 8 7 6 6 6 8 3 0 3 6 6 6 8 6 6 6 6 7 
 49      10 13 6 6 3 0 3 8 8 8 8 10 8 8 8 6 6 6 6 6 6 6 
 50      <BLANKLINE> 
 51   
 52  Note these examples only show the first 50 bases to keep the output short. 
 53  """ 
 54   
 55  from __future__ import print_function 
 56   
 57  from Bio.SeqRecord import SeqRecord 
 58  from Bio.Sequencing import Phd 
 59  from Bio.SeqIO.Interfaces import SequentialSequenceWriter 
 60  from Bio.SeqIO import QualityIO 
 61   
 62  __docformat__ = "restructuredtext en" 
 63   
 64   
65 -def PhdIterator(handle):
66 """Returns SeqRecord objects from a PHD file. 67 68 This uses the Bio.Sequencing.Phd module to do the hard work. 69 """ 70 phd_records = Phd.parse(handle) 71 for phd_record in phd_records: 72 # Convert the PHY record into a SeqRecord... 73 # The "filename" can contain spaces, e.g. 'HWI-EAS94_4_1_1_602_99 1' 74 # from unit test example file phd_solexa. 75 # This will cause problems if used as the record identifier 76 # (e.g. output for FASTQ format). 77 name = phd_record.file_name.split(None, 1)[0] 78 seq_record = SeqRecord(phd_record.seq, 79 id=name, name=name, 80 description=phd_record.file_name) 81 # Just re-use the comments dictionary as the SeqRecord's annotations 82 seq_record.annotations = phd_record.comments 83 # And store the qualities and peak locations as per-letter-annotation 84 seq_record.letter_annotations["phred_quality"] = \ 85 [int(site[1]) for site in phd_record.sites] 86 try: 87 seq_record.letter_annotations["peak_location"] = \ 88 [int(site[2]) for site in phd_record.sites] 89 except IndexError: 90 # peak locations are not always there according to 91 # David Gordon (the Consed author) 92 pass 93 yield seq_record
94 # All done 95 96
97 -class PhdWriter(SequentialSequenceWriter):
98 """Class to write Phd format files""" 99
100 - def __init__(self, handle):
102
103 - def write_record(self, record):
104 """Write a single Phd record to the file.""" 105 assert record.seq, "No sequence present in SeqRecord" 106 # This method returns the 'phred_quality' scores or converted 107 # 'solexa_quality' scores if present, else raises a value error 108 phred_qualities = QualityIO._get_phred_quality(record) 109 peak_locations = record.letter_annotations.get("peak_location", None) 110 assert len(record.seq) == len(phred_qualities), "Number of " + \ 111 "phd quality scores does not match length of sequence" 112 if peak_locations: 113 assert len(record.seq) == len(peak_locations), "Number " + \ 114 "of peak location scores does not match length of sequence" 115 if None in phred_qualities: 116 raise ValueError("A quality value of None was found") 117 if record.description.startswith("%s " % record.id): 118 title = record.description 119 else: 120 title = "%s %s" % (record.id, record.description) 121 self.handle.write("BEGIN_SEQUENCE %s\nBEGIN_COMMENT\n" 122 % self.clean(title)) 123 for annot in [k.lower() for k in Phd.CKEYWORDS]: 124 value = None 125 if annot == "trim": 126 if record.annotations.get("trim", None): 127 value = "%s %s %.4f" % record.annotations["trim"] 128 elif annot == "trace_peak_area_ratio": 129 if record.annotations.get("trace_peak_area_ratio", None): 130 value = "%.4f" % record.annotations[ 131 "trace_peak_area_ratio"] 132 else: 133 value = record.annotations.get(annot, None) 134 if value or value == 0: 135 self.handle.write("%s: %s\n" % (annot.upper(), value)) 136 137 self.handle.write("END_COMMENT\nBEGIN_DNA\n") 138 for i, site in enumerate(record.seq): 139 if peak_locations: 140 self.handle.write("%s %i %i\n" % ( 141 site, 142 round(phred_qualities[i]), 143 peak_locations[i]) 144 ) 145 else: 146 self.handle.write("%s %i\n" % ( 147 site, 148 round(phred_qualities[i])) 149 ) 150 151 self.handle.write("END_DNA\nEND_SEQUENCE\n")
152 153 154 if __name__ == "__main__": 155 from Bio._utils import run_doctest 156 run_doctest() 157