--- title: "_systemPipeR_: NGS Workflow and Reporting Environment" subtitle: "Running on Single Machines and Compute Clusters. -- Today's Focus: VAR-seq" author: "Author: Thomas Girke" date: "Date: July 22, 2015" output: BiocStyle::html_document: toc: true toc_depth: 3 vignette: > % \VignetteIndexEntry{Bioconductor for High Throughput Sequence Analysis} % \VignetteEngine{knitr::rmarkdown} fontsize: 14pt bibliography: bibtex.bib --- ```{r style, echo = FALSE, results = 'asis'} BiocStyle::markdown() options(width=100, max.print=1000) knitr::opts_chunk$set( eval=as.logical(Sys.getenv("KNITR_EVAL", "TRUE")), cache=as.logical(Sys.getenv("KNITR_CACHE", "TRUE"))) ``` ```{r setup, echo=FALSE, messages=FALSE, warnings=FALSE} suppressPackageStartupMessages({ library(systemPipeR) library(systemPipeRdata) library(BiocParallel) library(Biostrings) library(Rsamtools) library(GenomicRanges) library(ggplot2) library(GenomicAlignments) library(ShortRead) }) ``` Note: the most recent version of this tutorial can be found here. # Introduction [_`systemPipeR`_](http://www.bioconductor.org/packages/devel/bioc/html/systemPipeR.html) provides utilities for building \textit{end-to-end} analysis workflows with automated report generation for next generation sequence (NGS) applications such as RNA-Seq, ChIP-Seq, VAR-Seq and many others [@Girke2014-oy]. 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` [@Li2013-oy; @Li2009-oc], `TopHat2` [@Kim2013-vg] and `Bowtie2` [@Langmead2012-bs], as well as the R-based NGS aligners [_`Rsubread`_](http://www.bioconductor.org/packages/devel/bioc/html/Rsubread.html) [@Liao2013-bn] and [_`gsnap (gmapR)`_](http://www.bioconductor.org/packages/devel/bioc/html/gmapR.html) [@Wu2010-iq]. 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 [@Lawrence2013-kt]. 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](figures/SystemPipeR_Workflow.png) 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`_](https://github.com/tgirke/systemPipeR/blob/master/inst/extdata/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`_](https://github.com/tgirke/systemPipeR/blob/master/vignettes/systemPipeRNAseq.pdf?raw=true), [_`systemPipeChIPseq`_](https://github.com/tgirke/systemPipeR/blob/master/vignettes/systemPipeChIPseq.pdf?raw=true) and [_`systemPipeVARseq`_](https://github.com/tgirke/systemPipeR/blob/master/vignettes/systemPipeVARseq.pdf?raw=true). To work with _`*.Rnw`_ or _`*.Rmd`_ files efficiently, basic knowledge of [_`Sweave`_](https://www.stat.uni-muenchen.de/~leisch/Sweave/) or [_`knitr`_](http://yihui.name/knitr/) and [_`Latex`_](http://www.latex-project.org/) or [_`R Markdown v2`_](http://rmarkdown.rstudio.com/) is required. # Getting Started ## Installation The R software for running [_`systemPipeR`_](http://www.bioconductor.org/packages/devel/bioc/html/systemPipeR.html) and [_`systemPipeRdata`_](https://github.com/tgirke/systemPipeRdata) can be downloaded from [_CRAN_](http://cran.at.r-project.org/). The _`systemPipeR`_ environment can be installed from R using the `biocLite` install command. ```{r install, eval=FALSE} 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 ``` ## Loading package and documentation ```{r documentation, eval=FALSE} library("systemPipeR") # Loads the package library(help="systemPipeR") # Lists package info vignette("systemPipeR") # Opens vignette ``` ## 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`](http://biocluster.ucr.edu/~tgirke/projects/systemPipeR_test_data.zip). The chosen data set [`SRP010938`](http://www.ncbi.nlm.nih.gov/sra/?term=SRP010938) contains 18 paired-end (PE) read sets from _Arabidposis thaliana_ [@Howard2013-fq]. 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. ## 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. ```{r targetsSE, eval=TRUE} library(systemPipeR) targetspath <- system.file("extdata", "targets.txt", package="systemPipeR") read.delim(targetspath, comment.char = "#") ``` ## Structure of _`targets`_ file for paired end (PE) samples ```{r targetsPE, eval=TRUE} targetspath <- system.file("extdata", "targetsPE.txt", package="systemPipeR") read.delim(targetspath, comment.char = "#")[1:2,1:6] ``` ## Sample comparisons Sample comparisons are defined in the header lines of the _`targets`_ file starting with '``# ``'. 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. ```{r targetscomp, eval=TRUE} readComp(file=targetspath, format="vector", delim="-") ``` ## 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. ```{r param_structure, eval=TRUE} parampath <- system.file("extdata", "tophat.param", package="systemPipeR") read.delim(parampath, comment.char = "#") ``` 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`_. ```{r param_import, eval=TRUE} args <- suppressWarnings(systemArgs(sysma=parampath, mytargets=targetspath)) args ``` Several accessor functions are available that are named after the slot names of the _`SYSargs`_ object class. ```{r sysarg_access, eval=TRUE} names(args) modules(args) cores(args) outpaths(args)[1] sysargs(args)[1] ``` # Workflow overview ## Define environment settings and samples Load package ```{r load_package, eval=FALSE} library(systemPipeR) ``` Construct _`SYSargs`_ object from _`param`_ and _`targets`_ files. ```{r construct_sysargs, eval=FALSE} args <- systemArgs(sysma="trim.param", mytargets="targets.txt") ``` ## 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, \textit{e.g.} running the NGS alignments using the trimmed FASTQ files. ```{r preprocessing, eval=FALSE} 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). ```{r custom_preprocessing, eval=FALSE} 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") ``` ## 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. ```{r fastq_quality, eval=FALSE} fqlist <- seeFastq(fastq=infile1(args), batchsize=10000, klength=8) pdf("./results/fastqReport.pdf", height=18, width=4*length(fqlist)) seeFastqPlot(fqlist) dev.off() ``` ![fastqReport](figures/fastqReport.png) Parallelization of QC report on single machine with multiple cores ```{r fastq_quality_parallel_single, eval=FALSE} 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 ```{r fastq_quality_parallel_cluster, eval=FALSE} 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)) ``` ## Alignment with _`Tophat2`_ Build _`Bowtie2`_ index. ```{r bowtie_index, eval=FALSE} 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. ```{r run_bowtie_single, eval=FALSE} 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](\href{https://code.google.com/p/batchjobs/wiki/DortmundUsage)) and a template file (see _`*.tmpl`_ samples [here](\href{https://github.com/tudo-r/BatchJobs/tree/master/examples)) 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. ```{r run_bowtie_parallel, eval=FALSE} 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 ```{r process_monitoring, eval=FALSE} showStatus(reg) file.exists(outpaths(args)) sapply(1:length(args), function(x) loadResult(reg, x)) # Works after job completion ``` ## Read and alignment count stats Generate table of read and alignment counts for all samples. ```{r align_stats1, eval=FALSE} 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. ```{r align_stats2, eval=TRUE} read.table(system.file("extdata", "alignStats.xls", package="systemPipeR"), header=TRUE)[1:4,] ``` Parallelization of read/alignment stats on single machine with multiple cores ```{r align_stats_parallel, eval=FALSE} 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 ```{r align_stats_parallel_cluster, eval=FALSE} 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) ``` ## Create symbolic links for viewing BAM files in IGV The genome browser IGV supports reading of indexed/sorted BAM files via web URLs. This way it can be avoided to create unnecessary copies of these large files. To enable this approach, an HTML directory with http access needs to be available in the user account (_e.g._ _`home/publichtml`_) of a system. If this is not the case then the BAM files need to be moved or copied to the system where IGV runs. In the following, _`htmldir`_ defines the path to the HTML directory with http access where the symbolic links to the BAM files will be stored. The corresponding URLs will be written to a text file specified under the `_urlfile`_ argument. ```{r igv, eval=FALSE} symLink2bam(sysargs=args, htmldir=c("~/.html/", "somedir/"), urlbase="http://myserver.edu/~username/", urlfile="IGVurl.txt") ``` ## Alternative NGS Aligners ### Alignment with _`Bowtie2`_ (_e.g._ for miRNA profiling) The following example runs _`Bowtie2`_ as a single process without submitting it to a cluster. ```{r bowtie2, eval=FALSE} 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. ```{r bowtie2_cluster, eval=FALSE} 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")) ``` ### 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. ```{r bwamem_cluster, eval=FALSE} 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) ``` ### Alignment with Rsubread (_e.g._ for RNA-Seq) The following example shows how one can use within the \Rpackage{systemPipeR} environment the R-based aligner \Rpackage{Rsubread} or other R-based functions that read from input files and write to output files. ```{r rsubread, eval=FALSE} 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) ``` ### Alignment with _`gsnap`_ Another R-based short read aligner is _`gsnap`_ from the _`gmapR`_ package [@Wu2010-iq]. The code sample below introduces how to run this aligner on multiple nodes of a compute cluster. ```{r gsnap, eval=FALSE} 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) ``` # VAR-Seq workflow for single machine ## Generate workflow template Load one of the available NGS workflows into your current working directory (here for varseq). ```{r genVar_workflow_single, eval=FALSE} genWorkenvir(workflow="varseq") setwd("varseq") ``` ## Run workflow Next, run the chosen sample workflow _`systemPipeVARseq_single`_ ([PDF](https://github.com/tgirke/systemPipeRdata/blob/master/inst/extdata/workflows/varseq/systemPipeVARseq_single.pdf?raw=true), [Rnw](https://github.com/tgirke/systemPipeRdata/blob/master/inst/extdata/workflows/varseq/systemPipeVARseq_single.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](http://www.bioconductor.org/packages/devel/bioc/html/systemPipeR.html). # 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](https://github.com/tgirke/systemPipeRdata/blob/master/inst/extdata/workflows/varseq/systemPipeVARseq.pdf?raw=true), [Rnw](https://github.com/tgirke/systemPipeRdata/blob/master/inst/extdata/workflows/varseq/systemPipeVARseq.Rnw)) . # _`sessionInfo()`_ ```{r sessionInfo} sessionInfo() ``` # References