]> git.donarmstrong.com Git - rsem.git/blob - rsem-calculate-expression
Use PATH env to instead the variable
[rsem.git] / rsem-calculate-expression
1 #!/usr/bin/env perl
2
3 use Getopt::Long;
4 use Pod::Usage;
5
6 use FindBin;
7 use lib $FindBin::RealBin;
8 use rsem_perl_utils qw(runCommand collectResults showVersionInfo);
9
10 use Env qw(@PATH);
11 push(@PATH, $FindBin::RealBin, "$FindBin::RealBin/sam");
12
13 use strict;
14
15 #const
16 my $BURNIN = 200;
17 my $NCV = 1000;
18 my $SAMPLEGAP = 1;
19 my $CONFIDENCE = 0.95;
20 my $NSPC = 50;
21
22 my $NMB = 1024; # default
23 my $SortMem = "1G"; # default as 1G per thread
24
25 my $status = 0;
26
27 my $read_type = 1; # default, single end with qual
28
29 my $bowtie_path = "";
30 my $C = 2;
31 my $E = 99999999;
32 my $L = 25;
33 my $maxHits = 200;
34 my $chunkMbs = 0;       # 0 = use bowtie default
35 my $phred33 = 0;
36 my $phred64 = 0;
37 my $solexa = 0;
38
39 my $is_sam = 0;
40 my $is_bam = 0;
41 my $fn_list = "";
42 my $tagName = "XM";
43
44 my $probF = 0.5;
45
46 my $minL = 1;
47 my $maxL = 1000;
48 my $mean = -1;
49 my $sd = 0;
50
51 my $estRSPD = 0;
52 my $B = 20;
53
54 my $nThreads = 1;
55 my $genBamF = 1;  # default is generating transcript bam file
56 my $genGenomeBamF = 0;
57 my $sampling = 0;
58 my $calcPME = 0;
59 my $calcCI = 0;
60 my $var_opt = 0; # temporarily, only for internal use
61 my $quiet = 0;
62 my $help = 0;
63
64 my $paired_end = 0;
65 my $no_qual = 0;
66 my $keep_intermediate_files = 0;
67
68 my $strand_specific = 0;
69
70 my $bowtie2 = 0;
71 my $bowtie2_path = "";
72 my $bowtie2_mismatch_rate = 0.1;
73 my $bowtie2_k = 200;
74 my $bowtie2_sensitivity_level = "sensitive"; # must be one of "very_fast", "fast", "sensitive", "very_sensitive"
75
76 my $version = 0;
77
78 my $mTime = 0;
79 my ($time_start, $time_end, $time_alignment, $time_rsem, $time_ci) = (0, 0, 0, 0, 0);
80
81 my $mate1_list = "";
82 my $mate2_list = "";
83 my $inpF = "";
84
85 my ($refName, $sampleName, $sampleToken, $temp_dir, $stat_dir, $imdName, $statName) = ();
86 my $gap = 32;
87
88 my $alleleS = 0;
89
90 GetOptions("keep-intermediate-files" => \$keep_intermediate_files,
91            "temporary-folder=s" => \$temp_dir,
92            "no-qualities" => \$no_qual,
93            "paired-end" => \$paired_end,
94            "strand-specific" => \$strand_specific,
95            "sam" => \$is_sam,
96            "bam" => \$is_bam,
97            "sam-header-info=s" => \$fn_list,
98            "tag=s" => \$tagName,
99            "seed-length=i" => \$L,
100            "bowtie-path=s" => \$bowtie_path,
101            "bowtie-n=i" => \$C,
102            "bowtie-e=i" => \$E,
103            "bowtie-m=i" => \$maxHits,
104            "bowtie-chunkmbs=i" => \$chunkMbs,
105            "phred33-quals" => \$phred33,
106            "phred64-quals" => \$phred64, #solexa1.3-quals" => \$phred64,
107            "solexa-quals" => \$solexa,
108            "bowtie2" => \$bowtie2,
109            "bowtie2-path=s" => \$bowtie2_path,
110            "bowtie2-mismatch-rate=f" => \$bowtie2_mismatch_rate,
111            "bowtie2-k=i" => \$bowtie2_k,
112            "bowtie2-sensitivity-level=s" => \$bowtie2_sensitivity_level,
113            "forward-prob=f" => \$probF,
114            "fragment-length-min=i" => \$minL,
115            "fragment-length-max=i" => \$maxL,
116            "fragment-length-mean=f" => \$mean,
117            "fragment-length-sd=f" => \$sd,
118            "estimate-rspd" => \$estRSPD,
119            "num-rspd-bins=i" => \$B,
120            "p|num-threads=i" => \$nThreads,
121            "no-bam-output" => sub { $genBamF = 0; },
122            "output-genome-bam" => \$genGenomeBamF,
123            "sampling-for-bam" => \$sampling,
124            "calc-pme" => \$calcPME,
125            "var" => \$var_opt,
126            "calc-ci" => \$calcCI,
127            "ci-memory=i" => \$NMB,
128            "samtools-sort-mem=s" => \$SortMem,
129            "time" => \$mTime,
130            "version" => \$version,
131            "q|quiet" => \$quiet,
132            "h|help" => \$help) or pod2usage(-exitval => 2, -verbose => 2);
133
134 pod2usage(-verbose => 2) if ($help == 1);
135 &showVersionInfo($FindBin::RealBin) if ($version == 1);
136
137 #check parameters and options
138
139 if ($is_sam || $is_bam) {
140     pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (scalar(@ARGV) != 3);
141     pod2usage(-msg => "--sam and --bam cannot be active at the same time!", -exitval => 2, -verbose => 2) if ($is_sam == 1&& $is_bam == 1);
142     pod2usage(-msg => "--bowtie-path, --bowtie-n, --bowtie-e, --bowtie-m, --phred33-quals, --phred64-quals, --solexa-quals, --bowtie2, --bowtie2-path, --bowtie2-mismatch-rate, --bowtie2-k and --bowtie2-sensitivity-level cannot be set if input is SAM/BAM format!", -exitval => 2, -verbose => 2) if ($bowtie_path ne "" || $C != 2 || $E != 99999999 || $maxHits != 200 || $phred33 || $phred64 || $solexa || $bowtie2 || $bowtie2_path ne "" || $bowtie2_mismatch_rate != 0.1 || $bowtie2_k != 200 || $bowtie2_sensitivity_level ne "sensitive");
143 }
144 else {
145     pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (!$paired_end && scalar(@ARGV) != 3 || $paired_end && scalar(@ARGV) != 4);    
146     pod2usage(-msg => "If --no-qualities is set, neither --phred33-quals, --phred64-quals or --solexa-quals can be active!", -exitval => 2, -verbose => 2) if ($no_qual && ($phred33 + $phred64 + $solexa > 0));
147     pod2usage(-msg => "Only one of --phred33-quals, --phred64-quals, and --solexa-quals can be active!", -exitval => 2, -verbose => 2) if ($phred33 + $phred64 + $solexa > 1);    
148     pod2usage(-msg => "--sam , --bam or --sam-header-info cannot be set if use bowtie/bowtie2 aligner to produce alignments!", -exitval => 2, -verbose => 2) if ($is_sam || $is_bam || $fn_list ne "");
149     pod2usage(-msg => "--bowtie2-path, --bowtie2-mismatch-rate, --bowtie2-k and --bowtie2-sensitivity-level cannot be set if bowtie aligner is used!", -exitval => 2, -verbose => 2) if (!$bowtie2 && ($bowtie2_path ne "" || $bowtie2_mismatch_rate != 0.1 || $bowtie2_k != 200 || $bowtie2_sensitivity_level ne "sensitive"));
150     pod2usage(-msg => "--bowtie-path, --bowtie-n, --bowtie-e, --bowtie-m cannot be set if bowtie2 aligner is used!", -exitval => 2, -verbose => 2) if ($bowtie2 && ($bowtie_path ne "" || $C != 2 || $E != 99999999 || $maxHits != 200));
151     pod2usage(-msg => "Mismatch rate must be within [0, 1]!", -exitval => 2, -verbose => 2) if ($bowtie2 && ($bowtie2_mismatch_rate < 0.0 || $bowtie2_mismatch_rate > 1.0));
152     pod2usage(-msg => "Sensitivity level must be one of \"very_fast\", \"fast\", \"sensitive\", and \"very_sensitive\"!", -exitval => 2, -verbose => 2) if ($bowtie2 && (($bowtie2_sensitivity_level ne "very_fast") && ($bowtie2_sensitivity_level ne "fast") && ($bowtie2_sensitivity_level ne "sensitive") && ($bowtie2_sensitivity_level ne "very_sensitive"))); 
153 }
154
155 pod2usage(-msg => "Forward probability should be in [0, 1]!", -exitval => 2, -verbose => 2) if ($probF < 0 || $probF > 1);
156 pod2usage(-msg => "Min fragment length should be at least 1!", -exitval => 2, -verbose => 2) if ($minL < 1);
157 pod2usage(-msg => "Min fragment length should be smaller or equal to max fragment length!", -exitval => 2, -verbose => 2) if ($minL > $maxL);
158 pod2usage(-msg => "The memory allocated for calculating credibility intervals should be at least 1 MB!\n", -exitval => 2, -verbose => 2) if ($NMB < 1);
159 pod2usage(-msg => "Number of threads should be at least 1!\n", -exitval => 2, -verbose => 2) if ($nThreads < 1);
160 pod2usage(-msg => "Seed length should be at least 5!\n", -exitval => 2, -verbose => 2) if ($L < 5);
161 pod2usage(-msg => "--sampling-for-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($sampling && !$genBamF);
162 pod2usage(-msg => "--output-genome-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($genGenomeBamF && !$genBamF);
163
164 if ($L < 25) { print "Warning: the seed length set is less than 25! This is only allowed if the references are not added poly(A) tails.\n"; }
165
166 if ($strand_specific) { $probF = 1.0; }
167
168 if ($paired_end) {
169     if ($no_qual) { $read_type = 2; }
170     else { $read_type = 3; }
171 }
172 else {
173     if ($no_qual) { $read_type = 0; }
174     else { $read_type = 1; }
175 }
176
177 if (scalar(@ARGV) == 3) {
178     if ($is_sam || $is_bam) { $inpF = $ARGV[0]; } 
179     else {$mate1_list = $ARGV[0]; }
180     $refName = $ARGV[1];
181     $sampleName = $ARGV[2];
182 }
183 else {
184     $mate1_list = $ARGV[0];
185     $mate2_list = $ARGV[1];
186     $refName = $ARGV[2];
187     $sampleName = $ARGV[3];
188 }
189
190 if (((-e "$refName.ta") && !(-e "$refName.gt")) || (!(-e "$refName.ta") && (-e "$refName.gt"))) {
191     print "Allele-specific expression related reference files are corrupted!\n";
192     exit(-1);
193 }
194
195 $alleleS = (-e "$refName.ta") && (-e "$refName.gt");
196
197 if ($genGenomeBamF) {
198     open(INPUT, "$refName.ti");
199     my $line = <INPUT>; chomp($line);
200     close(INPUT);
201     my ($M, $type) = split(/ /, $line);
202     pod2usage(-msg => "No genome information provided, so genome bam file cannot be generated!\n", -exitval => 2, -verbose => 2) if ($type != 0);
203 }
204
205 my $pos = rindex($sampleName, '/');
206 if ($pos < 0) { $sampleToken = $sampleName; }
207 else { $sampleToken = substr($sampleName, $pos + 1); }
208
209 if ($temp_dir eq "") { $temp_dir = "$sampleName.temp"; }
210 $stat_dir = "$sampleName.stat";
211
212 if (!(-d $temp_dir) && !mkdir($temp_dir)) { print "Fail to create folder $temp_dir.\n"; exit(-1); }
213 if (!(-d $stat_dir) && !mkdir($stat_dir)) { print "Fail to create folder $stat_dir.\n"; exit(-1); }
214
215 $imdName = "$temp_dir/$sampleToken";
216 $statName = "$stat_dir/$sampleToken";
217
218 if (!$is_sam && !$is_bam && !$no_qual && ($phred33 + $phred64 + $solexa == 0)) { $phred33 = 1; }
219
220 my ($mate_minL, $mate_maxL) = (1, $maxL);
221
222 if ($bowtie_path ne "") { $bowtie_path .= "/"; }
223 if ($bowtie2_path ne "") { $bowtie2_path .= "/"; }
224
225 my $command = "";
226
227 if (!$is_sam && !$is_bam) {
228     if (!$bowtie2) {
229         $command = $bowtie_path."bowtie";
230         if ($no_qual) { $command .= " -f"; }
231         else { $command .= " -q"; }
232     
233         if ($phred33) { $command .= " --phred33-quals"; }
234         elsif ($phred64) { $command .= " --phred64-quals"; }
235         elsif ($solexa) { $command .= " --solexa-quals"; }
236     
237         $command .= " -n $C -e $E -l $L";
238         if ($read_type == 2 || $read_type == 3) { $command .= " -I $minL -X $maxL"; }
239         if ($chunkMbs > 0) { $command .= " --chunkmbs $chunkMbs"; }
240         
241         if ($strand_specific || $probF == 1.0) { $command .= " --norc"; }
242         elsif ($probF == 0.0) { $command .= " --nofw"; }
243         
244         $command .= " -p $nThreads -a -m $maxHits -S";
245         if ($quiet) { $command .= " --quiet"; }    
246
247         $command .= " $refName";
248         if ($read_type == 0 || $read_type == 1) {
249             $command .= " $mate1_list"; 
250         }
251         else {
252             $command .= " -1 $mate1_list -2 $mate2_list";
253         }
254
255         # pipe to samtools to generate a BAM file
256         $command .= " | samtools view -S -b -o $imdName.bam -";
257     }
258     else {
259         $command = $bowtie2_path."bowtie2";
260         if ($no_qual) { $command .= " -f"; }
261         else { $command .= " -q"; }
262
263         if ($phred33) { $command .= " --phred33"; }
264         elsif ($phred64) { $command .= " --phred64"; }
265         elsif ($solexa) { $command .= " --solexa-quals"; }
266
267         if ($bowtie2_sensitivity_level eq "very_fast") { $command .= " --very-fast"; }
268         elsif ($bowtie2_sensitivity_level eq "fast") { $command .= " --fast"; }
269         elsif ($bowtie2_sensitivity_level eq "sensitive") { $command .= " --sensitive"; }
270         else { $command .= " --very-sensitive"; }
271
272         $command .= " --dpad 0 --gbar 99999999 --mp 1,1 --np 1 --score-min L,0,-$bowtie2_mismatch_rate";  
273
274         if ($read_type == 2 || $read_type == 3) { $command .= " -I $minL -X $maxL --no-mixed --no-discordant"; }
275
276         if ($strand_specific || $probF == 1.0) { $command .= " --norc"; }
277         elsif ($probF == 0.0) { $command .= " --nofw"; }
278         
279         $command .= " -p $nThreads -k $bowtie2_k";
280         if ($quiet) { $command .= " --quiet"; }    
281
282         $command .= " -x $refName";
283         if ($read_type == 0 || $read_type == 1) {
284             $command .= " -U $mate1_list"; 
285         }
286         else {
287             $command .= " -1 $mate1_list -2 $mate2_list";
288         }
289
290         # pipe to samtools to generate a BAM file
291         $command .= " | samtools view -S -b -o $imdName.bam -";
292     }
293
294     if ($mTime) { $time_start = time(); }
295
296     &runCommand($command);
297
298     if ($mTime) { $time_end = time(); $time_alignment = $time_end - $time_start; }
299
300     $inpF = "$imdName.bam";
301     $is_bam = 1; # alignments are outputed as a BAM file
302 }
303
304 if ($mTime) { $time_start = time(); }
305
306 $command = "rsem-parse-alignments $refName $imdName $statName";
307
308 my $samInpType;
309 if ($is_sam) { $samInpType = "s"; } 
310 elsif ($is_bam) { $samInpType = "b"; }
311
312 $command .= " $samInpType $inpF -t $read_type";
313 if ($fn_list ne "") { $command .= " -l $fn_list"; }
314 if ($tagName ne "") { $command .= " -tag $tagName"; }
315 if ($quiet) { $command .= " -q"; }
316
317 &runCommand($command);
318
319 $command = "rsem-build-read-index $gap"; 
320 if ($read_type == 0) { $command .= " 0 $quiet $imdName\_alignable.fa"; }
321 elsif ($read_type == 1) { $command .= " 1 $quiet $imdName\_alignable.fq"; }
322 elsif ($read_type == 2) { $command .= " 0 $quiet $imdName\_alignable_1.fa $imdName\_alignable_2.fa"; }
323 elsif ($read_type == 3) { $command .= " 1 $quiet $imdName\_alignable_1.fq $imdName\_alignable_2.fq"; }
324 else { print "Impossible! read_type is not in [1,2,3,4]!\n"; exit(-1); }
325 &runCommand($command);
326
327 my $doesOpen = open(OUTPUT, ">$imdName.mparams");
328 if ($doesOpen == 0) { print "Cannot generate $imdName.mparams!\n"; exit(-1); }
329 print OUTPUT "$minL $maxL\n";
330 print OUTPUT "$probF\n";
331 print OUTPUT "$estRSPD\n";
332 print OUTPUT "$B\n";
333 print OUTPUT "$mate_minL $mate_maxL\n";
334 print OUTPUT "$mean $sd\n";
335 print OUTPUT "$L\n";
336 close(OUTPUT);  
337
338 $command = "rsem-run-em $refName $read_type $sampleName $imdName $statName -p $nThreads";
339 if ($genBamF) { 
340     $command .= " -b $samInpType $inpF";
341     if ($fn_list ne "") { $command .= " 1 $fn_list"; }
342     else { $command .= " 0"; }
343     if ($sampling) { $command .= " --sampling"; }
344 }
345 if ($calcPME || $var_opt || $calcCI) { $command .= " --gibbs-out"; }
346 if ($quiet) { $command .= " -q"; }
347
348 &runCommand($command);
349
350 if ($alleleS) {
351     &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level
352     &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
353     &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
354
355 else {
356     &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
357     &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
358 }
359
360 if ($genBamF) {
361     $command = "samtools sort -@ $nThreads -m $SortMem $sampleName.transcript.bam $sampleName.transcript.sorted";
362     &runCommand($command);
363     $command = "samtools index $sampleName.transcript.sorted.bam";
364     &runCommand($command);
365
366     if ($genGenomeBamF) {
367         $command = "rsem-tbam2gbam $refName $sampleName.transcript.bam $sampleName.genome.bam";
368         &runCommand($command);
369         $command = "samtools sort -@ $nThreads -m $SortMem $sampleName.genome.bam $sampleName.genome.sorted";
370         &runCommand($command);
371         $command = "samtools index $sampleName.genome.sorted.bam";
372         &runCommand($command);
373     }
374 }
375
376 if ($mTime) { $time_end = time(); $time_rsem = $time_end - $time_start; }
377
378 if ($mTime) { $time_start = time(); }
379
380 if ($calcPME || $var_opt || $calcCI ) {
381     $command = "rsem-run-gibbs $refName $imdName $statName $BURNIN $NCV $SAMPLEGAP";
382     $command .= " -p $nThreads";
383     if ($var_opt) { $command .= " --var"; }
384     if ($quiet) { $command .= " -q"; }
385     &runCommand($command);
386 }
387
388 if ($calcPME || $calcCI) {
389     if ($alleleS) {
390         system("mv $sampleName.alleles.results $imdName.alleles.results.bak1");
391         system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak1");
392         system("mv $sampleName.genes.results $imdName.genes.results.bak1");
393         &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level
394         &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
395         &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
396     }
397     else {
398         system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak1");
399         system("mv $sampleName.genes.results $imdName.genes.results.bak1");
400         &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
401         &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
402     }
403 }
404
405 if ($calcCI) {
406     $command = "rsem-calculate-credibility-intervals $refName $imdName $statName $CONFIDENCE $NCV $NSPC $NMB";
407     $command .= " -p $nThreads";
408     if ($quiet) { $command .= " -q"; }
409     &runCommand($command);
410
411     if ($alleleS) {
412         system("mv $sampleName.alleles.results $imdName.alleles.results.bak2");
413         system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak2");
414         system("mv $sampleName.genes.results $imdName.genes.results.bak2");
415         &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level
416         &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
417         &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
418     }
419     else {
420         system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak2");
421         system("mv $sampleName.genes.results $imdName.genes.results.bak2");
422         &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
423         &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
424     }
425 }
426
427 if ($mTime) { $time_end = time(); $time_ci = $time_end - $time_start; }
428
429 if ($mTime) { $time_start = time(); }
430
431 if (!$keep_intermediate_files) {
432     &runCommand("rm -rf $temp_dir", "Fail to delete the temporary folder!");
433 }
434
435 if ($mTime) { $time_end = time(); }
436
437 if ($mTime) { 
438     open(OUTPUT, ">$sampleName.time");
439     print OUTPUT "Aligning reads: $time_alignment s.\n";
440     print OUTPUT "Estimating expression levels: $time_rsem s.\n";
441     print OUTPUT "Calculating credibility intervals: $time_ci s.\n";
442 #    my $time_del = $time_end - $time_start;
443 #    print OUTPUT "Delete: $time_del s.\n";
444     close(OUTPUT);
445 }
446
447 __END__
448
449 =head1 NAME
450
451 rsem-calculate-expression
452
453 =head1 SYNOPSIS
454
455  rsem-calculate-expression [options] upstream_read_file(s) reference_name sample_name 
456  rsem-calculate-expression [options] --paired-end upstream_read_file(s) downstream_read_file(s) reference_name sample_name 
457  rsem-calculate-expression [options] --sam/--bam [--paired-end] input reference_name sample_name
458
459 =head1 ARGUMENTS
460
461 =over
462
463 =item B<upstream_read_files(s)>
464
465 Comma-separated list of files containing single-end reads or upstream reads for paired-end data.  By default, these files are assumed to be in FASTQ format.  If the --no-qualities option is specified, then FASTA format is expected.
466
467 =item B<downstream_read_file(s)>
468
469 Comma-separated list of files containing downstream reads which are paired with the upstream reads.  By default, these files are assumed to be in FASTQ format.  If the --no-qualities option is specified, then FASTA format is expected.
470
471 =item B<input>
472
473 SAM/BAM formatted input file.  If "-" is specified for the filename, SAM/BAM input is instead assumed to come from standard input. RSEM requires all alignments of the same read group together. For paired-end reads, RSEM also requires the two mates of any alignment be adjacent. See Description section for how to make input file obey RSEM's requirements.
474
475 =item B<reference_name>                        
476
477 The name of the reference used.  The user must have run 'rsem-prepare-reference' with this reference_name before running this program.
478
479 =item B<sample_name>
480
481 The name of the sample analyzed. All output files are prefixed by this name (e.g., sample_name.genes.results)
482
483 =back
484
485 =head1 OPTIONS
486
487 =over
488
489 =item B<--paired-end>
490
491 Input reads are paired-end reads. (Default: off)
492
493 =item B<--no-qualities>
494
495 Input reads do not contain quality scores. (Default: off)
496
497 =item B<--strand-specific>
498
499 The RNA-Seq protocol used to generate the reads is strand specific, i.e., all (upstream) reads are derived from the forward strand.  This option is equivalent to --forward-prob=1.0.  With this option set, if RSEM runs the Bowtie/Bowtie 2 aligner, the '--norc' Bowtie/Bowtie 2 option will be used, which disables alignment to the reverse strand of transcripts.  (Default: off)
500
501 =item B<--sam>
502
503 Input file is in SAM format. (Default: off)
504
505 =item B<--bam>
506
507 Input file is in BAM format. (Default: off)
508
509 =item B<--sam-header-info> <file>
510
511 RSEM reads header information from input by default. If this option is on, header information is read from the specified file. For the format of the file, please see SAM official website. (Default: "")
512
513 =item B<-p/--num-threads> <int>
514
515 Number of threads to use. Both Bowtie/Bowtie2, expression estimation and 'samtools sort' will use this many threads. (Default: 1)
516
517 =item B<--no-bam-output>
518
519 Do not output any BAM file. (Default: off)
520
521 =item B<--output-genome-bam>
522
523 Generate a BAM file, 'sample_name.genome.bam', with alignments mapped to genomic coordinates and annotated with their posterior probabilities. In addition, RSEM will call samtools (included in RSEM package) to sort and index the bam file. 'sample_name.genome.sorted.bam' and 'sample_name.genome.sorted.bam.bai' will be generated. (Default: off)
524
525 =item B<--sampling-for-bam>
526
527 When RSEM generates a BAM file, instead of outputing all alignments a read has with their posterior probabilities, one alignment is sampled according to the posterior probabilities. The sampling procedure includes the alignment to the "noise" transcript, which does not appear in the BAM file. Only the sampled alignment has a weight of 1. All other alignments have weight 0. If the "noise" transcript is sampled, all alignments appeared in the BAM file should have weight 0. (Default: off)
528
529 =item B<--calc-pme>
530
531 Run RSEM's collapsed Gibbs sampler to calculate posterior mean estimates. (Default: off) 
532
533 =item B<--calc-ci>
534
535 Calculate 95% credibility intervals and posterior mean estimates. (Default: off)
536
537 =item B<--seed-length> <int>
538
539 Seed length used by the read aligner.  Providing the correct value is important for RSEM. If RSEM runs Bowtie, it uses this value for Bowtie's seed length parameter. Any read with its or at least one of its mates' (for paired-end reads) length less than this value will be ignored. If the references are not added poly(A) tails, the minimum allowed value is 5, otherwise, the minimum allowed value is 25. Note that this script will only check if the value >= 5 and give a warning message if the value < 25 but >= 5. (Default: 25)
540
541 =item B<--tag> <string>
542
543 The name of the optional field used in the SAM input for identifying a read with too many valid alignments. The field should have the format <tagName>:i:<value>, where a <value> bigger than 0 indicates a read with too many alignments. (Default: "")
544
545 =item B<--bowtie-path> <path>
546
547 The path to the Bowtie executables. (Default: the path to the Bowtie executables is assumed to be in the user's PATH environment variable)
548
549 =item B<--bowtie-n> <int>
550
551 (Bowtie parameter) max # of mismatches in the seed. (Range: 0-3, Default: 2)
552
553 =item B<--bowtie-e> <int>
554
555 (Bowtie parameter) max sum of mismatch quality scores across the alignment. (Default: 99999999)
556
557 =item B<--bowtie-m> <int>
558
559 (Bowtie parameter) suppress all alignments for a read if > <int> valid alignments exist. (Default: 200)
560
561 =item B<--bowtie-chunkmbs> <int>
562
563 (Bowtie parameter) memory allocated for best first alignment calculation (Default: 0 - use Bowtie's default)
564
565 =item B<--phred33-quals>
566
567 Input quality scores are encoded as Phred+33. (Default: on)
568
569 =item B<--phred64-quals>
570
571 Input quality scores are encoded as Phred+64 (default for GA Pipeline ver. >= 1.3). (Default: off)
572
573 =item B<--solexa-quals>
574
575 Input quality scores are solexa encoded (from GA Pipeline ver. < 1.3). (Default: off)
576
577 =item B<--bowtie2>
578
579 Use Bowtie 2 instead of Bowtie to align reads. Since currently RSEM does not handle indel, local and discordant alignments, the Bowtie2 parameters are set in a way to avoid those alignments. In particular, we use options '--sensitive --dpad 0 --gbar 99999999 --mp 1,1 --np 1 --score-min L,0,-0.1' by default. "-0.1", the last parameter of '--score-min' is the negative value of the maximum mismatch rate allowed. This rate can be set by option '--bowtie2-mismatch-rate'. If reads are paired-end, we additionally use options '--no-mixed' and '--no-discordant'. (Default: off)
580
581 =item B<--bowtie2-path> <path>
582
583 (Bowtie 2 parameter) The path to the Bowtie 2 executables. (Default: the path to the Bowtie 2 executables is assumed to be in the user's PATH environment variable)
584
585 =item B<--bowtie2-mismatch-rate> <double>
586
587 (Bowtie 2 parameter) The maximum mismatch rate allowed. (Default: 0.1)
588
589 =item B<--bowtie2-k> <int>
590
591 (Bowtie 2 parameter) Find up to <int> alignments per read. (Default: 200)
592
593 =item B<--bowtie2-sensitivity-level> <string>
594
595 (Bowtie 2 parameter) Set Bowtie 2's preset options in --end-to-end mode. This option controls how hard Bowtie 2 tries to find alignments. <string> must be one of "very_fast", "fast", "sensitive" and "very_sensitive". The four candidates correspond to Bowtie 2's "--very-fast", "--fast", "--sensitive" and "--very-sensitive" options. (Default: "sensitive" - use Bowtie 2's default)
596
597 =item B<--forward-prob> <double>
598
599 Probability of generating a read from the forward strand of a transcript. Set to 1 for a strand-specific protocol where all (upstream) reads are derived from the forward strand, 0 for a strand-specific protocol where all (upstream) read are derived from the reverse strand, or 0.5 for a non-strand-specific protocol. (Default: 0.5)
600
601 =item B<--fragment-length-min> <int>
602
603 Minimum read/insert length allowed. This is also the value for the Bowtie/Bowtie2 -I option. (Default: 1)
604
605 =item B<--fragment-length-max> <int>
606
607 Maximum read/insert length allowed. This is also the value for the Bowtie/Bowtie 2 -X option. (Default: 1000)
608
609 =item B<--fragment-length-mean> <double>
610
611 (single-end data only) The mean of the fragment length distribution, which is assumed to be a Gaussian. (Default: -1, which disables use of the fragment length distribution)
612
613 =item B<--fragment-length-sd> <double>
614
615 (single-end data only) The standard deviation of the fragment length distribution, which is assumed to be a Gaussian.  (Default: 0, which assumes that all fragments are of the same length, given by the rounded value of B<--fragment-length-mean>)
616
617 =item B<--estimate-rspd>
618
619 Set this option if you want to estimate the read start position distribution (RSPD) from data. Otherwise, RSEM will use a uniform RSPD. (Default: off)
620
621 =item B<--num-rspd-bins> <int>
622
623 Number of bins in the RSPD. Only relevant when '--estimate-rspd' is specified.  Use of the default setting is recommended. (Default: 20)
624
625 =item B<--ci-memory> <int>
626
627 Maximum size (in memory, MB) of the auxiliary buffer used for computing credibility intervals (CI). Set it larger for a faster CI calculation. However, leaving 2 GB memory free for other usage is recommended. (Default: 1024)
628
629 =item B<--samtools-sort-mem> <string>
630
631 Set the maximum memory per thread that can be used by 'samtools sort'. <string> represents the memory and accepts suffices 'K/M/G'. RSEM will pass <string> to the '-m' option of 'samtools sort'.  Please note that the default used here is different from the default used by samtools. (Default: 1G)
632
633 =item B<--keep-intermediate-files>
634
635 Keep temporary files generated by RSEM.  RSEM creates a temporary directory, 'sample_name.temp', into which it puts all intermediate output files. If this directory already exists, RSEM overwrites all files generated by previous RSEM runs inside of it. By default, after RSEM finishes, the temporary directory is deleted.  Set this option to prevent the deletion of this directory and the intermediate files inside of it. (Default: off)
636
637 =item B<--temporary-folder> <string>
638
639 Set where to put the temporary files generated by RSEM. If the folder specified does not exist, RSEM will try to create it. (Default: sample_name.temp)
640
641 =item B<--time>
642
643 Output time consumed by each step of RSEM to 'sample_name.time'. (Default: off)
644
645 =item B<-q/--quiet>
646
647 Suppress the output of logging information. (Default: off)
648
649 =item B<-h/--help>
650
651 Show help information.
652
653 =item B<--version>
654
655 Show version information.
656
657 =back
658
659 =head1 DESCRIPTION
660
661 In its default mode, this program aligns input reads against a reference transcriptome with Bowtie and calculates expression values using the alignments.  RSEM assumes the data are single-end reads with quality scores, unless the '--paired-end' or '--no-qualities' options are specified.  Users may use an alternative aligner by specifying one of the --sam and --bam options, and providing an alignment file in the specified format. However, users should make sure that they align against the indices generated by 'rsem-prepare-reference' and the alignment file satisfies the requirements mentioned in ARGUMENTS section. 
662
663 One simple way to make the alignment file satisfying RSEM's requirements (assuming the aligner used put mates in a paired-end read adjacent) is to use 'convert-sam-for-rsem' script. This script only accept SAM format files as input. If a BAM format file is obtained, please use samtools to convert it to a SAM file first. For example, if '/ref/mouse_125' is the 'reference_name' and the SAM file is named 'input.sam', you can run the following command: 
664
665   convert-sam-for-rsem /ref/mouse_125 input.sam -o input_for_rsem.sam  
666
667 For details, please refer to 'convert-sam-for-rsem's documentation page.
668
669 The SAM/BAM format RSEM uses is v1.4. However, it is compatible with old SAM/BAM format. However, RSEM cannot recognize 0x100 in the FLAG field. In addition, RSEM requires SEQ and QUAL are not '*'. 
670
671 The user must run 'rsem-prepare-reference' with the appropriate reference before using this program.
672
673 For single-end data, it is strongly recommended that the user provide the fragment length distribution parameters (--fragment-length-mean and --fragment-length-sd).  For paired-end data, RSEM will automatically learn a fragment length distribution from the data.
674
675 Please note that some of the default values for the Bowtie parameters are not the same as those defined for Bowtie itself.
676
677 The temporary directory and all intermediate files will be removed when RSEM finishes unless '--keep-intermediate-files' is specified.
678
679 With the '--calc-ci' option, 95% credibility intervals and posterior mean estimates will be calculated in addition to maximum likelihood estimates.
680
681 =head1 OUTPUT
682
683 =over
684
685 =item B<sample_name.isoforms.results> 
686
687 File containing isoform level expression estimates. The first line
688 contains column names separated by the tab character. The format of
689 each line in the rest of this file is:
690
691 transcript_id gene_id length effective_length expected_count TPM FPKM IsoPct [pme_expected_count pme_TPM pme_FPKM IsoPct_from_pme_TPM TPM_ci_lower_bound TPM_ci_upper_bound FPKM_ci_lower_bound FPKM_ci_upper_bound]
692
693 Fields are separated by the tab character. Fields within "[]" are
694 optional. They will not be presented if neither '--calc-pme' nor
695 '--calc-ci' is set.
696
697 'transcript_id' is the transcript name of this transcript. 'gene_id'
698 is the gene name of the gene which this transcript belongs to (denote
699 this gene as its parent gene). If no gene information is provided,
700 'gene_id' and 'transcript_id' are the same.
701
702 'length' is this transcript's sequence length (poly(A) tail is not
703 counted). 'effective_length' counts only the positions that can
704 generate a valid fragment. If no poly(A) tail is added,
705 'effective_length' is equal to transcript length - mean fragment
706 length + 1. If one transcript's effective length is less than 1, this
707 transcript's both effective length and abundance estimates are set to
708 0.
709
710 'expected_count' is the sum of the posterior probability of each read
711 comes from this transcript over all reads. Because 1) each read
712 aligning to this transcript has a probability of being generated from
713 background noise; 2) RSEM may filter some alignable low quality reads,
714 the sum of expected counts for all transcript are generally less than
715 the total number of reads aligned.
716
717 'TPM' stands for Transcripts Per Million. It is a relative measure of
718 transcript abundance. The sum of all transcripts' TPM is 1
719 million. 'FPKM' stands for Fragments Per Kilobase of transcript per
720 Million mapped reads. It is another relative measure of transcript
721 abundance. If we define l_bar be the mean transcript length in a
722 sample, which can be calculated as
723
724 l_bar = \sum_i TPM_i / 10^6 * effective_length_i (i goes through every transcript), 
725
726 the following equation is hold:
727
728 FPKM_i = 10^3 / l_bar * TPM_i.
729
730 We can see that the sum of FPKM is not a constant across samples.
731
732 'IsoPct' stands for isoform percentage. It is the percentage of this
733 transcript's abandunce over its parent gene's abandunce. If its parent
734 gene has only one isoform or the gene information is not provided,
735 this field will be set to 100.
736
737 'pme_expected_count', 'pme_TPM', 'pme_FPKM' are posterior mean
738 estimates calculated by RSEM's Gibbs sampler. 'IsoPct_from_pme_TPM' is
739 the isoform percentage calculated from 'pme_TPM' values.
740
741 'TPM_ci_lower_bound', 'TPM_ci_upper_bound', 'FPKM_ci_lower_bound' and
742 'FPKM_ci_upper_bound' are lower(l) and upper(u) bounds of 95%
743 credibility intervals for TPM and FPKM values. The bounds are
744 inclusive (i.e. [l, u]).
745
746 =item B<sample_name.genes.results>
747
748 File containing gene level expression estimates. The first line
749 contains column names separated by the tab character. The format of
750 each line in the rest of this file is:
751
752 gene_id transcript_id(s) length effective_length expected_count TPM FPKM [pme_expected_count pme_TPM pme_FPKM TPM_ci_lower_bound TPM_ci_upper_bound FPKM_ci_lower_bound FPKM_ci_upper_bound]
753
754 Fields are separated by the tab character. Fields within "[]" are
755 optional. They will not be presented if neither '--calc-pme' nor
756 '--calc-ci' is set.
757
758 'transcript_id(s)' is a comma-separated list of transcript_ids
759 belonging to this gene. If no gene information is provided, 'gene_id'
760 and 'transcript_id(s)' are identical (the 'transcript_id').
761
762 A gene's 'length' and 'effective_length' are
763 defined as the weighted average of its transcripts' lengths and
764 effective lengths (weighted by 'IsoPct'). A gene's abundance estimates
765 are just the sum of its transcripts' abundance estimates.
766
767 =item B<sample_name.alleles.results>
768
769 Only generated when the RSEM references are built with allele-specific
770 transcripts.
771
772 This file contains allele level expression estimates for
773 allele-specific expression calculation. The first line
774 contains column names separated by the tab character. The format of
775 each line in the rest of this file is:
776
777 allele_id transcript_id gene_id length effective_length expected_count TPM FPKM AlleleIsoPct AlleleGenePct [pme_expected_count pme_TPM pme_FPKM AlleleIsoPct_from_pme_TPM AlleleGenePct_from_pme_TPM TPM_ci_lower_bound TPM_ci_upper_bound FPKM_ci_lower_bound FPKM_ci_upper_bound]
778
779 Fields are separated by the tab character. Fields within "[]" are
780 optional. They will not be presented if neither '--calc-pme' nor
781 '--calc-ci' is set.
782
783 'allele_id' is the allele-specific name of this allele-specific transcript.
784
785 'AlleleIsoPct' stands for allele-specific percentage on isoform
786 level. It is the percentage of this allele-specific transcript's
787 abundance over its parent transcript's abundance. If its parent
788 transcript has only one allele variant form, this field will be set to
789 100.
790
791 'AlleleGenePct' stands for allele-specific percentage on gene
792 level. It is the percentage of this allele-specific transcript's
793 abundance over its parent gene's abundance.
794
795 'AlleleIsoPct_from_pme_TPM' and 'AlleleGenePct_from_pme_TPM' have
796 similar meanings. They are calculated based on posterior mean
797 estimates.
798
799 Please note that if this file is present, the fields 'length' and
800 'effective_length' in 'sample_name.isoforms.results' should be
801 interpreted similarly as the corresponding definitions in
802 'sample_name.genes.results'.
803
804 =item B<sample_name.transcript.bam, sample_name.transcript.sorted.bam and sample_name.transcript.sorted.bam.bai>
805
806 Only generated when --no-bam-output is not specified.
807
808 'sample_name.transcript.bam' is a BAM-formatted file of read
809 alignments in transcript coordinates. The MAPQ field of each alignment
810 is set to min(100, floor(-10 * log10(1.0 - w) + 0.5)), where w is the
811 posterior probability of that alignment being the true mapping of a
812 read.  In addition, RSEM pads a new tag ZW:f:value, where value is a
813 single precision floating number representing the posterior
814 probability. Because this file contains all alignment lines produced
815 by bowtie or user-specified aligners, it can also be used as a
816 replacement of the aligner generated BAM/SAM file. For paired-end
817 reads, if one mate has alignments but the other does not, this file
818 marks the alignable mate as "unmappable" (flag bit 0x4) and appends an
819 optional field "Z0:A:!".
820
821 'sample_name.transcript.sorted.bam' and
822 'sample_name.transcript.sorted.bam.bai' are the sorted BAM file and
823 indices generated by samtools (included in RSEM package).
824
825 =item B<sample_name.genome.bam, sample_name.genome.sorted.bam and sample_name.genome.sorted.bam.bai>
826
827 Only generated when --no-bam-output is not specified and --output-genome-bam is specified.
828
829 'sample_name.genome.bam' is a BAM-formatted file of read alignments in
830 genomic coordinates. Alignments of reads that have identical genomic
831 coordinates (i.e., alignments to different isoforms that share the
832 same genomic region) are collapsed into one alignment.  The MAPQ field
833 of each alignment is set to min(100, floor(-10 * log10(1.0 - w) +
834 0.5)), where w is the posterior probability of that alignment being
835 the true mapping of a read.  In addition, RSEM pads a new tag
836 ZW:f:value, where value is a single precision floating number
837 representing the posterior probability. If an alignment is spliced, a
838 XS:A:value tag is also added, where value is either '+' or '-'
839 indicating the strand of the transcript it aligns to.
840
841 'sample_name.genome.sorted.bam' and 'sample_name.genome.sorted.bam.bai' are the
842 sorted BAM file and indices generated by samtools (included in RSEM package).
843
844 =item B<sample_name.time>
845
846 Only generated when --time is specified.
847
848 It contains time (in seconds) consumed by aligning reads, estimating expression levels and calculating credibility intervals.
849
850 =item B<sample_name.stat>
851
852 This is a folder instead of a file. All model related statistics are stored in this folder. Use 'rsem-plot-model' can generate plots using this folder.
853
854 =back
855
856 =head1 EXAMPLES
857
858 Assume the path to the bowtie executables is in the user's PATH environment variable. Reference files are under '/ref' with name 'mouse_125'. 
859
860 1) '/data/mmliver.fq', single-end reads with quality scores. Quality scores are encoded as for 'GA pipeline version >= 1.3'. We want to use 8 threads and generate a genome BAM file:
861
862  rsem-calculate-expression --phred64-quals \
863                            -p 8 \
864                            --output-genome-bam \
865                            /data/mmliver.fq \
866                            /ref/mouse_125 \
867                            mmliver_single_quals
868
869 2) '/data/mmliver_1.fq' and '/data/mmliver_2.fq', paired-end reads with quality scores. Quality scores are in SANGER format. We want to use 8 threads and do not generate a genome BAM file:
870
871  rsem-calculate-expression -p 8 \
872                            --paired-end \
873                            /data/mmliver_1.fq \
874                            /data/mmliver_2.fq \
875                            /ref/mouse_125 \
876                            mmliver_paired_end_quals
877
878 3) '/data/mmliver.fa', single-end reads without quality scores. We want to use 8 threads:
879
880  rsem-calculate-expression -p 8 \
881                            --no-qualities \
882                            /data/mmliver.fa \
883                            /ref/mouse_125 \
884                            mmliver_single_without_quals
885
886 4) Data are the same as 1). This time we assume the bowtie executables are under '/sw/bowtie'. We want to take a fragment length distribution into consideration. We set the fragment length mean to 150 and the standard deviation to 35. In addition to a BAM file, we also want to generate credibility intervals. We allow RSEM to use 1GB of memory for CI calculation:
887
888  rsem-calculate-expression --bowtie-path /sw/bowtie \
889                            --phred64-quals \
890                            --fragment-length-mean 150.0 \
891                            --fragment-length-sd 35.0 \
892                            -p 8 \
893                            --output-genome-bam \
894                            --calc-ci \
895                            --ci-memory 1024 \
896                            /data/mmliver.fq \
897                            /ref/mouse_125 \
898                            mmliver_single_quals
899
900 5) '/data/mmliver_paired_end_quals.bam', paired-end reads with quality scores. We want to use 8 threads:
901
902  rsem-calculate-expression --paired-end \
903                            --bam \
904                            -p 8 \
905                            /data/mmliver_paired_end_quals.bam \
906                            /ref/mouse_125 \
907                            mmliver_paired_end_quals
908
909 =cut