Contents

Note: the most recent version of this tutorial can be found here.

1 Introduction

systemPipeR provides utilities for building analysis workflows with automated report generation for next generation sequence (NGS) applications such as RNA-Seq, ChIP-Seq, VAR-Seq and many others (Girke 2014). An important feature is support for running command-line software, such as NGS aligners, on both single machines or compute clusters. This includes both interactive job submissions or batch submissions to queuing systems of clusters. For instance, systemPipeR can be used with most command-line aligners such as BWA (Heng Li 2013; H Li and Durbin 2009), TopHat2 (Kim et al. 2013) and Bowtie2 (Langmead and Salzberg 2012), as well as the R-based NGS aligners Rsubread (Liao, Smyth, and Shi 2013) and gsnap (gmapR) (Wu and Nacu 2010). Efficient handling of complex sample sets and experimental designs is facilitated by a well-defined sample annotation infrastructure which improves reproducibility and user-friendliness of many typical analysis workflows in the NGS area (Lawrence et al. 2013).

A central concept for designing workflows within the sytemPipeR environment is the use of sample management containers called SYSargs. Instances of this S4 object class are constructed by the systemArgs function from two simple tablular files: a targets file and a param file. The latter is optional for workflow steps lacking command-line software. Typically, a SYSargs instance stores all sample-level inputs as well as the paths to the corresponding outputs generated by command-line- or R-based software generating sample-level output files, such as read preprocessors (trimmed/filtered FASTQ files), aligners (SAM/BAM files), variant callers (VCF/BCF files) or peak callers (BED/WIG files). Each sample level input/outfile operation uses its own SYSargs instance. The outpaths of SYSargs usually define the sample inputs for the next SYSargs instance. This connectivity is established by writing the outpaths with the writeTargetsout function to a new targets file that serves as input to the next systemArgs call. By chaining several SYSargs steps together one can construct complex workflows involving many sample-level input/output file operations with any combinaton of command-line or R-based software.

SystemPipeR_Workflow

The intended way of running sytemPipeR workflows is via *.Rnw or *.Rmd files, which can be executed either line-wise in interactive mode or with a single command from R or the command-line using a Makefile. This way comprehensive and reproducible analysis reports in PDF or HTML format can be generated in a fully automated manner. Templates for setting up custom project reports are provided as *.Rnw files in the vignettes subdirectory of this package. The corresponding PDFs of these report templates are linked here: systemPipeRNAseq, systemPipeChIPseq and systemPipeVARseq. To work with *.Rnw or *.Rmd files efficiently, basic knowledge of Sweave or knitr and Latex or R Markdown v2 is required.

2 Getting Started

2.1 Installation

The R software for running systemPipeR and systemPipeRdata can be downloaded from CRAN. The systemPipeR environment can be installed from R using the biocLite install command.

source("http://bioconductor.org/biocLite.R") # Sources the biocLite.R installation script 
biocLite("systemPipeR") # Installs systemPipeR from Bioconductor 
biocLite("tgirke/systemPipeRdata", build_vignettes=TRUE, dependencies=TRUE) # From github

2.2 Loading package and documentation

library("systemPipeR") # Loads the package
library(help="systemPipeR") # Lists package info
vignette("systemPipeR") # Opens vignette

2.3 Sample FASTQ files

The mini sample FASTQ files used by this overview vignette as well as the associated workflow reporting vignettes can be downloaded from here. The chosen data set SRP010938 contains 18 paired-end (PE) read sets from Arabidposis thaliana (Howard et al. 2013). To minimize processing time during testing, each FASTQ file has been subsetted to 90,000-100,000 randomly sampled PE reads that map to the first 100,000 nucleotides of each chromosome of the A. thalina genome. The corresponding reference genome sequence (FASTA) and its GFF annotion files (provided in the same download) have been truncated accordingly. This way the entire test sample data set is less than 200MB in storage space. A PE read set has been chosen for this test data set for flexibility, because it can be used for testing both types of analysis routines requiring either SE (single end) reads or PE reads.

2.4 Structure of targets file

