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