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