The targets file defines all input files (e.g. FASTQ, BAM, BCF) and sample comparisons of an analysis workflow. The following shows the format of a sample targets file provided by this package. In a target file with a single type of input files, here FASTQ files of single end (SE) reads, the first three columns are mandatory including their column names, while it is four mandatory columns for FASTQ files for PE reads. All subsequent columns are optional and any number of additional columns can be added as needed.

library(systemPipeR)
targetspath <- system.file("extdata", "targets.txt", package="systemPipeR")
read.delim(targetspath, comment.char = "#")
##                    FileName SampleName Factor SampleLong Experiment        Date
## 1  ./data/SRR446027_1.fastq        M1A     M1  Mock.1h.A          1 23-Mar-2012
## 2  ./data/SRR446028_1.fastq        M1B     M1  Mock.1h.B          1 23-Mar-2012
## 3  ./data/SRR446029_1.fastq        A1A     A1   Avr.1h.A          1 23-Mar-2012
## 4  ./data/SRR446030_1.fastq        A1B     A1   Avr.1h.B          1 23-Mar-2012
## 5  ./data/SRR446031_1.fastq        V1A     V1   Vir.1h.A          1 23-Mar-2012
## 6  ./data/SRR446032_1.fastq        V1B     V1   Vir.1h.B          1 23-Mar-2012
## 7  ./data/SRR446033_1.fastq        M6A     M6  Mock.6h.A          1 23-Mar-2012
## 8  ./data/SRR446034_1.fastq        M6B     M6  Mock.6h.B          1 23-Mar-2012
## 9  ./data/SRR446035_1.fastq        A6A     A6   Avr.6h.A          1 23-Mar-2012
## 10 ./data/SRR446036_1.fastq        A6B     A6   Avr.6h.B          1 23-Mar-2012
## 11 ./data/SRR446037_1.fastq        V6A     V6   Vir.6h.A          1 23-Mar-2012
## 12 ./data/SRR446038_1.fastq        V6B     V6   Vir.6h.B          1 23-Mar-2012
## 13 ./data/SRR446039_1.fastq       M12A    M12 Mock.12h.A          1 23-Mar-2012
## 14 ./data/SRR446040_1.fastq       M12B    M12 Mock.12h.B          1 23-Mar-2012
## 15 ./data/SRR446041_1.fastq       A12A    A12  Avr.12h.A          1 23-Mar-2012
## 16 ./data/SRR446042_1.fastq       A12B    A12  Avr.12h.B          1 23-Mar-2012
## 17 ./data/SRR446043_1.fastq       V12A    V12  Vir.12h.A          1 23-Mar-2012
## 18 ./data/SRR446044_1.fastq       V12B    V12  Vir.12h.B          1 23-Mar-2012

2.5 Structure of targets file for paired end (PE) samples

targetspath <- system.file("extdata", "targetsPE.txt", package="systemPipeR")
read.delim(targetspath, comment.char = "#")[1:2,1:6]
##                  FileName1                FileName2 SampleName Factor SampleLong Experiment
## 1 ./data/SRR446027_1.fastq ./data/SRR446027_2.fastq        M1A     M1  Mock.1h.A          1
## 2 ./data/SRR446028_1.fastq ./data/SRR446028_2.fastq        M1B     M1  Mock.1h.B          1

2.6 Sample comparisons

Sample comparisons are defined in the header lines of the targets file starting with ‘# <CMP>’. The function readComp imports the comparison and stores them in a list. Alternatively, readComp can obtain the comparison information from the corresponding SYSargs object (see below). Note, the header lines are optional. They are mainly useful for controlling comparative analysis according to certain biological expectations, such as simple pairwise comparisons in RNA-Seq experiments.

