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