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