readComp(file=targetspath, format="vector", delim="-")
## $CMPset1
## [1] "M1-A1"   "M1-V1"   "A1-V1"   "M6-A6"   "M6-V6"   "A6-V6"   "M12-A12" "M12-V12" "A12-V12"
## 
## $CMPset2
##  [1] "M1-A1"   "M1-V1"   "M1-M6"   "M1-A6"   "M1-V6"   "M1-M12"  "M1-A12"  "M1-V12"  "A1-V1"  
## [10] "A1-M6"   "A1-A6"   "A1-V6"   "A1-M12"  "A1-A12"  "A1-V12"  "V1-M6"   "V1-A6"   "V1-V6"  
## [19] "V1-M12"  "V1-A12"  "V1-V12"  "M6-A6"   "M6-V6"   "M6-M12"  "M6-A12"  "M6-V12"  "A6-V6"  
## [28] "A6-M12"  "A6-A12"  "A6-V12"  "V6-M12"  "V6-A12"  "V6-V12"  "M12-A12" "M12-V12" "A12-V12"

2.7 Structure of param file and SYSargs container

The param file defines the parameters of the command-line software. The following shows the format of a sample param file provided by this package.

parampath <- system.file("extdata", "tophat.param", package="systemPipeR")
read.delim(parampath, comment.char = "#")
##      PairSet         Name                                  Value
## 1    modules         <NA>                          bowtie2/2.1.0
## 2    modules         <NA>                          tophat/2.0.8b
## 3   software         <NA>                                 tophat
## 4      cores           -p                                      4
## 5      other         <NA> -g 1 --segment-length 25 -i 30 -I 3000
## 6   outfile1           -o                            <FileName1>
## 7   outfile1         path                             ./results/
## 8   outfile1       remove                                   <NA>
## 9   outfile1       append                                .tophat
## 10  outfile1 outextension              .tophat/accepted_hits.bam
## 11 reference         <NA>                    ./data/tair10.fasta
## 12   infile1         <NA>                            <FileName1>
## 13   infile1         path                                   <NA>
## 14   infile2         <NA>                            <FileName2>
## 15   infile2         path                                   <NA>

The systemArgs function imports the definitions of both the param file and the targets file, and stores all relevant information as SYSargs object. To run the pipeline without command-line software, one can assign NULL to sysma instead of a param file. In addition, one can start the systemPipeR workflow with pre-generated BAM files by providing a targets file where the FileName column gives the paths to the BAM files and sysma is assigned NULL.

args <- suppressWarnings(systemArgs(sysma=parampath, mytargets=targetspath))
args
## An instance of 'SYSargs' for running 'tophat' on 18 samples

Several accessor functions are available that are named after the slot names of the SYSargs object class.

names(args)
##  [1] "targetsin"     "targetsout"    "targetsheader" "modules"       "software"      "cores"        
##  [7] "other"         "reference"     "results"       "infile1"       "infile2"       "outfile1"     
## [13] "sysargs"       "outpaths"
modules(args)
## [1] "bowtie2/2.1.0" "tophat/2.0.8b"
cores(args)
## [1] 4
outpaths(args)[1]
##                                                                                                               M1A 
## "/tmp/RtmpJpRaUD/Rbuild2fe51fd5c47e/systemPipeRdata/vignettes/results/SRR446027_1.fastq.tophat/accepted_hits.bam"
sysargs(args)[1]
##                                                                                                                                                                                                                                                                                    M1A 
## "tophat -p 4 -g 1 --segment-length 25 -i 30 -I 3000 -o /tmp/RtmpJpRaUD/Rbuild2fe51fd5c47e/systemPipeRdata/vignettes/results/SRR446027_1.fastq.tophat /tmp/RtmpJpRaUD/Rbuild2fe51fd5c47e/systemPipeRdata/vignettes/data/tair10.fasta ./data/SRR446027_1.fastq ./data/SRR446027_2.fastq"

3 Workflow overview

3.1 Define environment settings and samples

Load package

library(systemPipeR)

Construct SYSargs object from param and targets files.

args <- systemArgs(sysma="trim.param", mytargets="targets.txt")

3.2 Read Preprocessing

The function preprocessReads allows to apply predefined or custom read preprocessing functions to all FASTQ files referenced in a SYSargs container, such as quality filtering or adaptor trimming routines. The paths to the resulting output FASTQ files are stored in the outpaths slot of the SYSargs object. Internally, preprocessReads uses the FastqStreamer function from the ShortRead package to stream through large FASTQ files in a memory-efficient manner. The following example performs adaptor trimming with the trimLRPatterns function from the Biostrings package. After the trimming step a new targets file is generated (here targets_trim.txt) containing the paths to the trimmed FASTQ files. The new targets file can be used for the next workflow step with an updated SYSargs instance, running the NGS alignments using the trimmed FASTQ files.

