]> git.donarmstrong.com Git - biopieces.git/blob - code_ruby/lib/maasha/seq.rb
fixed wrap bug in seq.rb
[biopieces.git] / code_ruby / lib / maasha / seq.rb
1 # Copyright (C) 2007-2012 Martin A. Hansen.
2
3 # This program is free software; you can redistribute it and/or
4 # modify it under the terms of the GNU General Public License
5 # as published by the Free Software Foundation; either version 2
6 # of the License, or (at your option) any later version.
7
8 # This program is distributed in the hope that it will be useful,
9 # but WITHOUT ANY WARRANTY; without even the implied warranty of
10 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
11 # GNU General Public License for more details.
12
13 # You should have received a copy of the GNU General Public License
14 # along with this program; if not, write to the Free Software
15 # Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
16
17 # http://www.gnu.org/copyleft/gpl.html
18
19 # >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>><<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
20
21 # This software is part of the Biopieces framework (www.biopieces.org).
22
23 # >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>><<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
24
25 require 'maasha/bits'
26 require 'maasha/seq/digest'
27 require 'maasha/seq/trim'
28 require 'narray'
29
30 autoload :BackTrack,   'maasha/seq/backtrack'
31 autoload :Dynamic,     'maasha/seq/dynamic'
32 autoload :Homopolymer, 'maasha/seq/homopolymer'
33 autoload :Hamming,     'maasha/seq/hamming'
34 autoload :Levenshtein, 'maasha/seq/levenshtein'
35 autoload :Ambiguity,   'maasha/seq/ambiguity'
36
37 # Residue alphabets
38 DNA     = %w[a t c g]
39 RNA     = %w[a u c g]
40 PROTEIN = %w[f l s y c w p h q r i m t n k v a d e g]
41 INDELS  = %w[. - _ ~]
42
43 # Translation table 11
44 # (http://www.ncbi.nlm.nih.gov/Taxonomy/taxonomyhome.html/index.cgi?chapter=cgencodes#SG11)
45 #   AAs  = FFLLSSSSYY**CC*WLLLLPPPPHHQQRRRRIIIMTTTTNNKKSSRRVVVVAAAADDEEGGGG
46 # Starts = ---M---------------M------------MMMM---------------M------------
47 # Base1  = TTTTTTTTTTTTTTTTCCCCCCCCCCCCCCCCAAAAAAAAAAAAAAAAGGGGGGGGGGGGGGGG
48 # Base2  = TTTTCCCCAAAAGGGGTTTTCCCCAAAAGGGGTTTTCCCCAAAAGGGGTTTTCCCCAAAAGGGG
49 # Base3  = TCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAGTCAG
50 TRANS_TAB11_START = {
51   "TTG" => "M", "CTG" => "M", "ATT" => "M", "ATC" => "M",
52   "ATA" => "M", "ATG" => "M", "GTG" => "M"
53 }
54
55 TRANS_TAB11 = {
56   "TTT" => "F", "TCT" => "S", "TAT" => "Y", "TGT" => "C",
57   "TTC" => "F", "TCC" => "S", "TAC" => "Y", "TGC" => "C",
58   "TTA" => "L", "TCA" => "S", "TAA" => "*", "TGA" => "*",
59   "TTG" => "L", "TCG" => "S", "TAG" => "*", "TGG" => "W",
60   "CTT" => "L", "CCT" => "P", "CAT" => "H", "CGT" => "R",
61   "CTC" => "L", "CCC" => "P", "CAC" => "H", "CGC" => "R",
62   "CTA" => "L", "CCA" => "P", "CAA" => "Q", "CGA" => "R",
63   "CTG" => "L", "CCG" => "P", "CAG" => "Q", "CGG" => "R",
64   "ATT" => "I", "ACT" => "T", "AAT" => "N", "AGT" => "S",
65   "ATC" => "I", "ACC" => "T", "AAC" => "N", "AGC" => "S",
66   "ATA" => "I", "ACA" => "T", "AAA" => "K", "AGA" => "R",
67   "ATG" => "M", "ACG" => "T", "AAG" => "K", "AGG" => "R",
68   "GTT" => "V", "GCT" => "A", "GAT" => "D", "GGT" => "G",
69   "GTC" => "V", "GCC" => "A", "GAC" => "D", "GGC" => "G",
70   "GTA" => "V", "GCA" => "A", "GAA" => "E", "GGA" => "G",
71   "GTG" => "V", "GCG" => "A", "GAG" => "E", "GGG" => "G"
72 }
73
74 # Error class for all exceptions to do with Seq.
75 class SeqError < StandardError; end
76
77 class Seq
78   # Quality scores bases
79   SCORE_BASE = 33
80   SCORE_MIN  = 0
81   SCORE_MAX  = 40
82
83   include Digest
84   include Trim
85
86   attr_accessor :seq_name, :seq, :type, :qual
87
88   # Class method to instantiate a new Sequence object given
89   # a Biopiece record.
90   def self.new_bp(record)
91     seq_name = record[:SEQ_NAME]
92     seq      = record[:SEQ]
93     type     = record[:SEQ_TYPE].to_sym if record[:SEQ_TYPE]
94     qual     = record[:SCORES]
95
96     self.new(seq_name: seq_name, seq: seq, type: type, qual: qual)
97   end
98
99   # Class method that generates all possible oligos of a specifed length and type.
100   def self.generate_oligos(length, type)
101     raise SeqError, "Cannot generate negative oligo length: #{length}" if length <= 0
102
103     case type.downcase
104     when :dna     then alph = DNA
105     when :rna     then alph = RNA
106     when :protein then alph = PROTEIN
107     else
108       raise SeqError, "Unknown sequence type: #{type}"
109     end
110
111     oligos = [""]
112
113     (1 .. length).each do
114       list = []
115
116       oligos.each do |oligo|
117         alph.each do |char|
118           list << oligo + char
119         end
120       end
121
122       oligos = list
123     end
124
125     oligos
126   end
127
128   # Initialize a sequence object with the following options:
129   # - :seq_name   Name of the sequence.
130   # - :seq        The sequence.
131   # - :type       The sequence type - DNA, RNA, or protein
132   # - :qual       An Illumina type quality scores string.
133   def initialize(options = {})
134     @seq_name = options[:seq_name]
135     @seq      = options[:seq]
136     @type     = options[:type]
137     @qual     = options[:qual]
138
139     if @seq and @qual and @seq.length != @qual.length
140       raise SeqError, "Sequence length and score length mismatch: #{@seq.length} != #{@qual.length}"
141     end
142   end
143
144   # Method that guesses and returns the sequence type
145   # by inspecting the first 100 residues.
146   def type_guess
147     raise SeqError, "Guess failed: sequence is nil" if self.seq.nil?
148
149     case self.seq[0 ... 100].downcase
150     when /[flpqie]/ then return :protein
151     when /[u]/      then return :rna
152     else                 return :dna
153     end
154   end
155
156   # Method that guesses and sets the sequence type
157   # by inspecting the first 100 residues.
158   def type_guess!
159     self.type = self.type_guess
160     self
161   end
162
163   # Returns the length of a sequence.
164   def length
165     self.seq.nil? ? 0 : self.seq.length
166   end
167
168   alias :len :length
169
170   # Return the number indels in a sequence.
171   def indels
172     regex = Regexp.new(/[#{Regexp.escape(INDELS.join(""))}]/)
173     self.seq.scan(regex).size
174   end
175
176   # Method to remove indels from seq and qual if qual.
177   def indels_remove
178     if self.qual.nil?
179       self.seq.delete!(Regexp.escape(INDELS.join('')))
180     else
181       na_seq  = NArray.to_na(self.seq, "byte")
182       na_qual = NArray.to_na(self.qual, "byte")
183       mask    = NArray.byte(self.length)
184
185       INDELS.each do |c|
186         mask += na_seq.eq(c.ord)
187       end
188
189       mask = mask.eq(0)
190
191       self.seq  = na_seq[mask].to_s
192       self.qual = na_qual[mask].to_s
193     end
194
195     self
196   end
197
198   # Method that returns true is a given sequence type is DNA.
199   def is_dna?
200     self.type == :dna
201   end
202
203   # Method that returns true is a given sequence type is RNA.
204   def is_rna?
205     self.type == :rna
206   end
207
208   # Method that returns true is a given sequence type is protein.
209   def is_protein?
210     self.type == :protein
211   end
212
213   # Method to transcribe DNA to RNA.
214   def to_rna
215     raise SeqError, "Cannot transcribe 0 length sequence" if self.length == 0
216     raise SeqError, "Cannot transcribe sequence type: #{self.type}" unless self.is_dna?
217     self.type = :rna
218     self.seq.tr!('Tt','Uu')
219   end
220
221   # Method to reverse-transcribe RNA to DNA.
222   def to_dna
223     raise SeqError, "Cannot reverse-transcribe 0 length sequence" if self.length == 0
224     raise SeqError, "Cannot reverse-transcribe sequence type: #{self.type}" unless self.is_rna?
225     self.type = :dna
226     self.seq.tr!('Uu','Tt')
227   end
228
229   # Method to translate a DNA sequence to protein.
230   def translate!(trans_tab = 11)
231     raise SeqError, "Sequence type must be 'dna' - not #{self.type}" unless self.type == :dna
232     raise SeqError, "Sequence length must be a multiplum of 3 - was: #{self.length}" unless (self.length % 3) == 0
233
234     case trans_tab
235     when 11
236       codon_start_hash = TRANS_TAB11_START
237       codon_hash       = TRANS_TAB11
238     else
239       raise SeqError, "Unknown translation table: #{trans_tab}"
240     end
241
242     codon  = self.seq[0 ... 3].upcase
243
244     aa = codon_start_hash[codon]
245
246     raise SeqError, "Unknown start codon: #{codon}" if aa.nil?
247
248     protein = aa
249
250     i = 3
251
252     while i < self.length
253       codon = self.seq[i ... i + 3].upcase
254
255       aa = codon_hash[codon]
256
257       raise SeqError, "Unknown codon: #{codon}" if aa.nil?
258
259       protein << aa
260
261       i += 3
262     end
263
264     self.seq  = protein
265     self.qual = nil
266     self.type = :protein
267
268     self
269   end
270
271   alias :to_protein! :translate!
272
273   def translate(trans_tab = 11)
274     self.dup.translate!(trans_tab)
275   end
276
277   alias :to_protein :translate
278
279   # Method that given a Seq entry returns a Biopieces record (a hash).
280   def to_bp
281     raise SeqError, "Missing seq_name" if self.seq_name.nil?
282     raise SeqError, "Missing seq"      if self.seq.nil?
283
284     record             = {}
285     record[:SEQ_NAME] = self.seq_name
286     record[:SEQ]      = self.seq
287     record[:SEQ_LEN]  = self.length
288     record[:SCORES]   = self.qual if self.qual
289     record
290   end
291
292   # Method that given a Seq entry returns a FASTA entry (a string).
293   def to_fasta(wrap = nil)
294     raise SeqError, "Missing seq_name" if self.seq_name.nil? or self.seq_name == ''
295     raise SeqError, "Missing seq"      if self.seq.nil?      or self.seq.empty?
296
297     seq_name = self.seq_name
298     seq      = self.seq.dup
299
300     unless wrap.nil?
301       seq.gsub!(/(.{#{wrap}})/) do |match|
302         match << $/
303       end
304
305       seq.chomp!
306     end
307
308     ">" + seq_name + $/ + seq + $/
309   end
310
311   # Method that given a Seq entry returns a FASTQ entry (a string).
312   def to_fastq
313     raise SeqError, "Missing seq_name" if self.seq_name.nil?
314     raise SeqError, "Missing seq"      if self.seq.nil?
315     raise SeqError, "Missing qual"     if self.qual.nil?
316
317     seq_name = self.seq_name.to_s
318     seq      = self.seq.to_s
319     qual     = self.qual.to_s
320
321     "@" + seq_name + $/ + seq + $/ + "+" + $/ + qual + $/
322   end
323
324   # Method that generates a unique key for a
325   # DNA sequence and return this key as a Fixnum.
326   def to_key
327     key = 0
328     
329     self.seq.upcase.each_char do |char|
330       key <<= 2
331       
332       case char
333       when 'A' then key |= 0
334       when 'C' then key |= 1
335       when 'G' then key |= 2
336       when 'T' then key |= 3
337       else raise SeqError, "Bad residue: #{char}"
338       end
339     end
340     
341     key
342   end
343
344   # Method to reverse the sequence.
345   def reverse
346     entry = Seq.new(
347       seq_name: self.seq_name,
348       seq:      self.seq.reverse,
349       type:     self.type,
350       qual:     (self.qual ? self.qual.reverse : self.qual)
351     )
352
353     entry
354   end
355
356   # Method to reverse the sequence.
357   def reverse!
358     self.seq.reverse!
359     self.qual.reverse! if self.qual
360     self
361   end
362
363   # Method that complements sequence including ambiguity codes.
364   def complement
365     raise SeqError, "Cannot complement 0 length sequence" if self.length == 0
366
367     entry = Seq.new(
368       seq_name: self.seq_name,
369       type:     self.type,
370       qual:     self.qual
371     )
372
373     if self.is_dna?
374       entry.seq = self.seq.tr('AGCUTRYWSMKHDVBNagcutrywsmkhdvbn', 'TCGAAYRWSKMDHBVNtcgaayrwskmdhbvn')
375     elsif self.is_rna?
376       entry.seq = self.seq.tr('AGCUTRYWSMKHDVBNagcutrywsmkhdvbn', 'UCGAAYRWSKMDHBVNucgaayrwskmdhbvn')
377     else
378       raise SeqError, "Cannot complement sequence type: #{self.type}"
379     end
380
381     entry
382   end
383
384   # Method that complements sequence including ambiguity codes.
385   def complement!
386     raise SeqError, "Cannot complement 0 length sequence" if self.length == 0
387
388     if self.is_dna?
389       self.seq.tr!('AGCUTRYWSMKHDVBNagcutrywsmkhdvbn', 'TCGAAYRWSKMDHBVNtcgaayrwskmdhbvn')
390     elsif self.is_rna?
391       self.seq.tr!('AGCUTRYWSMKHDVBNagcutrywsmkhdvbn', 'UCGAAYRWSKMDHBVNucgaayrwskmdhbvn')
392     else
393       raise SeqError, "Cannot complement sequence type: #{self.type}"
394     end
395
396     self
397   end
398
399   # Method to determine the Hamming Distance between
400   # two Sequence objects (case insensitive).
401   def hamming_distance(entry, options = nil)
402     if options and options[:ambiguity]
403       Hamming.distance(self.seq, entry.seq)
404     else
405       self.seq.upcase.hamming_distance(entry.seq.upcase)
406     end
407   end
408
409   # Method to determine the Edit Distance between
410   # two Sequence objects (case insensitive).
411   def edit_distance(entry)
412     Levenshtein.distance(self.seq, entry.seq)
413   end
414
415   # Method that generates a random sequence of a given length and type.
416   def generate(length, type)
417     raise SeqError, "Cannot generate sequence length < 1: #{length}" if length <= 0
418
419     case type
420     when :dna     then alph = DNA
421     when :rna     then alph = RNA
422     when :protein then alph = PROTEIN
423     else
424       raise SeqError, "Unknown sequence type: #{type}"
425     end
426
427     seq_new   = Array.new(length) { alph[rand(alph.size)] }.join("")
428     self.seq  = seq_new
429     self.type = type
430     seq_new
431   end
432
433   # Method to return a new Seq object with shuffled sequence.
434   def shuffle
435     Seq.new(
436       seq_name: self.seq_name,
437       seq:      self.seq.split('').shuffle!.join,
438       type:     self.type,
439       qual:     self.qual
440     )
441   end
442
443   # Method to shuffle a sequence randomly inline.
444   def shuffle!
445     self.seq = self.seq.split('').shuffle!.join
446     self
447   end
448
449   # Method to add two Seq objects.
450   def +(entry)
451     new_entry = Seq.new()
452     new_entry.seq  = self.seq  + entry.seq
453     new_entry.type = self.type              if self.type == entry.type
454     new_entry.qual = self.qual + entry.qual if self.qual and entry.qual
455     new_entry
456   end
457
458   # Method to concatenate sequence entries.
459   def <<(entry)
460     raise SeqError, "sequences of different types" unless self.type == entry.type
461     raise SeqError, "qual is missing in one entry" unless self.qual.class == entry.qual.class
462
463     self.seq  << entry.seq
464     self.qual << entry.qual unless entry.qual.nil?
465
466     self
467   end
468
469   # Index method for Seq objects.
470   def [](*args)
471     entry = Seq.new
472     entry.seq_name = self.seq_name
473     entry.seq      = self.seq[*args]
474     entry.type     = self.type
475     entry.qual     = self.qual[*args] unless self.qual.nil?
476
477     entry
478   end
479
480   # Index assignment method for Seq objects.
481   def []=(*args, entry)
482     self.seq[*args]  = entry.seq[*args]
483     self.qual[*args] = entry.qual[*args] unless self.qual.nil?
484
485     self
486   end
487
488   # Method that returns a subsequence of from a given start position
489   # and of a given length.
490   def subseq(start, length = self.length - start)
491     raise SeqError, "subsequence start: #{start} < 0"                                                if start  < 0
492     raise SeqError, "subsequence length: #{length} < 0"                                              if length < 0
493     raise SeqError, "subsequence start + length > Seq.length: #{start} + #{length} > #{self.length}" if start + length > self.length
494
495     if length == 0
496       seq  = ""
497       qual = "" unless self.qual.nil?
498     else
499       stop = start + length - 1
500
501       seq  = self.seq[start .. stop]
502       qual = self.qual[start .. stop] unless self.qual.nil?
503     end
504
505     seq_name = self.seq_name.nil? ? nil : self.seq_name.dup
506
507     Seq.new(seq_name: seq_name, seq: seq, type: self.type, qual: qual)
508   end
509
510   # Method that replaces a sequence with a subsequence from a given start position
511   # and of a given length.
512   def subseq!(start, length = self.length - start)
513     s = subseq(start, length)
514
515     self.seq_name = s.seq_name
516     self.seq      = s.seq
517     self.type     = s.type
518     self.qual     = s.qual
519
520     self
521   end
522
523   # Method that returns a subsequence of a given length
524   # beginning at a random position.
525   def subseq_rand(length)
526     if self.length - length + 1 == 0
527       start = 0
528     else
529       start = rand(self.length - length + 1)
530     end
531
532     self.subseq(start, length)
533   end
534
535   # Method that returns the residue compositions of a sequence in
536   # a hash where the key is the residue and the value is the residue
537   # count.
538   def composition
539     comp = Hash.new(0);
540
541     self.seq.upcase.each_char do |char|
542       comp[char] += 1
543     end
544
545     comp
546   end
547
548   # Method that returns the percentage of hard masked residues
549   # or N's in a sequence.
550   def hard_mask
551     ((self.seq.upcase.scan("N").size.to_f / (self.len - self.indels).to_f) * 100).round(2)
552   end
553
554   # Method that returns the percentage of soft masked residues
555   # or lower cased residues in a sequence.
556   def soft_mask
557     ((self.seq.scan(/[a-z]/).size.to_f / (self.len - self.indels).to_f) * 100).round(2)
558   end
559
560   # Hard masks sequence residues where the corresponding quality score
561   # is below a given cutoff.
562   def mask_seq_hard!(cutoff)
563     raise SeqError, "seq is nil"  if self.seq.nil?
564     raise SeqError, "qual is nil" if self.qual.nil?
565     raise SeqError, "cufoff value: #{cutoff} out of range #{SCORE_MIN} .. #{SCORE_MAX}" unless (SCORE_MIN .. SCORE_MAX).include? cutoff
566
567     na_seq  = NArray.to_na(self.seq, "byte")
568     na_qual = NArray.to_na(self.qual, "byte")
569     mask    = (na_qual - SCORE_BASE) < cutoff
570     mask   *= na_seq.ne("-".ord)
571
572     na_seq[mask] = 'N'.ord
573
574     self.seq = na_seq.to_s
575
576     self
577   end
578
579   # Soft masks sequence residues where the corresponding quality score
580   # is below a given cutoff.
581   def mask_seq_soft!(cutoff)
582     raise SeqError, "seq is nil"  if self.seq.nil?
583     raise SeqError, "qual is nil" if self.qual.nil?
584     raise SeqError, "cufoff value: #{cutoff} out of range #{SCORE_MIN} .. #{SCORE_MAX}" unless (SCORE_MIN .. SCORE_MAX).include? cutoff
585
586     na_seq  = NArray.to_na(self.seq, "byte")
587     na_qual = NArray.to_na(self.qual, "byte")
588     mask    = (na_qual - SCORE_BASE) < cutoff
589     mask   *= na_seq.ne("-".ord)
590
591     na_seq[mask] ^= ' '.ord
592
593     self.seq = na_seq.to_s
594
595     self
596   end
597
598   # Method that determines if a quality score string can be
599   # absolutely identified as base 33.
600   def qual_base33?
601     self.qual.match(/[!-:]/) ? true : false
602   end
603  
604   # Method that determines if a quality score string may be base 64.
605   def qual_base64?
606     self.qual.match(/[K-h]/) ? true : false
607   end
608
609   # Method to determine if a quality score is valid accepting only 0-40 range.
610   def qual_valid?(encoding)
611     raise SeqError, "Missing qual" if self.qual.nil?
612
613     case encoding
614     when :base_33 then return true if self.qual.match(/^[!-I]*$/)
615     when :base_64 then return true if self.qual.match(/^[@-h]*$/)
616     else raise SeqError, "unknown quality score encoding: #{encoding}"
617     end
618
619     false
620   end
621
622   # Method to coerce quality scores to be within the 0-40 range.
623   def qual_coerce!(encoding)
624     raise SeqError, "Missing qual" if self.qual.nil?
625
626     case encoding
627     when :base_33 then self.qual.tr!("[J-~]", "I")
628     when :base_64 then self.qual.tr!("[i-~]", "h")
629     else raise SeqError, "unknown quality score encoding: #{encoding}"
630     end 
631
632     self
633   end
634
635   # Method to convert quality scores.
636   def qual_convert!(from, to)
637     raise SeqError, "unknown quality score encoding: #{from}" unless from == :base_33 or from == :base_64
638     raise SeqError, "unknown quality score encoding: #{to}"   unless to   == :base_33 or to   == :base_64
639
640     if from == :base_33 and to == :base_64
641       na_qual   = NArray.to_na(self.qual, "byte")
642       na_qual  += 64 - 33
643       self.qual = na_qual.to_s
644     elsif from == :base_64 and to == :base_33
645       self.qual.tr!("[;-?]", "@")  # Handle negative Solexa values from -5 to -1 (set these to 0).
646       na_qual   = NArray.to_na(self.qual, "byte")
647       na_qual  -= 64 - 33
648       self.qual = na_qual.to_s
649     end
650
651     self
652   end
653
654   # Method to calculate and return the mean quality score.
655   def scores_mean
656     raise SeqError, "Missing qual in entry" if self.qual.nil?
657
658     na_qual = NArray.to_na(self.qual, "byte")
659     na_qual -= SCORE_BASE
660
661     na_qual.mean
662   end
663
664   # Method to find open reading frames (ORFs).
665   def each_orf(size_min, size_max, start_codons, stop_codons, pick_longest = false)
666     orfs    = []
667     pos_beg = 0
668
669     regex_start = Regexp.new(start_codons.join('|'), true)
670     regex_stop  = Regexp.new(stop_codons.join('|'), true)
671
672     while pos_beg = self.seq.index(regex_start, pos_beg)
673       pos_end = pos_beg + 3
674        
675       while pos_end = self.seq.index(regex_stop, pos_end)
676         length = (pos_end - pos_beg) + 3
677
678         if (length % 3) == 0
679           if size_min <= length and length <= size_max
680             subseq = self.subseq(pos_beg, length)
681
682             orfs << [subseq, pos_beg, pos_end + 3]
683           end
684
685           break
686         end
687
688         pos_end += 1
689       end
690
691       pos_beg += 1
692     end
693
694     if pick_longest
695       orf_hash = {}
696
697       orfs.each { |orf| orf_hash[orf.last] = orf unless orf_hash[orf.last] }
698
699       orfs = orf_hash.values
700     end
701
702     if block_given?
703       orfs.each { |orf| yield orf }
704     else
705       return orfs
706     end
707   end
708 end
709
710 __END__