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