preprocessReads(args=args, Fct="trimLRPatterns(Rpattern='GCCCGGGTAA', subject=fq)", 
                batchsize=100000, overwrite=TRUE, compress=TRUE)
writeTargetsout(x=args, file="targets_trim.txt")

The following example shows how one can design a custom read preprocessing function using utilities provided by the ShortRead package, and then run it in batch mode with the ‘preprocessReads’ function (here on paired-end reads).

args <- systemArgs(sysma="trimPE.param", mytargets="targetsPE.txt")
filterFct <- function(fq, cutoff=20, Nexceptions=0) {
    qcount <- rowSums(as(quality(fq), "matrix") <= cutoff)
    fq[qcount <= Nexceptions] # Retains reads where Phred scores are >= cutoff with N exceptions
}
preprocessReads(args=args, Fct="filterFct(fq, cutoff=20, Nexceptions=0)", batchsize=100000)
writeTargetsout(x=args, file="targets_PEtrim.txt")

3.3 FASTQ quality report

The following seeFastq and seeFastqPlot functions generate and plot a series of useful quality statistics for a set of FASTQ files including per cycle quality box plots, base proportions, base-level quality trends, relative k-mer diversity, length and occurrence distribution of reads, number of reads above quality cutoffs and mean quality distribution.

fqlist <- seeFastq(fastq=infile1(args), batchsize=10000, klength=8)
pdf("./results/fastqReport.pdf", height=18, width=4*length(fqlist))
seeFastqPlot(fqlist)
dev.off()

fastqReport

Parallelization of QC report on single machine with multiple cores

args <- systemArgs(sysma="tophat.param", mytargets="targets.txt")
f <- function(x) seeFastq(fastq=infile1(args)[x], batchsize=100000, klength=8)
fqlist <- bplapply(seq(along=args), f, BPPARAM = MulticoreParam(workers=8))
seeFastqPlot(unlist(fqlist, recursive=FALSE))

Parallelization of QC report via scheduler (e.g. Torque) across several compute nodes

library(BiocParallel); library(BatchJobs)
f <- function(x) {
    library(systemPipeR)
    args <- systemArgs(sysma="tophat.param", mytargets="targets.txt")
    seeFastq(fastq=infile1(args)[x], batchsize=100000, klength=8)
}
funs <- makeClusterFunctionsTorque("torque.tmpl")
param <- BatchJobsParam(length(args), resources=list(walltime="20:00:00", nodes="1:ppn=1", memory="6gb"), cluster.functions=funs)
register(param)
fqlist <- bplapply(seq(along=args), f)
seeFastqPlot(unlist(fqlist, recursive=FALSE))

3.4 Alignment with Tophat2

Build Bowtie2 index.

args <- systemArgs(sysma="tophat.param", mytargets="targets.txt")
moduleload(modules(args)) # Skip if module system is not available
system("bowtie2-build ./data/tair10.fasta ./data/tair10.fasta")

Execute SYSargs on a single machine without submitting to a queuing system of a compute cluster. This way the input FASTQ files will be processed sequentially. If available, multiple CPU cores can be used for processing each file. The number of CPU cores (here 4) to use for each process is defined in the *.param file. With cores(args) one can return this value from the SYSargs object. Note, if a module system is not installed or used, then the corresponding *.param file needs to be edited accordingly by either providing an empty field in the line(s) starting with module or by deleting these lines.

bampaths <- runCommandline(args=args)

