]> git.donarmstrong.com Git - rsem.git/blob - rsem-calculate-expression
changed output format to contain FPKM etc. ; fixed a bug for paired-end reads
[rsem.git] / rsem-calculate-expression
1 #!/usr/bin/perl
2
3 use Getopt::Long;
4 use Pod::Usage;
5 use File::Basename;
6 use Switch;
7 use strict;
8
9 #const
10 my $BURNIN = 200;
11 my $NCV = 1000;
12 my $SAMPLEGAP = 1;
13 my $CONFIDENCE = 0.95;
14 my $NSPC = 50;
15
16 my $NMB = 1024; # default
17
18 my $status = 0;
19
20 my $read_type = 1; # default, single end with qual
21
22 my @transcript_title = ("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");
23
24 my @gene_title = ("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");
25
26 my $bowtie_path = "";
27 my $C = 2;
28 my $E = 99999999;
29 my $L = 25;
30 my $maxHits = 200;
31 my $chunkMbs = 0;       # 0 = use bowtie default
32 my $phred33 = 0;
33 my $phred64 = 0;
34 my $solexa = 0;
35
36 my $is_sam = 0;
37 my $is_bam = 0;
38 my $fn_list = "";
39 my $tagName = "XM";
40
41 my $probF = 0.5;
42
43 my $minL = 1;
44 my $maxL = 1000;
45 my $mean = -1;
46 my $sd = 0;
47
48 my $estRSPD = 0;
49 my $B = 20;
50
51 my $nThreads = 1;
52 my $genBamF = 1;  # default is generating transcript bam file
53 my $genGenomeBamF = 0;
54 my $sampling = 0;
55 my $calcCI = 0;
56 my $var_opt = 0; # temporarily, only for internal use
57 my $quiet = 0;
58 my $help = 0;
59
60 my $paired_end = 0;
61 my $no_qual = 0;
62 my $keep_intermediate_files = 0;
63
64 my $strand_specific = 0;
65
66 my $mTime = 0;
67 my ($time_start, $time_end, $time_alignment, $time_rsem, $time_ci) = (0, 0, 0, 0, 0);
68
69 my $mate1_list = "";
70 my $mate2_list = "";
71 my $inpF = "";
72
73 my ($refName, $sampleName, $sampleToken, $temp_dir, $stat_dir, $imdName, $statName) = ();
74 my $gap = 32;
75
76 GetOptions("keep-intermediate-files" => \$keep_intermediate_files,
77            "temporary-folder=s" => \$temp_dir,
78            "no-qualities" => \$no_qual,
79            "paired-end" => \$paired_end,
80            "strand-specific" => \$strand_specific,
81            "sam" => \$is_sam,
82            "bam" => \$is_bam,
83            "sam-header-info=s" => \$fn_list,
84            "tag=s" => \$tagName,
85            "seed-length=i" => \$L,
86            "bowtie-path=s" => \$bowtie_path,
87            "bowtie-n=i" => \$C,
88            "bowtie-e=i" => \$E,
89            "bowtie-m=i" => \$maxHits,
90            "bowtie-chunkmbs=i" => \$chunkMbs,
91            "phred33-quals" => \$phred33,
92            "phred64-quals" => \$phred64, #solexa1.3-quals" => \$phred64,
93            "solexa-quals" => \$solexa,
94            "forward-prob=f" => \$probF,
95            "fragment-length-min=i" => \$minL,
96            "fragment-length-max=i" => \$maxL,
97            "fragment-length-mean=f" => \$mean,
98            "fragment-length-sd=f" => \$sd,
99            "estimate-rspd" => \$estRSPD,
100            "num-rspd-bins=i" => \$B,
101            "p|num-threads=i" => \$nThreads,
102            "no-bam-output" => sub { $genBamF = 0; },
103            "output-genome-bam" => \$genGenomeBamF,
104            "sampling-for-bam" => \$sampling,
105            "var" => \$var_opt,
106            "calc-ci" => \$calcCI,
107            "ci-memory=i" => \$NMB,
108            "time" => \$mTime,
109            "q|quiet" => \$quiet,
110            "h|help" => \$help) or pod2usage(-exitval => 2, -verbose => 2);
111
112 pod2usage(-verbose => 2) if ($help == 1);
113
114
115 #check parameters and options
116
117 if ($is_sam || $is_bam) {
118     pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (scalar(@ARGV) != 3);
119     pod2usage(-msg => "--sam and --bam cannot be active at the same time!", -exitval => 2, -verbose => 2) if ($is_sam == 1&& $is_bam == 1);
120     pod2usage(-msg => "--bowtie-path, --bowtie-n, --bowtie-e, --bowtie-m, --phred33-quals, --phred64-quals or --solexa-quals 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);
121 }
122 else {
123     pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (!$paired_end && scalar(@ARGV) != 3 || $paired_end && scalar(@ARGV) != 4);    
124     pod2usage(-msg => "Only one of --phred33-quals --phred64-quals/--solexa1.3-quals --solexa-suqls can be active!", -exitval => 2, -verbose => 2) if ($phred33 + $phred64 + $solexa > 1);    
125     podwusage(-msg => "--sam , --bam or --sam-header-info cannot be set if use bowtie aligner to produce alignments!", -exitval => 2, -verbose => 2) if ($is_sam || $is_bam || $fn_list ne "");
126 }
127
128 pod2usage(-msg => "Forward probability should be in [0, 1]!", -exitval => 2, -verbose => 2) if ($probF < 0 || $probF > 1);
129 pod2usage(-msg => "Min fragment length should be at least 1!", -exitval => 2, -verbose => 2) if ($minL < 1);
130 pod2usage(-msg => "Min fragment length should be smaller or equal to max fragment length!", -exitval => 2, -verbose => 2) if ($minL > $maxL);
131 pod2usage(-msg => "The memory allocated for calculating credibility intervals should be at least 1 MB!\n", -exitval => 2, -verbose => 2) if ($NMB < 1);
132 pod2usage(-msg => "Number of threads should be at least 1!\n", -exitval => 2, -verbose => 2) if ($nThreads < 1);
133 pod2usage(-msg => "Seed length should be at least 5!\n", -exitval => 2, -verbose => 2) if ($L < 5);
134 pod2usage(-msg => "--sampling-for-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($sampling && !$genBamF);
135 pod2usage(-msg => "--output-genome-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($genGenomeBamF && !$genBamF);
136
137 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"; }
138
139 if ($strand_specific) { $probF = 1.0; }
140
141 if ($paired_end) {
142     if ($no_qual) { $read_type = 2; }
143     else { $read_type = 3; }
144 }
145 else {
146     if ($no_qual) { $read_type = 0; }
147     else { $read_type = 1; }
148 }
149
150 if (scalar(@ARGV) == 3) {
151     if ($is_sam || $is_bam) { $inpF = $ARGV[0]; } 
152     else {$mate1_list = $ARGV[0]; }
153     $refName = $ARGV[1];
154     $sampleName = $ARGV[2];
155 }
156 else {
157     $mate1_list = $ARGV[0];
158     $mate2_list = $ARGV[1];
159     $refName = $ARGV[2];
160     $sampleName = $ARGV[3];
161 }
162
163 if ($genGenomeBamF) {
164     open(INPUT, "$refName.ti");
165     my $line = <INPUT>; chomp($line);
166     close(INPUT);
167     my ($M, $type) = split(/ /, $line);
168     pod2usage(-msg => "No genome information provided, so genome bam file cannot be generated!\n", -exitval => 2, -verbose => 2) if ($type != 0);
169 }
170
171 my $pos = rindex($sampleName, '/');
172 if ($pos < 0) { $sampleToken = $sampleName; }
173 else { $sampleToken = substr($sampleName, $pos + 1); }
174
175 if ($temp_dir eq "") { $temp_dir = "$sampleName.temp"; }
176 $stat_dir = "$sampleName.stat";
177
178 if (!(-d $temp_dir) && !mkdir($temp_dir)) { print "Fail to create folder $temp_dir.\n"; exit(-1); }
179 if (!(-d $stat_dir) && !mkdir($stat_dir)) { print "Fail to create folder $stat_dir.\n"; exit(-1); }
180
181 $imdName = "$temp_dir/$sampleToken";
182 $statName = "$stat_dir/$sampleToken";
183
184 if (!$is_sam && !$is_bam && !$no_qual && ($phred33 + $phred64 + $solexa == 0)) { $phred33 = 1; }
185
186 my ($mate_minL, $mate_maxL) = (1, $maxL);
187
188 if ($bowtie_path ne "") { $bowtie_path .= "/"; }
189
190 my ($fn, $dir, $suf) = fileparse($0);
191 my $command = "";
192
193 if (!$is_sam && !$is_bam) {
194     $command = $bowtie_path."bowtie";
195     if ($no_qual) { $command .= " -f"; }
196     else { $command .= " -q"; }
197     
198     if ($phred33) { $command .= " --phred33-quals"; }
199     elsif ($phred64) { $command .= " --phred64-quals"; }
200     elsif ($solexa) { $command .= " --solexa-quals"; }
201     
202     $command .= " -n $C -e $E -l $L";
203     if ($read_type == 2 || $read_type == 3) { $command .= " -I $minL -X $maxL"; }
204     if ($chunkMbs > 0) { $command .= " --chunkmbs $chunkMbs"; }
205     
206     if ($strand_specific || $probF == 1.0) { $command .= " --norc"; }
207     elsif ($probF == 0.0) { $command .= " --nofw"; }
208
209     $command .= " -p $nThreads -a -m $maxHits -S";
210     if ($quiet) { $command .= " --quiet"; }    
211
212     $command .= " $refName";
213     if ($read_type == 0 || $read_type == 1) {
214         $command .= " $mate1_list"; 
215     }
216     else {
217         $command .= " -1 $mate1_list -2 $mate2_list";
218     }
219
220     # pipe to samtools to generate a BAM file
221     $command .= " | $dir\sam/samtools view -S -b -o $imdName.bam -";
222
223     if ($mTime) { $time_start = time(); }
224
225     &runCommand($command);
226
227     if ($mTime) { $time_end = time(); $time_alignment = $time_end - $time_start; }
228
229     $inpF = "$imdName.bam";
230     $is_bam = 1; # alignments are outputed as a BAM file
231 }
232
233 if ($mTime) { $time_start = time(); }
234
235 $command = $dir."rsem-parse-alignments $refName $imdName $statName";
236
237 my $samInpType;
238 if ($is_sam) { $samInpType = "s"; } 
239 elsif ($is_bam) { $samInpType = "b"; }
240
241 $command .= " $samInpType $inpF -t $read_type";
242 if ($fn_list ne "") { $command .= " -l $fn_list"; }
243 if ($tagName ne "") { $command .= " -tag $tagName"; }
244 if ($quiet) { $command .= " -q"; }
245
246 &runCommand($command);
247
248 $command = $dir."rsem-build-read-index $gap"; 
249 switch($read_type) {
250     case 0  { $command .= " 0 $quiet $imdName\_alignable.fa"; }
251     case 1  { $command .= " 1 $quiet $imdName\_alignable.fq"; }
252     case 2  { $command .= " 0 $quiet $imdName\_alignable_1.fa $imdName\_alignable_2.fa"; }
253     case 3  { $command .= " 1 $quiet $imdName\_alignable_1.fq $imdName\_alignable_2.fq"; }
254 }
255 &runCommand($command);
256
257 my $doesOpen = open(OUTPUT, ">$imdName.mparams");
258 if ($doesOpen == 0) { print "Cannot generate $imdName.mparams!\n"; exit(-1); }
259 print OUTPUT "$minL $maxL\n";
260 print OUTPUT "$probF\n";
261 print OUTPUT "$estRSPD\n";
262 print OUTPUT "$B\n";
263 print OUTPUT "$mate_minL $mate_maxL\n";
264 print OUTPUT "$mean $sd\n";
265 print OUTPUT "$L\n";
266 close(OUTPUT);  
267
268 $command = $dir."rsem-run-em $refName $read_type $sampleName $imdName $statName -p $nThreads";
269 if ($genBamF) { 
270     $command .= " -b $samInpType $inpF";
271     if ($fn_list ne "") { $command .= " 1 $fn_list"; }
272     else { $command .= " 0"; }
273     if ($sampling) { $command .= " --sampling"; }
274 }
275 if ($calcCI || $var_opt) { $command .= " --gibbs-out"; }
276 if ($quiet) { $command .= " -q"; }
277
278 &runCommand($command);
279
280 &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
281 &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
282
283 if ($genBamF) {
284     $command = $dir."sam/samtools sort $sampleName.transcript.bam $sampleName.transcript.sorted";
285     &runCommand($command);
286     $command = $dir."sam/samtools index $sampleName.transcript.sorted.bam";
287     &runCommand($command);
288
289     if ($genGenomeBamF) {
290         $command = $dir."rsem-tbam2gbam $refName $sampleName.transcript.bam $sampleName.genome.bam";
291         &runCommand($command);
292         $command = $dir."sam/samtools sort $sampleName.genome.bam $sampleName.genome.sorted";
293         &runCommand($command);
294         $command = $dir."sam/samtools index $sampleName.genome.sorted.bam";
295         &runCommand($command);
296     }
297 }
298
299 if ($mTime) { $time_end = time(); $time_rsem = $time_end - $time_start; }
300
301 if ($mTime) { $time_start = time(); }
302
303 if ($calcCI || $var_opt) {
304     $command = $dir."rsem-run-gibbs $refName $imdName $statName $BURNIN $NCV $SAMPLEGAP";
305     $command .= " -p $nThreads";
306     if ($var_opt) { $command .= " --var"; }
307     if ($quiet) { $command .= " -q"; }
308     &runCommand($command);
309 }
310
311 if ($calcCI) {
312     system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak1");
313     system("mv $sampleName.genes.results $imdName.genes.results.bak1");
314     &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
315     &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
316
317     $command = $dir."rsem-calculate-credibility-intervals $refName $imdName $statName $CONFIDENCE $NCV $NSPC $NMB";
318     $command .= " -p $nThreads";
319     if ($quiet) { $command .= " -q"; }
320     &runCommand($command);
321
322     system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak2");
323     system("mv $sampleName.genes.results $imdName.genes.results.bak2");
324     &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level
325     &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level
326 }
327
328 if ($mTime) { $time_end = time(); $time_ci = $time_end - $time_start; }
329
330 if ($mTime) { $time_start = time(); }
331
332 if (!$keep_intermediate_files) {
333     &runCommand("rm -rf $temp_dir", "Fail to delete the temporary folder!");
334 }
335
336 if ($mTime) { $time_end = time(); }
337
338 if ($mTime) { 
339     open(OUTPUT, ">$sampleName.time");
340     print OUTPUT "Aligning reads: $time_alignment s.\n";
341     print OUTPUT "Estimating expression levels: $time_rsem s.\n";
342     print OUTPUT "Calculating credibility intervals: $time_ci s.\n";
343     my $time_del = $time_end - $time_start;
344 #    print OUTPUT "Delete: $time_del s.\n";
345     close(OUTPUT);
346 }
347
348 # command, {err_msg}
349 sub runCommand {
350     print $_[0]."\n";
351     my $status = system($_[0]);
352     if ($status != 0) {
353         my $errmsg = "";
354         if (scalar(@_) > 1) { $errmsg .= $_[1]."\n"; }
355         $errmsg .= "\"$_[0]\" failed! Plase check if you provide correct parameters/options for the pipeline!\n";
356         print $errmsg;
357         exit(-1);
358     }
359     print "\n";
360 }
361  
362 # inpF, outF
363 sub collectResults {
364     my $local_status;
365     my ($inpF, $outF);
366     my @results = ();
367     my $line;
368
369     $inpF = $_[1];
370     $outF = $_[2];
371
372     $local_status = open(INPUT, $inpF);
373     if ($local_status == 0) { print "Fail to open file $inpF!\n"; exit(-1); }
374     
375     @results = ();
376     
377     while ($line = <INPUT>) {
378         chomp($line);
379         my @local_arr = split(/\t/, $line);
380         push(@results, \@local_arr); 
381     }
382
383     close(INPUT);
384
385     $local_status = open(OUTPUT, ">$outF");
386     if ($local_status == 0) { print "Fail to create file $outF!\n"; exit(-1); }
387
388     my $n = scalar(@results);
389     my $m = scalar(@{$results[0]});
390
391     $" = "\t";
392
393     my @out_arr = ();
394     for (my $i = 0; $i < $n; $i++) {
395         if ($_[0] eq "isoform") { push(@out_arr, $transcript_title[$i]); }
396         elsif ($_[0] eq "gene") { push(@out_arr, $gene_title[$i]); }
397         else { print "A bug on 'collectResults' is detected!\n"; exit(-1); }
398     }
399     print OUTPUT "@out_arr\n";
400
401     for (my $i = 0; $i < $m; $i++) {
402         @out_arr = ();
403         for (my $j = 0; $j < $n; $j++) { push(@out_arr, $results[$j][$i]); }
404         print OUTPUT "@out_arr\n"; 
405     }
406
407     close(OUTPUT);
408 }
409
410 __END__
411
412 =head1 NAME
413
414 rsem-calculate-expression
415
416 =head1 SYNOPSIS
417
418 =over
419
420  rsem-calculate-expression [options] upstream_read_file(s) reference_name sample_name
421  rsem-calculate-expression [options] --paired-end upstream_read_file(s) downstream_read_file(s) reference_name sample_name
422  rsem-calculate-expression [options] --sam/--bam [--paired-end] input reference_name sample_name
423
424 =back
425
426 =head1 ARGUMENTS
427
428 =over
429
430 =item B<upstream_read_files(s)>
431
432 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.
433
434 =item B<downstream_read_file(s)>
435
436 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.
437
438 =item B<input>
439
440 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.
441
442 =item B<reference_name>                        
443
444 The name of the reference used.  The user must have run 'rsem-prepare-reference' with this reference_name before running this program.
445
446 =item B<sample_name>
447
448 The name of the sample analyzed. All output files are prefixed by this name (e.g., sample_name.genes.results)
449
450 =back
451
452 =head1 OPTIONS
453
454 =over
455
456 =item B<--paired-end>
457
458 Input reads are paired-end reads. (Default: off)
459
460 =item B<--no-qualities>
461
462 Input reads do not contain quality scores. (Default: off)
463
464 =item B<--strand-specific>
465
466 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 aligner, the '--norc' Bowtie option will be used, which disables alignment to the reverse strand of transcripts.  (Default: off)
467
468 =item B<--sam>
469
470 Input file is in SAM format. (Default: off)
471
472 =item B<--bam>
473
474 Input file is in BAM format. (Default: off)
475
476 =item B<--sam-header-info> <file>
477
478 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: "")
479
480 =item B<-p/--num-threads> <int>
481
482 Number of threads to use. Both Bowtie and expression estimation will use this many threads. (Default: 1)
483
484 =item B<--no-bam-output>
485
486 Do not output any BAM file. (Default: off)
487
488 =item B<--output-genome-bam>
489
490 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)
491
492 =item B<--sampling-for-bam>
493
494 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)
495
496 =item B<--calc-ci>
497
498 Calculate 95% credibility intervals and posterior mean estimates.  (Default: off)
499
500 =item B<--seed-length> <int>
501
502 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)
503
504 =item B<--tag> <string>
505
506 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: "")
507
508 =item B<--bowtie-path> <path>
509
510 The path to the bowtie executables. (Default: the path to the bowtie executables is assumed to be in the user's PATH environment variable)
511
512 =item B<--bowtie-n> <int>
513
514 (Bowtie parameter) max # of mismatches in the seed. (Range: 0-3, Default: 2)
515
516 =item B<--bowtie-e> <int>
517
518 (Bowtie parameter) max sum of mismatch quality scores across the alignment. (Default: 99999999)
519
520 =item B<--bowtie-m> <int>
521
522 (Bowtie parameter) suppress all alignments for a read if > <int> valid alignments exist. (Default: 200)
523
524 =item B<--bowtie-chunkmbs> <int>
525
526 (Bowtie parameter) memory allocated for best first alignment calculation (Default: 0 - use bowtie's default)
527
528 =item B<--phred33-quals>
529
530 Input quality scores are encoded as Phred+33. (Default: on)
531
532 =item B<--phred64-quals>
533
534 Input quality scores are encoded as Phred+64 (default for GA Pipeline ver. >= 1.3). (Default: off)
535
536 =item B<--solexa-quals>
537
538 Input quality scores are solexa encoded (from GA Pipeline ver. < 1.3). (Default: off)
539
540 =item B<--forward-prob> <double>
541
542 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)
543
544 =item B<--fragment-length-min> <int>
545
546 Minimum read/insert length allowed. This is also the value for the bowtie -I option. (Default: 1)
547
548 =item B<--fragment-length-max> <int>
549
550 Maximum read/insert length allowed. This is also the value for the bowtie -X option. (Default: 1000)
551
552 =item B<--fragment-length-mean> <double>
553
554 (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)
555
556 =item B<--fragment-length-sd> <double>
557
558 (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>)
559
560 =item B<--estimate-rspd>
561
562 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)
563
564 =item B<--num-rspd-bins> <int>
565
566 Number of bins in the RSPD. Only relevant when '--estimate-rspd' is specified.  Use of the default setting is recommended. (Default: 20)
567
568 =item B<--ci-memory> <int>
569
570 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)
571
572 =item B<--keep-intermediate-files>
573
574 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)
575
576 =item B<--temporary-folder> <string>
577
578 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)
579
580 =item B<--time>
581
582 Output time consumed by each step of RSEM to 'sample_name.time'. (Default: off)
583
584 =item B<-q/--quiet>
585
586 Suppress the output of logging information. (Default: off)
587
588 =item B<-h/--help>
589
590 Show help information.
591
592 =back
593
594 =head1 DESCRIPTION
595
596 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. 
597
598 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: 
599
600   convert-sam-for-rsem /ref/mouse_125 input.sam -o input_for_rsem.sam  
601
602 For details, please refer to 'convert-sam-for-rsem's documentation page.
603
604 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 '*'. 
605
606 The user must run 'rsem-prepare-reference' with the appropriate reference before using this program.
607
608 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.
609
610 Please note that some of the default values for the Bowtie parameters are not the same as those defined for Bowtie itself.
611
612 The temporary directory and all intermediate files will be removed when RSEM finishes unless '--keep-intermediate-files' is specified.
613
614 With the '--calc-ci' option, 95% credibility intervals and posterior mean estimates will be calculated in addition to maximum likelihood estimates.
615
616 =head1 OUTPUT
617
618 =over
619
620 =item B<sample_name.isoforms.results> 
621
622 File containing isoform level expression estimates. The first line
623 contains column names separated by the tab character. The format of
624 each line in the rest of this file is:
625
626 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]
627
628 Fields are separated by the tab character. Fields within "[]" are only
629 presented if '--calc-ci' is set.
630
631 'transcript_id' is the transcript name of this transcript. 'gene_id'
632 is the gene name of the gene which this transcript belongs to (denote
633 this gene as its parent gene). If no gene information is provided,
634 'gene_id' and 'transcript_id' are the same.
635
636 'length' is this transcript's sequence length (poly(A) tail is not
637 counted). 'effective_length' counts only the positions that can
638 generate a valid fragment. If no poly(A) tail is added,
639 'effective_length' is equal to transcript length - mean fragment
640 length + 1. If one transcript's effective length is less than 1, this
641 transcript's both effective length and abundance estimates are set to
642 0.
643
644 'expected_count' is the sum of the posterior probability of each read
645 comes from this transcript over all reads. Because 1) each read
646 aligning to this transcript has a probability of being generated from
647 background noise; 2) RSEM may filter some alignable low quality reads,
648 the sum of expected counts for all transcript are generally less than
649 the total number of reads aligned.
650
651 'TPM' stands for Transcripts Per Million. It is a relative measure of
652 transcript abundance. The sum of all transcripts' TPM is 1
653 million. 'FPKM' stands for Fragments Per Kilobase of transcript per
654 Million mapped reads. It is another relative measure of transcript
655 abundance. If we define l_bar be the mean transcript length in a
656 sample, which can be calculated as
657
658 l_bar = \sum_i TPM_i / 10^6 * effective_length_i (i goes through every transcript), 
659
660 the following equation is hold:
661
662 FPKM_i = 10^3 / l_bar * TPM_i.
663
664 We can see that the sum of FPKM is not a constant across samples.
665
666 'IsoPct' stands for isoform percentage. It is the percentage of this
667 transcript's abandunce over its parent gene's abandunce. If its parent
668 gene has only one isoform or the gene information is not provided,
669 this field will be set to 100.
670
671 'pme_expected_count', 'pme_TPM', 'pme_FPKM' are posterior mean
672 estimates calculated by RSEM's Gibbs sampler. 'IsoPct_from_pme_TPM' is
673 the isoform percentage calculated from 'pme_TPM' values.
674
675 'TPM_ci_lower_bound', 'TPM_ci_upper_bound', 'FPKM_ci_lower_bound' and
676 'FPKM_ci_upper_bound' are lower(l) and upper(u) bounds of 95%
677 credibility intervals for TPM and FPKM values. The bounds are
678 inclusive (i.e. [l, u]).
679
680 =item B<sample_name.genes.results>
681
682 File containing gene level expression estimates. The first line
683 contains column names separated by the tab character. The format of
684 each line in the rest of this file is:
685
686 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]
687
688 Fields are separated by the tab character. Fields within "[]" are only
689 presented if '--calc-ci' is set. 
690
691 'transcript_id(s)' is a comma-separated list of transcript_ids
692 belonging to this gene. If no gene information is provided, 'gene_id'
693 and 'transcript_id(s)' are identical (the 'transcript_id').
694
695 A gene's 'length' and 'effective_length' are
696 defined as the weighted average of its transcripts' lengths and
697 effective lengths (weighted by 'IsoPct'). A gene's abundance estimates
698 are just the sum of its transcripts' abundance estimates.
699
700 =item B<sample_name.transcript.bam, sample_name.transcript.sorted.bam and sample_name.transcript.sorted.bam.bai>
701
702 Only generated when --no-bam-output is not specified.
703
704 'sample_name.transcript.bam' is a BAM-formatted file of read
705 alignments in transcript coordinates. The MAPQ field of each alignment
706 is set to min(100, floor(-10 * log10(1.0 - w) + 0.5)), where w is the
707 posterior probability of that alignment being the true mapping of a
708 read.  In addition, RSEM pads a new tag ZW:f:value, where value is a
709 single precision floating number representing the posterior
710 probability. Because this file contains all alignment lines produced
711 by bowtie or user-specified aligners, it can also be used as a
712 replacement of the aligner generated BAM/SAM file. For paired-end
713 reads, if one mate has alignments but the other does not, this file
714 marks the alignable mate as "unmappable" (flag bit 0x4) and appends an
715 optional field "Z0:A:!".
716
717 'sample_name.transcript.sorted.bam' and
718 'sample_name.transcript.sorted.bam.bai' are the sorted BAM file and
719 indices generated by samtools (included in RSEM package).
720
721 =item B<sample_name.genome.bam, sample_name.genome.sorted.bam and sample_name.genome.sorted.bam.bai>
722
723 Only generated when --no-bam-output is not specified and --output-genome-bam is specified.
724
725 'sample_name.genome.bam' is a BAM-formatted file of read alignments in
726 genomic coordinates. Alignments of reads that have identical genomic
727 coordinates (i.e., alignments to different isoforms that share the
728 same genomic region) are collapsed into one alignment.  The MAPQ field
729 of each alignment is set to min(100, floor(-10 * log10(1.0 - w) +
730 0.5)), where w is the posterior probability of that alignment being
731 the true mapping of a read.  In addition, RSEM pads a new tag
732 ZW:f:value, where value is a single precision floating number
733 representing the posterior probability. If an alignment is spliced, a
734 XS:A:value tag is also added, where value is either '+' or '-'
735 indicating the strand of the transcript it aligns to.
736
737 'sample_name.genome.sorted.bam' and 'sample_name.genome.sorted.bam.bai' are the
738 sorted BAM file and indices generated by samtools (included in RSEM package).
739
740 =item B<sample_name.time>
741
742 Only generated when --time is specified.
743
744 It contains time (in seconds) consumed by aligning reads, estimating expression levels and calculating credibility intervals.
745
746 =item B<sample_name.stat>
747
748 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.
749
750 =back
751
752 =head1 EXAMPLES
753
754 Assume the path to the bowtie executables is in the user's PATH environment variable. Reference files are under '/ref' with name 'mouse_125'. 
755
756 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:
757
758  rsem-calculate-expression --phred64-quals \
759                            -p 8 \
760                            --output-genome-bam \
761                            /data/mmliver.fq \
762                            /ref/mouse_125 \
763                            mmliver_single_quals
764
765 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:
766
767  rsem-calculate-expression -p 8 \
768                            --paired-end \
769                            /data/mmliver_1.fq \
770                            /data/mmliver_2.fq \
771                            /ref/mouse_125 \
772                            mmliver_paired_end_quals
773
774 3) '/data/mmliver.fa', single-end reads without quality scores. We want to use 8 threads:
775
776  rsem-calculate-expression -p 8 \
777                            --no-qualities \
778                            /data/mmliver.fa \
779                            /ref/mouse_125 \
780                            mmliver_single_without_quals
781
782 4) Data are the same as 1). 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:
783
784  rsem-calculate-expression --bowtie-path /sw/bowtie \
785                            --phred64-quals \
786                            --fragment-length-mean 150.0 \
787                            --fragment-length-sd 35.0 \
788                            -p 8 \
789                            --output-genome-bam \
790                            --calc-ci \
791                            --ci-memory 1024 \
792                            /data/mmliver.fq \
793                            /ref/mouse_125 \
794                            mmliver_single_quals
795
796 5) '/data/mmliver_paired_end_quals.bam', paired-end reads with quality scores.  We want to use 8 threads:
797
798  rsem-calculate-expression --paired-end \
799                            --bam \
800                            -p 8 \
801                            /data/mmliver_paired_end_quals.bam \
802                            /ref/mouse_125 \
803                            mmliver_paired_end_quals
804
805 =cut