Class: Bio::DB::Sam
- Inherits:
-
Object
- Object
- Bio::DB::Sam
- Defined in:
- lib/bio/db/sam.rb
Instance Attribute Summary collapse
-
#sam_file ⇒ Object
readonly
Returns the value of attribute sam_file.
Class Method Summary collapse
-
.finalize(id) ⇒ Object
Destructor method that closes the file before letting the object be garbage collected.
-
.merge(files, merged_file, headers, add_RG, by_qname) ⇒ Object
Merges n BAM files.
Instance Method Summary collapse
-
#average_coverage(chromosome, qstart, len) ⇒ Object
Returns the average coverage of a region in a bam file.
-
#chromosome_coverage(chromosome, qstart, len) ⇒ Object
Returns an array with the coverage at each possition in the queried region This is a simple average coverage just calculated with the first and last possition of the alignment, ignoring the gaps.
-
#close ⇒ Object
Closes the sam file and destroys the C pointers using the functions provided by libbam.
-
#deprecated_pileup(cmd) ⇒ Object
utility method that does not use the samtools API, it calls samtools directly as if on the command line and catches the output, to use this method you must have a version of samtools that supports the pileup command (< 0.1.17) otherwise the command will fail.
-
#each_reference ⇒ Object
yields each reference name and its length.
-
#fetch(chromosome, qstart, qend) ⇒ Object
Returns an array of Alignments on a given region.
-
#fetch_reference(chromosome, qstart, qend) ⇒ Object
Returns the sequence for a given region.
-
#fetch_with_function(chromosome, qstart, qend, function) ⇒ Object
Executes a function on each Alignment inside the queried region of the chromosome.
- #index_stats ⇒ Object
-
#initialize(optsa = {}) ⇒ Sam
constructor
To make a new sam object.
-
#load_index ⇒ Object
Loads the bam index to be used for fetching.
-
#load_reference ⇒ Object
Loads the reference file to be able to query regions of it.
-
#mpileup(opts = {}) ⇒ Object
calls the mpileup function, opts is a hash of options identical to the command line options for mpileup.
-
#mpileup_plus(opts) ⇒ Object
experimental method that spawns a samtools mpileup | bcftools view process and supports returning of pileup vcf otherwise works like mpileup.
-
#open ⇒ Object
Function that actually opens the sam file Throws a SAMException if the file can’t be open.
-
#query_string(chromosome, qstart, qend) ⇒ Object
Generates a query sting to be used by the region parser in samtools.
-
#to_s ⇒ Object
Prints a description of the sam file in a text format containg if it is binary or text, the path and the fasta file of the reference.
Constructor Details
#initialize(optsa = {}) ⇒ Sam
To make a new sam object. Initialize expects a hash optsa with the following elemets:
- fasta
-
The fasta file with the reference. (nil)
- bam
-
path to a binary SAM file (nil)
- tam
-
path to a text SAM file (nil)
- compressed
-
If the binary file is compressed (true)
- write
-
If the file is to be writen (false). Not supported yet.
NOTE: you can’t use binary and text formats simultaneusly. To make queries, the file has to be a sorted binary. This function doesn’t actually open the file, it just prepares the object to be opened in a later stage.
29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 |
# File 'lib/bio/db/sam.rb', line 29 def initialize(optsa={}) opts = { :fasta => nil, :bam => nil,:tam => nil, :compressed => true, :write => false }.merge!(optsa) @fasta_path = opts[:fasta] @compressed = opts[:compressed] @write = opts[:write] bam = opts[:bam] tam = opts[:tam] if bam == nil && tam == nil && @fasta_path == nil then raise SAMException.new(), "No alignment or reference file" elsif bam != nil && tam != nil then raise SAMException.new(), "Alignment has to be in either text or binary format, not both" elsif bam != nil then @binary = true @sam = bam elsif tam != nil then @sam = tam @binary = false end @fasta_file = nil @sam_file = nil ObjectSpace.define_finalizer(self, self.class.method(:finalize).to_proc) end |
Instance Attribute Details
#sam_file ⇒ Object (readonly)
Returns the value of attribute sam_file.
18 19 20 |
# File 'lib/bio/db/sam.rb', line 18 def sam_file @sam_file end |
Class Method Details
.finalize(id) ⇒ Object
Destructor method that closes the file before letting the object be garbage collected.
104 105 106 107 |
# File 'lib/bio/db/sam.rb', line 104 def Sam.finalize(id) id.close() puts "Finalizing #{id} at #{Time.new}" end |
.merge(files, merged_file, headers, add_RG, by_qname) ⇒ Object
Merges n BAM files. This doesn’t require to create a SAM object
- files
-
An array with the paths to the files.
- merged_file
-
The path to the merged file
- headers
-
The BAM file containing the header
- add_RG
-
If true, the RG tag is added (infered from the filenames)
- by_qname
-
If true, the bamfiles should by ordered by query name, if false, by coordinates.
259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 |
# File 'lib/bio/db/sam.rb', line 259 def self.merge(files, merged_file, headers, add_RG, by_qname) strptrs = [] strptrs << FFI::MemoryPointer.from_string("merge") files.each do |file| strptrs << FFI::MemoryPointer.from_string(file) end strptrs << nil # Now load all the pointers into a native memory block argv = FFI::MemoryPointer.new(:pointer, strptrs.length) strptrs.each_with_index do |p, i| argv[i].put_pointer(0, p) end #void bam_merge_core(int by_qname, const char *out, const char *headers, int n, char * const *fn, int add_RG) Bio::DB::SAM::Tools.bam_merge_core(by_qname, merged_file, headers, strptrs.length, argv, add_RG) end |
Instance Method Details
#average_coverage(chromosome, qstart, len) ⇒ Object
Returns the average coverage of a region in a bam file.
143 144 145 146 147 148 149 150 151 152 153 154 155 |
# File 'lib/bio/db/sam.rb', line 143 def average_coverage(chromosome, qstart, len) #reference = fetch_reference(chromosome, qstart,len) # len = reference.length if len > reference.length coverages = chromosome_coverage(chromosome, qstart, len) total = 0 len.times{ |i| total= total + coverages[i] } avg_cov = total.to_f / len #LibC.free reference avg_cov end |
#chromosome_coverage(chromosome, qstart, len) ⇒ Object
Returns an array with the coverage at each possition in the queried region This is a simple average coverage just calculated with the first and last possition of the alignment, ignoring the gaps.
160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 |
# File 'lib/bio/db/sam.rb', line 160 def chromosome_coverage(chromosome, qstart, len) # reference = fetch_reference(chromosome, qstart,len) # len = reference.length if len > reference.length #p qend.to_s + "-" + qstart.to_s + "framesize " + (qend - qstart).to_s coverages = Array.new(len, 0) chr_cov_proc = Proc.new do |alignment| #last = qstart + len #first = qstart #last = alignment.calend if last > alignment.calend #first = alignment.pos if first < alignment.pos # p first last = alignment.calend - qstart first = alignment.pos - qstart if last < first tmp = last last = first first = last end # STDERR.puts "#{first} #{last}\n" first.upto(last-1) { |i| coverages[i-1] = 1 + coverages[i-1] if i-1 < len && i > 0 } end fetch_with_function(chromosome, qstart, qstart+len, chr_cov_proc) #p coverages coverages end |
#close ⇒ Object
Closes the sam file and destroys the C pointers using the functions provided by libbam
95 96 97 98 99 100 101 |
# File 'lib/bio/db/sam.rb', line 95 def close() Bio::DB::SAM::Tools.fai_destroy(@fasta_index) unless @fasta_index.nil? || @fasta_index.null? Bio::DB::SAM::Tools.bam_index_destroy(@sam_index) unless @sam_index.nil? || @sam_index.null? Bio::DB::SAM::Tools.samclose(@sam_file) unless @sam_file.nil? @sam_file = nil @fasta_index = nil end |
#deprecated_pileup(cmd) ⇒ Object
utility method that does not use the samtools API, it calls samtools directly as if on the command line and catches the output, to use this method you must have a version of samtools that supports the pileup command (< 0.1.17) otherwise the command will fail. mpileup is the preferred method for getting pileups. With this method the sam object should be created as usual, but you need to pass this method a string of options for samtools you don’t need to provide the call to samtools pileup itself or -f <fasta file> or the bam file itself, these are taken from the sam object
473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 |
# File 'lib/bio/db/sam.rb', line 473 def deprecated_pileup( cmd ) system('samtools pileup > /dev/null 2>&1') ##assumes samtools is in the path... if $?.exitstatus > 1 raise RuntimeError, "samtools is required on the path. A version of samtools with the pileup function is required" end raise SAMException.new(), "No BAMFile provided" unless @sam and @binary raise SAMException.new(), "No FastA provided" unless @fasta_path command = 'samtools pileup ' + cmd + " -f #{@fasta_path}" + " #{@sam}" pipe = IO.popen(command) while line = pipe.gets yield Pileup.new(line) end pipe.close end |
#each_reference ⇒ Object
yields each reference name and its length
535 536 537 538 539 540 |
# File 'lib/bio/db/sam.rb', line 535 def each_reference refs = index_stats refs.each_pair do |k, v| yield k, v[:length] end end |
#fetch(chromosome, qstart, qend) ⇒ Object
Returns an array of Alignments on a given region.
211 212 213 214 215 216 217 218 219 |
# File 'lib/bio/db/sam.rb', line 211 def fetch(chromosome, qstart, qend) als = Array.new fetchAlignment = Proc.new do |alignment| als.push(alignment.clone) 0 end fetch_with_function(chromosome, qstart, qend, fetchAlignment) als end |
#fetch_reference(chromosome, qstart, qend) ⇒ Object
Returns the sequence for a given region.
193 194 195 196 197 198 199 200 201 |
# File 'lib/bio/db/sam.rb', line 193 def fetch_reference(chromosome, qstart,qend) load_reference if @fasta_index.nil? || @fasta_index.null? query = query_string(chromosome, qstart,qend) len = FFI::MemoryPointer.new :int reference = Bio::DB::SAM::Tools.fai_fetch(@fasta_index, query, len) raise SAMException.new(), "Unable to get sequence for reference: "+query if reference.nil? reference end |
#fetch_with_function(chromosome, qstart, qend, function) ⇒ Object
Executes a function on each Alignment inside the queried region of the chromosome. The chromosome can be either the textual name or a FixNum with the internal index. However, you need to get the internal index with the provided API, otherwise the pointer is outside the scope of the C library. Returns the count of alignments in the region. WARNING: Accepts an index already parsed by the library. It fails when you use your own FixNum (FFI-bug?)
226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 |
# File 'lib/bio/db/sam.rb', line 226 def fetch_with_function(chromosome, qstart, qend, function) load_index if @sam_index.nil? || @sam_index.null? chr = FFI::MemoryPointer.new :int beg = FFI::MemoryPointer.new :int last = FFI::MemoryPointer.new :int query = query_string(chromosome, qstart,qend) qpointer = FFI::MemoryPointer.from_string(query) header = @sam_file[:header] Bio::DB::SAM::Tools.bam_parse_region(header,qpointer, chr, beg, last) #raise SAMException.new(), "invalid query: " + query if(chr.read_int < 0) count = 0; fetchAlignment = Proc.new do |bam_alignment, data| alignment = Alignment.new alignment.set(bam_alignment, header) function.call(alignment) count = count + 1 0 end Bio::DB::SAM::Tools.bam_fetch(@sam_file[:x][:bam], @sam_index,chr.read_int,beg.read_int, last.read_int, nil, fetchAlignment) #LibC.free chr #LibC.free beg #LibC.free last #LibC.free qpointer count end |
#index_stats ⇒ Object
494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 |
# File 'lib/bio/db/sam.rb', line 494 def index_stats raise SAMException.new(), "No BAMFile provided" unless @sam and @binary raise SAMException.new(), "No FastA provided" unless @fasta_path strptrs = [] strptrs << FFI::MemoryPointer.from_string("idxstats") strptrs << FFI::MemoryPointer.from_string(@sam) strptrs << nil # Now load all the pointers into a native memory block argv = FFI::MemoryPointer.new(:pointer, strptrs.length) strptrs.each_with_index do |p, i| argv[i].put_pointer(0, p) end index_stats = {} old_stdout = STDOUT.clone read_pipe, write_pipe = IO.pipe() STDOUT.reopen(write_pipe) #int bam_idxstats(int argc, char *argv[]) Bio::DB::SAM::Tools.bam_idxstats(strptrs.length - 1,argv) if fork write_pipe.close STDOUT.reopen(old_stdout) #beware .. stdout from other processes eg tests calling this method can get mixed in... begin while line = read_pipe.readline #TAB delimited with each line consisting of reference sequence name, sequence length, # mapped reads and # unmapped reads. info = line.split(/\t/) next unless info.length == 4 index_stats[ info[0] ] = {:length => info[1].to_i, :mapped_reads => info[2].to_i, :unmapped_reads => info[3].to_i } end rescue EOFError read_pipe.close Process.wait end end #fork index_stats end |
#load_index ⇒ Object
Loads the bam index to be used for fetching. If the index doesn’t exists the index is built provided that the user has writing access to the folder where the BAM file is located. If the creation of the file fails a SAMException is thrown. If the index doesn’t exist, loading it will take more time. It is suggested to generate the index separatedly if the bam file sits on a server where the executing user may not have writing permissions in the server.
114 115 116 117 118 119 120 121 122 123 |
# File 'lib/bio/db/sam.rb', line 114 def load_index() raise SAMException.new(), "Indexes are only supported by BAM files, please use samtools to convert your SAM file" unless @binary @sam_index = Bio::DB::SAM::Tools.bam_index_load(@sam) if @sam_index.null? then p "Generating index for: " + @sam Bio::DB::SAM::Tools.bam_index_build(@sam) @sam_index = Bio::DB::SAM::Tools.bam_index_load(@sam) raise SAMException.new(), "Unable to generate bam index for: " + @sam if @sam_index.nil? || @sam_index.null? end end |
#load_reference ⇒ Object
Loads the reference file to be able to query regions of it. This requires the fai index to exist in the same folder than the reference. If it doesn’t exisits, this functions attempts to generate it. If user doesn’t have writing permissions on the folder, or the creation of the fai fails for any reason, a SAMException is thrown.
128 129 130 131 132 133 134 135 136 137 138 139 140 |
# File 'lib/bio/db/sam.rb', line 128 def load_reference() raise SAMException.new(), "No path for the refernce fasta file. " if @fasta_path.nil? @fasta_index = Bio::DB::SAM::Tools.fai_load(@fasta_path) if @fasta_index.null? then p "Generating index for: " + @fasta_path Bio::DB::SAM::Tools.fai_build(@fasta_path) @fasta_index = Bio::DB::SAM::Tools.fai_load(@fasta_path) raise SAMException.new(), "Unable to generate fasta index for: " + @fasta_path if @fasta_index.nil? || @fasta_index.null? end end |
#mpileup(opts = {}) ⇒ Object
calls the mpileup function, opts is a hash of options identical to the command line options for mpileup. is an iterator that yields a Pileup object for each postion the command line options that generate/affect BCF/VCF are ignored ie (g,u,e,h,I,L,o,p) call the option as a symbol of the flag, eg -r for region is called :r => “some SAM compatible region” eg bam.mpileup(:r => “chr1:1000-2000”, :q => 50) gets the bases with quality > 50 on chr1 between 1000-5000
281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 |
# File 'lib/bio/db/sam.rb', line 281 def mpileup( opts={}) raise SAMException.new(), "No BAMFile provided" unless @sam and @binary raise SAMException.new(), "No FastA provided" unless @fasta_path #long option form to short samtools form.. long_opts = { :region => :r, :illumina_quals => :six, :count_anomalous => :A, :no_baq => :B, :adjust_mapq => :C, :max_per_bam_depth => :d, :extended_baq => :E, :exclude_reads_file => :G, :list_of_positions => :l, :mapping_quality_cap => :M, :ignore_rg => :R, :min_mapping_quality => :q, :min_base_quality => :Q } ##convert any long_opts to short opts temp_opts = opts.dup opts.each_pair do |k,v| if long_opts[k] temp_opts[long_opts[k]] = v temp_opts.delete(k) end end opts = temp_opts ##remove any calls to -g or -u for mpileup, bcf output is not yet supported ##and also associated output options [:g, :u, :e, :h, :I, :L, :o, :p].each {|x| opts.delete(x) } sam_opts = [] #strptrs << FFI::MemoryPointer.from_string("mpileup") opts.each do |k,v| next unless opts[k] ##dont bother unless the values provided are true.. k = '6' if k == :six k = '-' + k.to_s sam_opts << k #strptrs << FFI::MemoryPointer.from_string(k) sam_opts << v.to_s unless ["-R", "-B", "-E", "-6", "-A"].include?(k) #these are just flags so don't pass a value... strptrs << FFI::MemoryPointer.from_string(v.to_s) end sam_opts = sam_opts + ['-f', @fasta_path, @sam] sam_command = "#{File.join(File.(File.dirname(__FILE__)),'sam','external','samtools')} mpileup #{sam_opts.join(' ')} 2> /dev/null" sam_pipe = IO.popen(sam_command) while line = sam_pipe.gets yield Bio::DB::Pileup.new(line) end sam_pipe.close #strptrs << FFI::MemoryPointer.from_string('-f') #strptrs << FFI::MemoryPointer.from_string(@fasta_path) #strptrs << FFI::MemoryPointer.from_string(@sam) #strptrs << nil # Now load all the pointers into a native memory block #argv = FFI::MemoryPointer.new(:pointer, strptrs.length) #strptrs.each_with_index do |p, i| # argv[i].put_pointer(0, p) #end #old_stdout = STDOUT.clone #read_pipe, write_pipe = IO.pipe() #STDOUT.reopen(write_pipe) #int bam_mpileup(int argc, char *argv[]) # Bio::DB::SAM::Tools.bam_mpileup(strptrs.length - 1,argv) #if fork # write_pipe.close # STDOUT.reopen(old_stdout) #beware .. stdout from other processes eg tests calling this method can get mixed in... # begin # while line = read_pipe.readline # yield Pileup.new(line) # end # rescue EOFError # read_pipe.close # Process.wait # end #end end |
#mpileup_plus(opts) ⇒ Object
experimental method that spawns a samtools mpileup | bcftools view process and supports returning of pileup vcf otherwise works like mpileup
363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 |
# File 'lib/bio/db/sam.rb', line 363 def mpileup_plus( opts ) raise SAMException.new(), "No BAMFile provided" unless @sam and @binary raise SAMException.new(), "No FastA provided" unless @fasta_path #long option form to short samtools form.. long_opts = { :region => :r, :illumina_quals => :six, :count_anomalous => :A, :no_baq => :B, :adjust_mapq => :C, :max_per_bam_depth => :d, :extended_baq => :E, :exclude_reads_file => :G, :list_of_positions => :l, :mapping_quality_cap => :M, :ignore_rg => :R, :min_mapping_quality => :q, :min_base_quality => :Q, ###following options are for the -g -u option :genotype_calling => :g, :uncompressed_bcf => :u, :extension_sequencing_probability => :e, :homopolymer_error_coefficient => :h, :no_indels => :I, :skip_indel_over_average_depth => :L, :gap_open_sequencing_error_probability => :o, :platforms => :P } ##convert any long_opts to short opts temp_opts = opts.dup opts.each_pair do |k,v| if long_opts[k] temp_opts[long_opts[k]] = v temp_opts.delete(k) end end opts = temp_opts ##remove any calls to -g or -u for mpileup, bcf output is not yet supported ##and also associated output options #[:g, :u, :e, :h, :I, :L, :o, :p].each {|x| opts.delete(x) } opts[:u] = true if opts[:g] #so that we always get uncompressed output opts.delete(:g) sam_opts = [] #strptrs << FFI::MemoryPointer.from_string("mpileup") opts.each do |k,v| next unless opts[k] ##dont bother unless the values provided are true.. k = '6' if k == :six k = '-' + k.to_s sam_opts << k #strptrs << FFI::MemoryPointer.from_string(k) sam_opts << v.to_s unless ["-R", "-B", "-E", "-6", "-A", "-g", "-u", "-I"].include?(k) #these are just flags so don't pass a value... strptrs << FFI::MemoryPointer.from_string(v.to_s) end sam_opts = sam_opts + ['-f', @fasta_path, @sam] command = "#{File.join(File.(File.dirname(__FILE__)),'sam','external','samtools')} mpileup #{sam_opts.join(' ')} 2> /dev/null" if opts[:u] command = command + " | #{File.join(File.(File.dirname(__FILE__)),'sam','external','bcftools')} view -cg -" end pipe = IO.popen(command) $stderr.puts command if opts[:u] while line = pipe.gets next if line[0,1] == '#' #skip any header or meta-lines, we dont do anything with those yield Bio::DB::Vcf.new(line) end else while line = pipe.gets yield Bio::DB::Pileup.new(line) end end pipe.close #strptrs << FFI::MemoryPointer.from_string('-f') #strptrs << FFI::MemoryPointer.from_string(@fasta_path) #strptrs << FFI::MemoryPointer.from_string(@sam) #strptrs << nil # Now load all the pointers into a native memory block #argv = FFI::MemoryPointer.new(:pointer, strptrs.length) #strptrs.each_with_index do |p, i| # argv[i].put_pointer(0, p) #end #old_stdout = STDOUT.clone #read_pipe, write_pipe = IO.pipe() #STDOUT.reopen(write_pipe) #int bam_mpileup(int argc, char *argv[]) # Bio::DB::SAM::Tools.bam_mpileup(strptrs.length - 1,argv) #if fork # write_pipe.close # STDOUT.reopen(old_stdout) #beware .. stdout from other processes eg tests calling this method can get mixed in... # begin # while line = read_pipe.readline # yield Pileup.new(line) # end # rescue EOFError # read_pipe.close # Process.wait # end #end end |
#open ⇒ Object
Function that actually opens the sam file Throws a SAMException if the file can’t be open.
60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 |
# File 'lib/bio/db/sam.rb', line 60 def open() raise SAMException.new(), "Writing not supported yet" if @write raise SAMException.new(), "No SAM file specified" unless @sam opts = @write ? "w" : "r" if @binary then opts += "b" if @write then unless @compressed then opts += "u" end end end valid = ["r", "w", "wh", "rb", "wb" , "wbu"] unless valid.include?(opts) then raise SAMException.new(), "Invalid options for samopen: " + opts end samFile = Bio::DB::SAM::Tools.samopen(@sam, opts, nil) if samFile.null? then @sam_file = nil raise SAMException.new(), "File not opened: " + @sam end @sam_file = Bio::DB::SAM::Tools::SamfileT.new(samFile) end |
#query_string(chromosome, qstart, qend) ⇒ Object
Generates a query sting to be used by the region parser in samtools. In principle, you shouldn’t need to use this function.
205 206 207 208 |
# File 'lib/bio/db/sam.rb', line 205 def query_string(chromosome, qstart,qend) query = chromosome + ":" + qstart.to_s + "-" + qend.to_s query end |
#to_s ⇒ Object
Prints a description of the sam file in a text format containg if it is binary or text, the path and the fasta file of the reference
90 91 92 |
# File 'lib/bio/db/sam.rb', line 90 def to_s() (@binary ? "Binary" : "Text") + " file: " + @sam + " with fasta: " + @fasta_path end |