Alternatively, the computation can be greatly accelerated by processing many files in parallel using several compute nodes of a cluster, where a scheduling/queuing system is used for load balancing. To avoid over-subscription of CPU cores on the compute nodes, the value from cores(args) is passed on to the submission command, here nodes in the resources list object. The number of independent parallel cluster processes is defined under the Njobs argument. The following example will run 18 processes in parallel using for each 4 CPU cores. If the resources available on a cluster allow to run all 18 processes at the same time then the shown sample submission will utilize in total 72 CPU cores. Note, runCluster can be used with most queueing systems as it is based on utilities from the BatchJobs package which supports the use of template files *.tmpl) for defining the run parameters of different schedulers. To run the following code, one needs to have both a conf file (see .BatchJob samples here) and a template file (see *.tmpl samples here) for the queueing available on a system. The following example uses the sample conf and template files for the Torque scheduler provided by this package.

file.copy(system.file("extdata", ".BatchJobs.R", package="systemPipeR"), ".")
file.copy(system.file("extdata", "torque.tmpl", package="systemPipeR"), ".")
resources <- list(walltime="20:00:00", nodes=paste0("1:ppn=", cores(args)), memory="10gb")
reg <- clusterRun(args, conffile=".BatchJobs.R", template="torque.tmpl", Njobs=18, runid="01", 
                  resourceList=resources)
waitForJobs(reg)

Useful commands for monitoring progress of submitted jobs

showStatus(reg)
file.exists(outpaths(args))
sapply(1:length(args), function(x) loadResult(reg, x)) # Works after job completion

3.5 Read and alignment count stats

Generate table of read and alignment counts for all samples.

read_statsDF <- alignStats(args) 
write.table(read_statsDF, "results/alignStats.xls", row.names=FALSE, quote=FALSE, sep="\t")

The following shows the first four lines of the sample alignment stats file provided by the systemPipeR package. For simplicity the number of PE reads is multiplied here by 2 to approximate proper alignment frequencies where each read in a pair is counted.

read.table(system.file("extdata", "alignStats.xls", package="systemPipeR"), header=TRUE)[1:4,]
##   FileName Nreads2x Nalign Perc_Aligned Nalign_Primary Perc_Aligned_Primary
## 1      M1A   192918 177961     92.24697         177961             92.24697
## 2      M1B   197484 159378     80.70426         159378             80.70426
## 3      A1A   189870 176055     92.72397         176055             92.72397
## 4      A1B   188854 147768     78.24457         147768             78.24457

Parallelization of read/alignment stats on single machine with multiple cores

f <- function(x) alignStats(args[x])
read_statsList <- bplapply(seq(along=args), f, BPPARAM = MulticoreParam(workers=8))
read_statsDF <- do.call("rbind", read_statsList)

Parallelization of read/alignment stats via scheduler (e.g. Torque) across several compute nodes

library(BiocParallel); library(BatchJobs)
f <- function(x) {
    library(systemPipeR)
    args <- systemArgs(sysma="tophat.param", mytargets="targets.txt")
    alignStats(args[x])
}
funs <- makeClusterFunctionsTorque("torque.tmpl")
param <- BatchJobsParam(length(args), resources=list(walltime="20:00:00", nodes="1:ppn=1", memory="6gb"), cluster.functions=funs)
register(param)
read_statsList <- bplapply(seq(along=args), f)
read_statsDF <- do.call("rbind", read_statsList)

3.7 Alternative NGS Aligners

3.7.1 Alignment with Bowtie2 (e.g. for miRNA profiling)

The following example runs Bowtie2 as a single process without submitting it to a cluster.

args <- systemArgs(sysma="bowtieSE.param", mytargets="targets.txt")
moduleload(modules(args)) # Skip if module system is not available
bampaths <- runCommandline(args=args)

Alternatively, submit the job to compute nodes.

qsubargs <- getQsubargs(queue="batch", cores=cores(args), memory="mem=10gb", time="walltime=20:00:00")
(joblist <- qsubRun(args=args, qsubargs=qsubargs, Nqsubs=18, package="systemPipeR"))

3.7.2 Alignment with BWA-MEM (e.g. for VAR-Seq)

The following example runs BWA-MEM as a single process without submitting it to a cluster.

args <- systemArgs(sysma="bwa.param", mytargets="targets.txt")
moduleload(modules(args)) # Skip if module system is not available
system("bwa index -a bwtsw ./data/tair10.fasta") # Indexes reference genome
bampaths <- runCommandline(args=args)

