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