#!/usr/bin/env perl use Getopt::Long; use Pod::Usage; use FindBin; use lib $FindBin::RealBin; use rsem_perl_utils qw(runCommand collectResults showVersionInfo); use Env qw(@PATH); @PATH = ($FindBin::RealBin, "$FindBin::RealBin/sam", @PATH); use strict; #const my $BURNIN = 200; my $NCV = 1000; my $SAMPLEGAP = 1; my $CONFIDENCE = 0.95; my $NSPC = 50; my $NMB = 1024; # default my $SortMem = "1G"; # default as 1G per thread my $status = 0; my $read_type = 1; # default, single end with qual my $bowtie_path = ""; my $C = 2; my $E = 99999999; my $L = 25; my $maxHits = 200; my $chunkMbs = 0; # 0 = use bowtie default my $phred33 = 0; my $phred64 = 0; my $solexa = 0; my $is_sam = 0; my $is_bam = 0; my $fn_list = ""; my $tagName = "XM"; my $probF = 0.5; my $minL = 1; my $maxL = 1000; my $mean = -1; my $sd = 0; my $estRSPD = 0; my $B = 20; my $nThreads = 1; my $genBamF = 1; # default is generating transcript bam file my $genGenomeBamF = 0; my $sampling = 0; my $calcPME = 0; my $calcCI = 0; my $quiet = 0; my $help = 0; my $paired_end = 0; my $no_qual = 0; my $keep_intermediate_files = 0; my $strand_specific = 0; my $bowtie2 = 0; my $bowtie2_path = ""; my $bowtie2_mismatch_rate = 0.1; my $bowtie2_k = 200; my $bowtie2_sensitivity_level = "sensitive"; # must be one of "very_fast", "fast", "sensitive", "very_sensitive" my $seed = "NULL"; my $version = 0; my $mTime = 0; my ($time_start, $time_end, $time_alignment, $time_rsem, $time_ci) = (0, 0, 0, 0, 0); my $mate1_list = ""; my $mate2_list = ""; my $inpF = ""; my ($refName, $sampleName, $sampleToken, $temp_dir, $stat_dir, $imdName, $statName) = (); my $gap = 32; my $alleleS = 0; GetOptions("keep-intermediate-files" => \$keep_intermediate_files, "temporary-folder=s" => \$temp_dir, "no-qualities" => \$no_qual, "paired-end" => \$paired_end, "strand-specific" => \$strand_specific, "sam" => \$is_sam, "bam" => \$is_bam, "sam-header-info=s" => \$fn_list, "tag=s" => \$tagName, "seed-length=i" => \$L, "bowtie-path=s" => \$bowtie_path, "bowtie-n=i" => \$C, "bowtie-e=i" => \$E, "bowtie-m=i" => \$maxHits, "bowtie-chunkmbs=i" => \$chunkMbs, "phred33-quals" => \$phred33, "phred64-quals" => \$phred64, #solexa1.3-quals" => \$phred64, "solexa-quals" => \$solexa, "bowtie2" => \$bowtie2, "bowtie2-path=s" => \$bowtie2_path, "bowtie2-mismatch-rate=f" => \$bowtie2_mismatch_rate, "bowtie2-k=i" => \$bowtie2_k, "bowtie2-sensitivity-level=s" => \$bowtie2_sensitivity_level, "forward-prob=f" => \$probF, "fragment-length-min=i" => \$minL, "fragment-length-max=i" => \$maxL, "fragment-length-mean=f" => \$mean, "fragment-length-sd=f" => \$sd, "estimate-rspd" => \$estRSPD, "num-rspd-bins=i" => \$B, "p|num-threads=i" => \$nThreads, "no-bam-output" => sub { $genBamF = 0; }, "output-genome-bam" => \$genGenomeBamF, "sampling-for-bam" => \$sampling, "calc-pme" => \$calcPME, "calc-ci" => \$calcCI, "ci-memory=i" => \$NMB, "samtools-sort-mem=s" => \$SortMem, "seed=i" => \$seed, "time" => \$mTime, "version" => \$version, "q|quiet" => \$quiet, "h|help" => \$help) or pod2usage(-exitval => 2, -verbose => 2); pod2usage(-verbose => 2) if ($help == 1); &showVersionInfo($FindBin::RealBin) if ($version == 1); #check parameters and options if ($is_sam || $is_bam) { pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (scalar(@ARGV) != 3); pod2usage(-msg => "--sam and --bam cannot be active at the same time!", -exitval => 2, -verbose => 2) if ($is_sam == 1&& $is_bam == 1); 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"); } else { pod2usage(-msg => "Invalid number of arguments!", -exitval => 2, -verbose => 2) if (!$paired_end && scalar(@ARGV) != 3 || $paired_end && scalar(@ARGV) != 4); 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)); pod2usage(-msg => "Only one of --phred33-quals, --phred64-quals, and --solexa-quals can be active!", -exitval => 2, -verbose => 2) if ($phred33 + $phred64 + $solexa > 1); 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 ""); 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")); 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)); 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)); 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"))); } pod2usage(-msg => "Forward probability should be in [0, 1]!", -exitval => 2, -verbose => 2) if ($probF < 0 || $probF > 1); pod2usage(-msg => "Min fragment length should be at least 1!", -exitval => 2, -verbose => 2) if ($minL < 1); pod2usage(-msg => "Min fragment length should be smaller or equal to max fragment length!", -exitval => 2, -verbose => 2) if ($minL > $maxL); pod2usage(-msg => "The memory allocated for calculating credibility intervals should be at least 1 MB!\n", -exitval => 2, -verbose => 2) if ($NMB < 1); pod2usage(-msg => "Number of threads should be at least 1!\n", -exitval => 2, -verbose => 2) if ($nThreads < 1); pod2usage(-msg => "Seed length should be at least 5!\n", -exitval => 2, -verbose => 2) if ($L < 5); pod2usage(-msg => "--sampling-for-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($sampling && !$genBamF); pod2usage(-msg => "--output-genome-bam cannot be specified if --no-bam-output is specified!\n", -exitval => 2, -verbose => 2) if ($genGenomeBamF && !$genBamF); 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)); 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"; } if ($strand_specific) { $probF = 1.0; } if ($paired_end) { if ($no_qual) { $read_type = 2; } else { $read_type = 3; } } else { if ($no_qual) { $read_type = 0; } else { $read_type = 1; } } if (scalar(@ARGV) == 3) { if ($is_sam || $is_bam) { $inpF = $ARGV[0]; } else {$mate1_list = $ARGV[0]; } $refName = $ARGV[1]; $sampleName = $ARGV[2]; } else { $mate1_list = $ARGV[0]; $mate2_list = $ARGV[1]; $refName = $ARGV[2]; $sampleName = $ARGV[3]; } if (((-e "$refName.ta") && !(-e "$refName.gt")) || (!(-e "$refName.ta") && (-e "$refName.gt"))) { print "Allele-specific expression related reference files are corrupted!\n"; exit(-1); } $alleleS = (-e "$refName.ta") && (-e "$refName.gt"); if ($genGenomeBamF) { open(INPUT, "$refName.ti"); my $line = ; chomp($line); close(INPUT); my ($M, $type) = split(/ /, $line); pod2usage(-msg => "No genome information provided, so genome bam file cannot be generated!\n", -exitval => 2, -verbose => 2) if ($type != 0); } my $pos = rindex($sampleName, '/'); if ($pos < 0) { $sampleToken = $sampleName; } else { $sampleToken = substr($sampleName, $pos + 1); } if ($temp_dir eq "") { $temp_dir = "$sampleName.temp"; } $stat_dir = "$sampleName.stat"; if (!(-d $temp_dir) && !mkdir($temp_dir)) { print "Fail to create folder $temp_dir.\n"; exit(-1); } if (!(-d $stat_dir) && !mkdir($stat_dir)) { print "Fail to create folder $stat_dir.\n"; exit(-1); } $imdName = "$temp_dir/$sampleToken"; $statName = "$stat_dir/$sampleToken"; if (!$is_sam && !$is_bam && !$no_qual && ($phred33 + $phred64 + $solexa == 0)) { $phred33 = 1; } my ($mate_minL, $mate_maxL) = (1, $maxL); if ($bowtie_path ne "") { $bowtie_path .= "/"; } if ($bowtie2_path ne "") { $bowtie2_path .= "/"; } my $command = ""; if (!$is_sam && !$is_bam) { if (!$bowtie2) { $command = $bowtie_path."bowtie"; if ($no_qual) { $command .= " -f"; } else { $command .= " -q"; } if ($phred33) { $command .= " --phred33-quals"; } elsif ($phred64) { $command .= " --phred64-quals"; } elsif ($solexa) { $command .= " --solexa-quals"; } $command .= " -n $C -e $E -l $L"; if ($read_type == 2 || $read_type == 3) { $command .= " -I $minL -X $maxL"; } if ($chunkMbs > 0) { $command .= " --chunkmbs $chunkMbs"; } if ($strand_specific || $probF == 1.0) { $command .= " --norc"; } elsif ($probF == 0.0) { $command .= " --nofw"; } $command .= " -p $nThreads -a -m $maxHits -S"; if ($quiet) { $command .= " --quiet"; } $command .= " $refName"; if ($read_type == 0 || $read_type == 1) { $command .= " $mate1_list"; } else { $command .= " -1 $mate1_list -2 $mate2_list"; } # pipe to samtools to generate a BAM file $command .= " | samtools view -S -b -o $imdName.bam -"; } else { $command = $bowtie2_path."bowtie2"; if ($no_qual) { $command .= " -f"; } else { $command .= " -q"; } if ($phred33) { $command .= " --phred33"; } elsif ($phred64) { $command .= " --phred64"; } elsif ($solexa) { $command .= " --solexa-quals"; } if ($bowtie2_sensitivity_level eq "very_fast") { $command .= " --very-fast"; } elsif ($bowtie2_sensitivity_level eq "fast") { $command .= " --fast"; } elsif ($bowtie2_sensitivity_level eq "sensitive") { $command .= " --sensitive"; } else { $command .= " --very-sensitive"; } $command .= " --dpad 0 --gbar 99999999 --mp 1,1 --np 1 --score-min L,0,-$bowtie2_mismatch_rate"; if ($read_type == 2 || $read_type == 3) { $command .= " -I $minL -X $maxL --no-mixed --no-discordant"; } if ($strand_specific || $probF == 1.0) { $command .= " --norc"; } elsif ($probF == 0.0) { $command .= " --nofw"; } $command .= " -p $nThreads -k $bowtie2_k"; if ($quiet) { $command .= " --quiet"; } $command .= " -x $refName"; if ($read_type == 0 || $read_type == 1) { $command .= " -U $mate1_list"; } else { $command .= " -1 $mate1_list -2 $mate2_list"; } # pipe to samtools to generate a BAM file $command .= " | samtools view -S -b -o $imdName.bam -"; } if ($mTime) { $time_start = time(); } &runCommand($command); if ($mTime) { $time_end = time(); $time_alignment = $time_end - $time_start; } $inpF = "$imdName.bam"; $is_bam = 1; # alignments are outputed as a BAM file } if ($mTime) { $time_start = time(); } $command = "rsem-parse-alignments $refName $imdName $statName"; my $samInpType; if ($is_sam) { $samInpType = "s"; } elsif ($is_bam) { $samInpType = "b"; } $command .= " $samInpType $inpF -t $read_type"; if ($fn_list ne "") { $command .= " -l $fn_list"; } if ($tagName ne "") { $command .= " -tag $tagName"; } if ($quiet) { $command .= " -q"; } &runCommand($command); $command = "rsem-build-read-index $gap"; if ($read_type == 0) { $command .= " 0 $quiet $imdName\_alignable.fa"; } elsif ($read_type == 1) { $command .= " 1 $quiet $imdName\_alignable.fq"; } elsif ($read_type == 2) { $command .= " 0 $quiet $imdName\_alignable_1.fa $imdName\_alignable_2.fa"; } elsif ($read_type == 3) { $command .= " 1 $quiet $imdName\_alignable_1.fq $imdName\_alignable_2.fq"; } else { print "Impossible! read_type is not in [1,2,3,4]!\n"; exit(-1); } &runCommand($command); my $doesOpen = open(OUTPUT, ">$imdName.mparams"); if ($doesOpen == 0) { print "Cannot generate $imdName.mparams!\n"; exit(-1); } print OUTPUT "$minL $maxL\n"; print OUTPUT "$probF\n"; print OUTPUT "$estRSPD\n"; print OUTPUT "$B\n"; print OUTPUT "$mate_minL $mate_maxL\n"; print OUTPUT "$mean $sd\n"; print OUTPUT "$L\n"; close(OUTPUT); my @seeds = (); if ($seed ne "NULL") { srand($seed); for (my $i = 0; $i < 3; $i++) { push(@seeds, int(rand(1 << 32))); } } $command = "rsem-run-em $refName $read_type $sampleName $imdName $statName -p $nThreads"; if ($genBamF) { $command .= " -b $samInpType $inpF"; if ($fn_list ne "") { $command .= " 1 $fn_list"; } else { $command .= " 0"; } if ($sampling) { $command .= " --sampling"; } if ($seed ne "NULL") { $command .= " --seed $seeds[0]"; } } if ($calcPME || $calcCI) { $command .= " --gibbs-out"; } if ($quiet) { $command .= " -q"; } &runCommand($command); if ($alleleS) { &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } else { &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } if ($genBamF) { $command = "samtools sort -@ $nThreads -m $SortMem $sampleName.transcript.bam $sampleName.transcript.sorted"; &runCommand($command); $command = "samtools index $sampleName.transcript.sorted.bam"; &runCommand($command); if ($genGenomeBamF) { $command = "rsem-tbam2gbam $refName $sampleName.transcript.bam $sampleName.genome.bam"; &runCommand($command); $command = "samtools sort -@ $nThreads -m $SortMem $sampleName.genome.bam $sampleName.genome.sorted"; &runCommand($command); $command = "samtools index $sampleName.genome.sorted.bam"; &runCommand($command); } } if ($mTime) { $time_end = time(); $time_rsem = $time_end - $time_start; } if ($mTime) { $time_start = time(); } if ($calcPME || $calcCI ) { $command = "rsem-run-gibbs $refName $imdName $statName $BURNIN $NCV $SAMPLEGAP"; $command .= " -p $nThreads"; if ($seed ne "NULL") { $command .= " --seed $seeds[1]"; } if ($quiet) { $command .= " -q"; } &runCommand($command); } if ($calcPME || $calcCI) { if ($alleleS) { system("mv $sampleName.alleles.results $imdName.alleles.results.bak1"); system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak1"); system("mv $sampleName.genes.results $imdName.genes.results.bak1"); &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } else { system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak1"); system("mv $sampleName.genes.results $imdName.genes.results.bak1"); &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } } if ($calcCI) { $command = "rsem-calculate-credibility-intervals $refName $imdName $statName $CONFIDENCE $NCV $NSPC $NMB"; $command .= " -p $nThreads"; if ($seed ne "NULL") { $command .= " --seed $seeds[2]"; } if ($quiet) { $command .= " -q"; } &runCommand($command); if ($alleleS) { system("mv $sampleName.alleles.results $imdName.alleles.results.bak2"); system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak2"); system("mv $sampleName.genes.results $imdName.genes.results.bak2"); &collectResults("allele", "$imdName.allele_res", "$sampleName.alleles.results"); # allele level &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } else { system("mv $sampleName.isoforms.results $imdName.isoforms.results.bak2"); system("mv $sampleName.genes.results $imdName.genes.results.bak2"); &collectResults("isoform", "$imdName.iso_res", "$sampleName.isoforms.results"); # isoform level &collectResults("gene", "$imdName.gene_res", "$sampleName.genes.results"); # gene level } } if ($mTime) { $time_end = time(); $time_ci = $time_end - $time_start; } if ($mTime) { $time_start = time(); } if (!$keep_intermediate_files) { &runCommand("rm -rf $temp_dir", "Fail to delete the temporary folder!"); } if ($mTime) { $time_end = time(); } if ($mTime) { open(OUTPUT, ">$sampleName.time"); print OUTPUT "Aligning reads: $time_alignment s.\n"; print OUTPUT "Estimating expression levels: $time_rsem s.\n"; print OUTPUT "Calculating credibility intervals: $time_ci s.\n"; # my $time_del = $time_end - $time_start; # print OUTPUT "Delete: $time_del s.\n"; close(OUTPUT); } __END__ =head1 NAME rsem-calculate-expression =head1 SYNOPSIS rsem-calculate-expression [options] upstream_read_file(s) reference_name sample_name rsem-calculate-expression [options] --paired-end upstream_read_file(s) downstream_read_file(s) reference_name sample_name rsem-calculate-expression [options] --sam/--bam [--paired-end] input reference_name sample_name =head1 ARGUMENTS =over =item B 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. =item B 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. =item B 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. =item B The name of the reference used. The user must have run 'rsem-prepare-reference' with this reference_name before running this program. =item B The name of the sample analyzed. All output files are prefixed by this name (e.g., sample_name.genes.results) =back =head1 OPTIONS =over =item B<--paired-end> Input reads are paired-end reads. (Default: off) =item B<--no-qualities> Input reads do not contain quality scores. (Default: off) =item B<--strand-specific> 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) =item B<--sam> Input file is in SAM format. (Default: off) =item B<--bam> Input file is in BAM format. (Default: off) =item B<--sam-header-info> 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: "") =item B<-p/--num-threads> Number of threads to use. Both Bowtie/Bowtie2, expression estimation and 'samtools sort' will use this many threads. (Default: 1) =item B<--no-bam-output> Do not output any BAM file. (Default: off) =item B<--output-genome-bam> 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) =item B<--sampling-for-bam> 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) =item B<--seed> 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) =item B<--calc-pme> Run RSEM's collapsed Gibbs sampler to calculate posterior mean estimates. (Default: off) =item B<--calc-ci> Calculate 95% credibility intervals and posterior mean estimates. (Default: off) =item B<--seed-length> 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) =item B<--tag> 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 :i:, where a bigger than 0 indicates a read with too many alignments. (Default: "") =item B<--bowtie-path> The path to the Bowtie executables. (Default: the path to the Bowtie executables is assumed to be in the user's PATH environment variable) =item B<--bowtie-n> (Bowtie parameter) max # of mismatches in the seed. (Range: 0-3, Default: 2) =item B<--bowtie-e> (Bowtie parameter) max sum of mismatch quality scores across the alignment. (Default: 99999999) =item B<--bowtie-m> (Bowtie parameter) suppress all alignments for a read if > valid alignments exist. (Default: 200) =item B<--bowtie-chunkmbs> (Bowtie parameter) memory allocated for best first alignment calculation (Default: 0 - use Bowtie's default) =item B<--phred33-quals> Input quality scores are encoded as Phred+33. (Default: on) =item B<--phred64-quals> Input quality scores are encoded as Phred+64 (default for GA Pipeline ver. >= 1.3). (Default: off) =item B<--solexa-quals> Input quality scores are solexa encoded (from GA Pipeline ver. < 1.3). (Default: off) =item B<--bowtie2> 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) =item B<--bowtie2-path> (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) =item B<--bowtie2-mismatch-rate> (Bowtie 2 parameter) The maximum mismatch rate allowed. (Default: 0.1) =item B<--bowtie2-k> (Bowtie 2 parameter) Find up to alignments per read. (Default: 200) =item B<--bowtie2-sensitivity-level> (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. 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) =item B<--forward-prob> 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) =item B<--fragment-length-min> Minimum read/insert length allowed. This is also the value for the Bowtie/Bowtie2 -I option. (Default: 1) =item B<--fragment-length-max> Maximum read/insert length allowed. This is also the value for the Bowtie/Bowtie 2 -X option. (Default: 1000) =item B<--fragment-length-mean> (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) =item B<--fragment-length-sd> (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>) =item B<--estimate-rspd> 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) =item B<--num-rspd-bins> Number of bins in the RSPD. Only relevant when '--estimate-rspd' is specified. Use of the default setting is recommended. (Default: 20) =item B<--ci-memory> 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) =item B<--samtools-sort-mem> Set the maximum memory per thread that can be used by 'samtools sort'. represents the memory and accepts suffices 'K/M/G'. RSEM will pass to the '-m' option of 'samtools sort'. Please note that the default used here is different from the default used by samtools. (Default: 1G) =item B<--keep-intermediate-files> 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) =item B<--temporary-folder> 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) =item B<--time> Output time consumed by each step of RSEM to 'sample_name.time'. (Default: off) =item B<-q/--quiet> Suppress the output of logging information. (Default: off) =item B<-h/--help> Show help information. =item B<--version> Show version information. =back =head1 DESCRIPTION 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. 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: convert-sam-for-rsem /ref/mouse_125 input.sam -o input_for_rsem.sam For details, please refer to 'convert-sam-for-rsem's documentation page. 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 '*'. The user must run 'rsem-prepare-reference' with the appropriate reference before using this program. 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. Please note that some of the default values for the Bowtie parameters are not the same as those defined for Bowtie itself. The temporary directory and all intermediate files will be removed when RSEM finishes unless '--keep-intermediate-files' is specified. With the '--calc-ci' option, 95% credibility intervals and posterior mean estimates will be calculated in addition to maximum likelihood estimates. =head1 OUTPUT =over =item B File containing isoform level expression estimates. The first line contains column names separated by the tab character. The format of each line in the rest of this file is: 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] Fields are separated by the tab character. Fields within "[]" are optional. They will not be presented if neither '--calc-pme' nor '--calc-ci' is set. 'transcript_id' is the transcript name of this transcript. 'gene_id' is the gene name of the gene which this transcript belongs to (denote this gene as its parent gene). If no gene information is provided, 'gene_id' and 'transcript_id' are the same. 'length' is this transcript's sequence length (poly(A) tail is not counted). 'effective_length' counts only the positions that can generate a valid fragment. If no poly(A) tail is added, 'effective_length' is equal to transcript length - mean fragment length + 1. If one transcript's effective length is less than 1, this transcript's both effective length and abundance estimates are set to 0. 'expected_count' is the sum of the posterior probability of each read comes from this transcript over all reads. Because 1) each read aligning to this transcript has a probability of being generated from background noise; 2) RSEM may filter some alignable low quality reads, the sum of expected counts for all transcript are generally less than the total number of reads aligned. 'TPM' stands for Transcripts Per Million. It is a relative measure of transcript abundance. The sum of all transcripts' TPM is 1 million. 'FPKM' stands for Fragments Per Kilobase of transcript per Million mapped reads. It is another relative measure of transcript abundance. If we define l_bar be the mean transcript length in a sample, which can be calculated as l_bar = \sum_i TPM_i / 10^6 * effective_length_i (i goes through every transcript), the following equation is hold: FPKM_i = 10^3 / l_bar * TPM_i. We can see that the sum of FPKM is not a constant across samples. 'IsoPct' stands for isoform percentage. It is the percentage of this transcript's abandunce over its parent gene's abandunce. If its parent gene has only one isoform or the gene information is not provided, this field will be set to 100. 'posterior_mean_count', 'pme_TPM', 'pme_FPKM' are posterior mean estimates calculated by RSEM's Gibbs sampler. 'posterior_standard_deviation_of_count' is the posterior standard deviation of counts. 'IsoPct_from_pme_TPM' is the isoform percentage calculated from 'pme_TPM' values. 'TPM_ci_lower_bound', 'TPM_ci_upper_bound', 'FPKM_ci_lower_bound' and 'FPKM_ci_upper_bound' are lower(l) and upper(u) bounds of 95% credibility intervals for TPM and FPKM values. The bounds are inclusive (i.e. [l, u]). =item B File containing gene level expression estimates. The first line contains column names separated by the tab character. The format of each line in the rest of this file is: 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] Fields are separated by the tab character. Fields within "[]" are optional. They will not be presented if neither '--calc-pme' nor '--calc-ci' is set. 'transcript_id(s)' is a comma-separated list of transcript_ids belonging to this gene. If no gene information is provided, 'gene_id' and 'transcript_id(s)' are identical (the 'transcript_id'). A gene's 'length' and 'effective_length' are defined as the weighted average of its transcripts' lengths and effective lengths (weighted by 'IsoPct'). A gene's abundance estimates are just the sum of its transcripts' abundance estimates. =item B Only generated when the RSEM references are built with allele-specific transcripts. This file contains allele level expression estimates for allele-specific expression calculation. The first line contains column names separated by the tab character. The format of each line in the rest of this file is: 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] Fields are separated by the tab character. Fields within "[]" are optional. They will not be presented if neither '--calc-pme' nor '--calc-ci' is set. 'allele_id' is the allele-specific name of this allele-specific transcript. 'AlleleIsoPct' stands for allele-specific percentage on isoform level. It is the percentage of this allele-specific transcript's abundance over its parent transcript's abundance. If its parent transcript has only one allele variant form, this field will be set to 100. 'AlleleGenePct' stands for allele-specific percentage on gene level. It is the percentage of this allele-specific transcript's abundance over its parent gene's abundance. 'AlleleIsoPct_from_pme_TPM' and 'AlleleGenePct_from_pme_TPM' have similar meanings. They are calculated based on posterior mean estimates. Please note that if this file is present, the fields 'length' and 'effective_length' in 'sample_name.isoforms.results' should be interpreted similarly as the corresponding definitions in 'sample_name.genes.results'. =item B Only generated when --no-bam-output is not specified. 'sample_name.transcript.bam' is a BAM-formatted file of read alignments in transcript coordinates. The MAPQ field of each alignment is set to min(100, floor(-10 * log10(1.0 - w) + 0.5)), where w is the posterior probability of that alignment being the true mapping of a read. In addition, RSEM pads a new tag ZW:f:value, where value is a single precision floating number representing the posterior probability. Because this file contains all alignment lines produced by bowtie or user-specified aligners, it can also be used as a replacement of the aligner generated BAM/SAM file. For paired-end reads, if one mate has alignments but the other does not, this file marks the alignable mate as "unmappable" (flag bit 0x4) and appends an optional field "Z0:A:!". 'sample_name.transcript.sorted.bam' and 'sample_name.transcript.sorted.bam.bai' are the sorted BAM file and indices generated by samtools (included in RSEM package). =item B Only generated when --no-bam-output is not specified and --output-genome-bam is specified. 'sample_name.genome.bam' is a BAM-formatted file of read alignments in genomic coordinates. Alignments of reads that have identical genomic coordinates (i.e., alignments to different isoforms that share the same genomic region) are collapsed into one alignment. The MAPQ field of each alignment is set to min(100, floor(-10 * log10(1.0 - w) + 0.5)), where w is the posterior probability of that alignment being the true mapping of a read. In addition, RSEM pads a new tag ZW:f:value, where value is a single precision floating number representing the posterior probability. If an alignment is spliced, a XS:A:value tag is also added, where value is either '+' or '-' indicating the strand of the transcript it aligns to. 'sample_name.genome.sorted.bam' and 'sample_name.genome.sorted.bam.bai' are the sorted BAM file and indices generated by samtools (included in RSEM package). =item B Only generated when --time is specified. It contains time (in seconds) consumed by aligning reads, estimating expression levels and calculating credibility intervals. =item B 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. =back =head1 EXAMPLES Assume the path to the bowtie executables is in the user's PATH environment variable. Reference files are under '/ref' with name 'mouse_125'. 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: rsem-calculate-expression --phred64-quals \ -p 8 \ --output-genome-bam \ /data/mmliver.fq \ /ref/mouse_125 \ mmliver_single_quals 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: rsem-calculate-expression -p 8 \ --paired-end \ /data/mmliver_1.fq \ /data/mmliver_2.fq \ /ref/mouse_125 \ mmliver_paired_end_quals 3) '/data/mmliver.fa', single-end reads without quality scores. We want to use 8 threads: rsem-calculate-expression -p 8 \ --no-qualities \ /data/mmliver.fa \ /ref/mouse_125 \ mmliver_single_without_quals 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: rsem-calculate-expression --bowtie-path /sw/bowtie \ --phred64-quals \ --fragment-length-mean 150.0 \ --fragment-length-sd 35.0 \ -p 8 \ --output-genome-bam \ --calc-ci \ --ci-memory 1024 \ /data/mmliver.fq \ /ref/mouse_125 \ mmliver_single_quals 5) '/data/mmliver_paired_end_quals.bam', paired-end reads with quality scores. We want to use 8 threads: rsem-calculate-expression --paired-end \ --bam \ -p 8 \ /data/mmliver_paired_end_quals.bam \ /ref/mouse_125 \ mmliver_paired_end_quals =cut