3.7.3 Alignment with Rsubread (e.g. for RNA-Seq)

The following example shows how one can use within the environment the R-based aligner or other R-based functions that read from input files and write to output files.

library(Rsubread)
args <- systemArgs(sysma="rsubread.param", mytargets="targets.txt")
buildindex(basename=reference(args), reference=reference(args)) # Build indexed reference genome
align(index=reference(args), readfile1=infile1(args), input_format="FASTQ", 
      output_file=outfile1(args), output_format="SAM", nthreads=8, indels=1, TH1=2)
for(i in seq(along=outfile1(args))) asBam(file=outfile1(args)[i], destination=gsub(".sam", "", outfile1(args)[i]), overwrite=TRUE, indexDestination=TRUE)

3.7.4 Alignment with gsnap

Another R-based short read aligner is gsnap from the gmapR package (Wu and Nacu 2010). The code sample below introduces how to run this aligner on multiple nodes of a compute cluster.

library(gmapR); library(BiocParallel); library(BatchJobs)
gmapGenome <- GmapGenome(reference(args), directory="data", name="gmap_tair10chr/", create=TRUE)
args <- systemArgs(sysma="gsnap.param", mytargets="targetsPE.txt")
f <- function(x) {
    library(gmapR); library(systemPipeR)
    args <- systemArgs(sysma="gsnap.param", mytargets="targetsPE.txt")
    gmapGenome <- GmapGenome(reference(args), directory="data", name="gmap_tair10chr/", create=FALSE)
    p <- GsnapParam(genome=gmapGenome, unique_only=TRUE, molecule="DNA", max_mismatches=3)
    o <- gsnap(input_a=infile1(args)[x], input_b=infile2(args)[x], params=p, output=outfile1(args)[x])
}
funs <- makeClusterFunctionsTorque("torque.tmpl")
param <- BatchJobsParam(length(args), resources=list(walltime="20:00:00", nodes="1:ppn=1", memory="6gb"), cluster.functions=funs)
register(param)
d <- bplapply(seq(along=args), f)

4 VAR-Seq workflow for single machine

4.1 Generate workflow template

Load one of the available NGS workflows into your current working directory (here for varseq).

genWorkenvir(workflow="varseq")
setwd("varseq")

4.2 Run workflow

Next, run the chosen sample workflow systemPipeVARseq_single (PDF, Rnw) by executing from the command-line make -B within the varseq directory. Alternatively, one can run the code from the provided *.Rnw template file from within R interactively. Much more detailed information is available in systemPipeR’s overview and workflow vignettes available here.

5 VAR-Seq workflow for computer cluster (demo)

This demonstration will run the above VAR-Seq workflow in parallel on multiple computer nodes of IIGB’s HPC cluster. The workflow template provided for this is called systemPipeVARseq.Rnw (PDF, Rnw) .

6 sessionInfo()

