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