sessionInfo()
## R version 3.2.1 (2015-06-18)
## Platform: x86_64-unknown-linux-gnu (64-bit)
## Running under: Ubuntu 14.04.2 LTS
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C               LC_TIME=en_US.UTF-8       
##  [4] LC_COLLATE=C               LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                  LC_ADDRESS=C              
## [10] LC_TELEPHONE=C             LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
##  [1] ggplot2_1.0.1              systemPipeRdata_0.99.2     systemPipeR_1.3.16        
##  [4] RSQLite_1.0.0              DBI_0.3.1                  ShortRead_1.27.5          
##  [7] GenomicAlignments_1.5.11   SummarizedExperiment_0.3.2 Biobase_2.29.1            
## [10] BiocParallel_1.3.34        Rsamtools_1.21.14          Biostrings_2.37.2         
## [13] XVector_0.9.1              GenomicRanges_1.21.16      GenomeInfoDb_1.5.8        
## [16] IRanges_2.3.14             S4Vectors_0.7.10           BiocGenerics_0.15.3       
## [19] BiocStyle_1.7.4           
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_0.11.6             lattice_0.20-33         GO.db_3.1.2             digest_0.6.8           
##  [5] plyr_1.8.3              futile.options_1.0.0    BatchJobs_1.6           evaluate_0.7           
##  [9] zlibbioc_1.15.0         annotate_1.47.1         Matrix_1.2-2            checkmate_1.6.1        
## [13] rmarkdown_0.7           proto_0.3-10            GOstats_2.35.1          splines_3.2.1          
## [17] stringr_1.0.0           pheatmap_1.0.7          munsell_0.4.2           sendmailR_1.2-1        
## [21] base64enc_0.1-2         BBmisc_1.9              htmltools_0.2.6         fail_1.2               
## [25] edgeR_3.11.2            codetools_0.2-14        XML_3.98-1.3            AnnotationForge_1.11.12
## [29] crayon_1.3.1            MASS_7.3-43             bitops_1.0-6            grid_3.2.1             
## [33] RBGL_1.45.1             xtable_1.7-4            GSEABase_1.31.3         gtable_0.1.2           
## [37] magrittr_1.5            formatR_1.2             scales_0.2.5            graph_1.47.2           
## [41] stringi_0.5-5           hwriter_1.3.2           reshape2_1.4.1          genefilter_1.51.0      
## [45] testthat_0.10.0         limma_3.25.13           latticeExtra_0.6-26     futile.logger_1.4.1    
## [49] brew_1.0-6              rjson_0.2.15            lambda.r_1.1.7          RColorBrewer_1.1-2     
## [53] tools_3.2.1             Category_2.35.1         survival_2.38-3         yaml_2.1.13            
## [57] AnnotationDbi_1.31.17   colorspace_1.2-6        memoise_0.2.1           knitr_1.10.5

References

Girke, Thomas. 2014. “systemPipeR: NGS Workflow and Report Generation Environment.” UC Riverside. https://github.com/tgirke/systemPipeR.

Howard, Brian E, Qiwen Hu, Ahmet Can Babaoglu, Manan Chandra, Monica Borghi, Xiaoping Tan, Luyan He, et al. 2013. “High-Throughput RNA Sequencing of Pseudomonas-Infected Arabidopsis Reveals Hidden Transcriptome Complexity and Novel Splice Variants.” PLoS One 8 (10): e74183. doi:10.1371/journal.pone.0074183.

Kim, Daehwan, Geo Pertea, Cole Trapnell, Harold Pimentel, Ryan Kelley, and Steven L Salzberg. 2013. “TopHat2: Accurate Alignment of Transcriptomes in the Presence of Insertions, Deletions and Gene Fusions.” Genome Biol. 14 (4): R36. doi:10.1186/gb-2013-14-4-r36.

Langmead, Ben, and Steven L Salzberg. 2012. “Fast Gapped-Read Alignment with Bowtie 2.” Nat. Methods 9 (4). Nature Publishing Group: 357–59. doi:10.1038/nmeth.1923.

Lawrence, Michael, Wolfgang Huber, Hervé Pagès, Patrick Aboyoun, Marc Carlson, Robert Gentleman, Martin T Morgan, and Vincent J Carey. 2013. “Software for Computing and Annotating Genomic Ranges.” PLoS Comput. Biol. 9 (8): e1003118. doi:10.1371/journal.pcbi.1003118.

Li, H, and R Durbin. 2009. “Fast and Accurate Short Read Alignment with Burrows-Wheeler Transform.” Bioinformatics 25 (14): 1754–60. doi:10.1093/bioinformatics/btp324.

Li, Heng. 2013. “Aligning Sequence Reads, Clone Sequences and Assembly Contigs with BWA-MEM.” ArXiv [Q-Bio.GN]. http://arxiv.org/abs/1303.3997.

Liao, Yang, Gordon K Smyth, and Wei Shi. 2013. “The Subread Aligner: Fast, Accurate and Scalable Read Mapping by Seed-and-Vote.” Nucleic Acids Res. 41 (10): e108. doi:10.1093/nar/gkt214.

Wu, T D, and S Nacu. 2010. “Fast and SNP-tolerant Detection of Complex Variants and Splicing in Short Reads.” Bioinformatics 26 (7): 873–81. doi:10.1093/bioinformatics/btq057.