%\VignetteIndexEntry{RNA-Seq Workflow Template} %\VignetteDepends{rjson, ggplot2, limma, edgeR, GOstats, GO.db, annotate, pheatmap} %\VignetteKeywords{compute cluster, pipeline, reports} %\VignetteEngine{knitr::knitr} %\VignettePackage{systemPipeR} % Generate vignette with knitr % R CMD Sweave --engine=knitr::knitr --pdf systemPipeRNAseq.Rnw \documentclass{article} %<>= %BiocStyle::latex(use.unsrturl=FALSE) %@ <>= BiocStyle::latex(use.unsrturl=FALSE) @ \usepackage[authoryear,round]{natbib} \bibliographystyle{plainnat} \def\bibsection{\section{References}} \usepackage{graphicx} \usepackage{color} \usepackage{hyperref} \usepackage{url} \usepackage{float} %\newcommand{\comment}[1]{} %\newcommand{\Rfunction}[1]{{\texttt{#1}}} %\newcommand{\Robject}[1]{{\texttt{#1}}} %\newcommand{\Rpackage}[1]{{\textit{#1}}} %\newcommand{\Rmethod}[1]{{\texttt{#1}}} \newcommand{\Rfunarg}[1]{{\texttt{#1}}} %\newcommand{\Rclass}[1]{{\textit{#1}}} % Define header and footer area with fandyhdr package (see: http://www.ctan.org/tex-archive/macros/latex/contrib/fancyhdr/fancyhdr.pdf) \usepackage{fancyhdr} \pagestyle{fancy} \fancyhead{} \fancyfoot{} \rhead{\nouppercase{\leftmark}} \lhead{\textit{systemPipeR RNA-Seq Workflow}} \rfoot{\thepage} \begin{document} <>= library(knitr) # set global chunk options for knitr opts_chunk$set(comment=NA, warning=FALSE, message=FALSE, fig.path='figure/systemPipeR-') options(formatR.arrow=TRUE, width=95) unlink("test.db") @ \title{RNA-Seq workflow template: Some Descriptive Title} \author{Project ID: RNAseq\_PI\_Name\_Organism\_Jun2014 \\ Project PI: First Last (first.last@inst.edu)\\ Author of Report: First Last (first.last@inst.edu)} \maketitle \tableofcontents \section{Introduction} This report describes the analysis of an RNA-Seq project from Dr. First Last's lab which studies the gene expression changes of ... in \textit{organism} .... The experimental design is as follows... \section{Sample definitions and environment settings} \subsection{Environment settings and input data} Typically, the user wants to record here the sources and versions of the reference genome sequence along with the corresponding annotations. In the provided sample data set all data inputs are stored in a \Robject{data} subdirectory and all results will be written to a separate \Robject{results} directory, while the \Robject{systemPipeRNAseq.Rnw} script and the \Robject{targets} file are expected to be located in the parent directory. The R session is expected to run from this parent directory. To run this sample report, mini sample FASTQ and reference genome files can be downloaded from \href{http://biocluster.ucr.edu/~tgirke/projects/systemPipeR_test_data.zip}{\textcolor{blue}{here}}. The chosen data set \href{http://www.ncbi.nlm.nih.gov/sra/?term=SRP010938}{\textcolor{blue}{SRP010938}} contains 18 paired-end (PE) read sets from \textit{Arabidposis thaliana} \citep{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 \textit{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. \subsection{Required packages and resources} The \Rpackage{systemPipeR} package needs to be loaded to perform the analysis steps shown in this report \citep{Girke2014-oy}. <>= library(systemPipeR) @ If applicable load custom functions not provided by \Rpackage{systemPipeR} <>= source("systemPipeRNAseq_Fct.R") @ \subsection{Experiment definition provided by \Robject{targets} file} The \href{run:targets.txt}{\Robject{targets}} file defines all FASTQ files and sample comparisons of the analysis workflow. <>= targetspath <- system.file("extdata", "targets.txt", package="systemPipeR") targets <- read.delim(targetspath, comment.char = "#")[,1:4] targets @ \section{Read preprocessing} \subsection{FASTQ quality report} The following \Rfunction{seeFastq} and \Rfunction{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. The results are written to a PDF file named \href{run:./results/fastqReport.pdf}{\Robject{fastqReport.pdf}}. <>= args <- systemArgs(sysma="tophat.param", mytargets="targets.txt") fqlist <- seeFastq(fastq=infile1(args), batchsize=100000, klength=8) pdf("./results/fastqReport.pdf", height=18, width=4*length(fqlist)) seeFastqPlot(fqlist) dev.off() @ \begin{figure}[H] \centering \includegraphics[width=18cm]{fastqReport.pdf} \caption{QC report for 18 FASTQ files.} \label{fig:fastqreport} \end{figure} \section{Alignments} \subsection{Read mapping with \Rfunction{Bowtie2/Tophat2}} The NGS reads of this project will be aligned against the reference genome sequence using \Robject{Bowtie2/TopHat2} \citep{Kim2013-vg, Langmead2012-bs}. The parameter settings of the aligner are defined in the \Robject{tophat.param} file. <>= args <- systemArgs(sysma="tophat.param", mytargets="targets.txt") sysargs(args)[1] # Command-line parameters for first FASTQ file @ Submission of alignment jobs to compute cluster, here using 72 CPU cores (18 \Robject{qsub} processes each with 4 CPU cores). <>= moduleload(modules(args)) system("bowtie2-build ./data/tair10.fasta ./data/tair10.fasta") 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) @ Check whether all BAM files have been created <>= file.exists(outpaths(args)) @ \subsection{Read and alignment stats} The following provides an overview of the number of reads in each sample and how many of them aligned to the reference. <>= read_statsDF <- alignStats(args=args) write.table(read_statsDF, "results/alignStats.xls", row.names=FALSE, quote=FALSE, sep="\t") @ <>= read.table(system.file("extdata", "alignStats.xls", package="systemPipeR"), header=TRUE)[1:4,] @ \subsection{Create symbolic links for viewing BAM files in IGV} The \Rfunction{symLink2bam} function creates symbolic links to view the BAM alignment files in a genome browser such as IGV. The corresponding URLs are written to a file with a path specified under \Robject{urlfile}, here \href{run:./results/IGVurl.txt}{IGVurl.txt}. <>= symLink2bam(sysargs=args, htmldir=c("~/.html/", "somedir/"), urlbase="http://biocluster.ucr.edu/~tgirke/", urlfile="./results/IGVurl.txt") @ \section{Read quantification per annotation range} \subsection{Read counting with \Rfunction{summarizeOverlaps} in parallel mode using multiple cores} Reads overlapping with annotation ranges of interest are counted for each sample using the \Rfunction{summarizeOverlaps} function \citep{Lawrence2013-kt}. The read counting is preformed for exonic gene regions in a non-strand-specific manner while ignoring overlaps among different genes. Subsequently, the expression count values are normalized by \textit{reads per kp per million mapped reads} (RPKM). The raw read count table (\href{run:./results/countDFeByg.xls}{countDFeByg.xls}) and the correspoding RPKM table (\href{run:./results/rpkmDFeByg.xls}{rpkmDFeByg.xls}) are written to separate files in the \Robject{results} directory of this project. Parallelization is achieved with the \Rpackage{BiocParallel} package, here using 8 CPU cores. <>= library("GenomicFeatures"); library(BiocParallel) txdb <- loadDb("./data/tair10.sqlite") eByg <- exonsBy(txdb, by=c("gene")) bfl <- BamFileList(outpaths(args), yieldSize=50000, index=character()) multicoreParam <- MulticoreParam(workers=8); register(multicoreParam); registered() counteByg <- bplapply(bfl, function(x) summarizeOverlaps(eByg, x, mode="Union", ignore.strand=TRUE, inter.feature=FALSE, singleEnd=TRUE)) countDFeByg <- sapply(seq(along=counteByg), function(x) assays(counteByg[[x]])$counts) rownames(countDFeByg) <- names(rowRanges(counteByg[[1]])); colnames(countDFeByg) <- names(bfl) rpkmDFeByg <- apply(countDFeByg, 2, function(x) returnRPKM(counts=x, ranges=eByg)) write.table(countDFeByg, "results/countDFeByg.xls", col.names=NA, quote=FALSE, sep="\t") write.table(rpkmDFeByg, "results/rpkmDFeByg.xls", col.names=NA, quote=FALSE, sep="\t") @ Sample of data slice of count table <>= read.delim("results/countDFeByg.xls", row.names=1, check.names=FALSE)[1:4,1:5] @ Sample of data slice of RPKM table <>= read.delim("results/rpkmDFeByg.xls", row.names=1, check.names=FALSE)[1:4,1:4] @ Note, for most statistical differential expression or abundance analysis methods, such as \Rpackage{edgeR} or \Rpackage{DESeq2}, the raw count values should be used as input. The usage of RPKM values should be restricted to specialty applications required by some users, \textit{e.g.} manually comparing the expression levels among different genes or features. \subsection{Sample-wise correlation analysis} The following computes the sample-wise Spearman correlation coefficients from the \Rfunarg{rlog} transformed expression values generated with the \Rpackage{DESeq2} package. After transformation to a distance matrix, hierarchical clustering is performed with the \Rfunction{hclust} function and the result is plotted as a dendrogram (\href{run:./results/sample_tree.pdf}{sample\_tree.pdf}). <>= library(DESeq2, quietly=TRUE); library(ape, warn.conflicts=FALSE) countDF <- as.matrix(read.table("./results/countDFeByg.xls")) colData <- data.frame(row.names=targetsin(args)$SampleName, condition=targetsin(args)$Factor) dds <- DESeqDataSetFromMatrix(countData = countDF, colData = colData, design = ~ condition) d <- cor(assay(rlog(dds)), method="spearman") hc <- hclust(dist(1-d)) pdf("results/sample_tree.pdf") plot.phylo(as.phylo(hc), type="p", edge.col="blue", edge.width=2, show.node.label=TRUE, no.margin=TRUE) dev.off() @ \begin{figure}[H] \centering \includegraphics[width=10cm]{sample_tree.pdf} \caption{Correlation dendrogram of samples.} \label{fig:sample_tree} \end{figure} \section{Analysis of differentially expressed genes with \Rpackage{edgeR}} The analysis of differentially expressed genes (DEGs) is performed with the glm method from the \Rpackage{edgeR} package \citep{Robinson2010-uk}. The sample comparisons used by this analysis are defined in the header lines of the \href{run:targets.txt}{\Robject{targets}} file starting with \texttt{}. <>= library(edgeR) countDF <- read.delim("countDFeByg.xls", row.names=1, check.names=FALSE) targets <- read.delim("targets.txt", comment="#") cmp <- readComp(file="targets.txt", format="matrix", delim="-") edgeDF <- run_edgeR(countDF=countDF, targets=targets, cmp=cmp[[1]], independent=FALSE, mdsplot="") @ Add custom functional descriptions. Skip this step if \Robject{desc.xls} is not available. <>= desc <- read.delim("data/desc.xls") desc <- desc[!duplicated(desc[,1]),] descv <- as.character(desc[,2]); names(descv) <- as.character(desc[,1]) edgeDF <- data.frame(edgeDF, Desc=descv[rownames(edgeDF)], check.names=FALSE) write.table(edgeDF, "./results/edgeRglm_allcomp.xls", quote=FALSE, sep="\t", col.names = NA) @ Filter and plot DEG results for up and down regulated genes. The definition of '\textit{up}' and '\textit{down}' is given in the corresponding help file. To open it, type \Rfunction{?filterDEGs} in the R console. <>= edgeDF <- read.delim("results/edgeRglm_allcomp.xls", row.names=1, check.names=FALSE) pdf("results/DEGcounts.pdf") DEG_list <- filterDEGs(degDF=edgeDF, filter=c(Fold=2, FDR=1)) dev.off() write.table(DEG_list$Summary, "./results/DEGcounts.xls", quote=FALSE, sep="\t", row.names=FALSE) @ \begin{figure}[H] \centering \includegraphics[width=10cm]{DEGcounts.pdf} \caption{Up and down regulated DEGs with FDR of 1\%.} \label{fig:DEGcounts} \end{figure} The function \Rfunction{overLapper} can compute Venn intersects for large numbers of sample sets (up to 20 or more) and \Rfunction{vennPlot} can plot 2-5 way Venn diagrams. A useful feature is the possiblity to combine the counts from several Venn comparisons with the same number of sample sets in a single Venn diagram (here for 4 up and down DEG sets). <>= vennsetup <- overLapper(DEG_list$Up[6:9], type="vennsets") vennsetdown <- overLapper(DEG_list$Down[6:9], type="vennsets") pdf("results/vennplot.pdf") vennPlot(list(vennsetup, vennsetdown), mymain="", mysub="", colmode=2, ccol=c("blue", "red")) dev.off() @ \begin{figure}[H] \centering \includegraphics[width=14cm]{vennplot.pdf} \caption{Venn Diagram for 4 Up and Down DEG Sets.} \label{fig:vennplot} \end{figure} \subsection{GO term enrichment analysis of DEGs} \subsubsection{Obtain gene-to-GO mappings} The following shows how to obtain gene-to-GO mappings from \Rpackage{biomaRt} (here for \textit{A. thaliana}) and how to organize them for the downstream GO term enrichment analysis. Alternatively, the gene-to-GO mappings can be obtained for many organisms from Bioconductor's \Robject{*.db} genome annotation packages or GO annotation files provided by various genome databases. For each annotation this relatively slow preprocessing step needs to be performed only once. Subsequently, the preprocessed data can be loaded with the \Rfunction{load} function as shown in the next subsection. <>= library("biomaRt") listMarts() # To choose BioMart database m <- useMart("ENSEMBL_MART_PLANT"); listDatasets(m) m <- useMart("ENSEMBL_MART_PLANT", dataset="athaliana_eg_gene") listAttributes(m) # Choose data types you want to download go <- getBM(attributes=c("go_accession", "tair_locus", "go_namespace_1003"), mart=m) go <- go[go[,3]!="",]; go[,3] <- as.character(go[,3]) go[go[,3]=="molecular_function", 3] <- "F"; go[go[,3]=="biological_process", 3] <- "P"; go[go[,3]=="cellular_component", 3] <- "C" go[1:4,] dir.create("./data/GO") write.table(go, "data/GO/GOannotationsBiomart_mod.txt", quote=FALSE, row.names=FALSE, col.names=FALSE, sep="\t") catdb <- makeCATdb(myfile="data/GO/GOannotationsBiomart_mod.txt", lib=NULL, org="", colno=c(1,2,3), idconv=NULL) save(catdb, file="data/GO/catdb.RData") @ \subsubsection{Batch GO term enrichment analysis} Apply the enrichment analysis to the DEG sets obtained the above differential expression analysis. Note, in the following example the \Rfunarg{FDR} filter is set here to an unreasonably high value, simply because of the small size of the toy data set used in this vignette. Batch enrichment analysis of many gene sets is performed with the \Rfunction{GOCluster\_Report} function. When \Rfunarg{method="all"}, it returns all GO terms passing the p-value cutoff specified under the \Rfunarg{cutoff} arguments. When \Rfunarg{method="slim"}, it returns only the GO terms specified under the \Rfunarg{myslimv} argument. The given example shows how a GO slim vector for a specific organism can be obtained from BioMart. <>= load("data/GO/catdb.RData") DEG_list <- filterDEGs(degDF=edgeDF, filter=c(Fold=2, FDR=50), plot=FALSE) up_down <- DEG_list$UporDown; names(up_down) <- paste(names(up_down), "_up_down", sep="") up <- DEG_list$Up; names(up) <- paste(names(up), "_up", sep="") down <- DEG_list$Down; names(down) <- paste(names(down), "_down", sep="") DEGlist <- c(up_down, up, down) DEGlist <- DEGlist[sapply(DEGlist, length) > 0] BatchResult <- GOCluster_Report(catdb=catdb, setlist=DEGlist, method="all", id_type="gene", CLSZ=2, cutoff=0.9, gocats=c("MF", "BP", "CC"), recordSpecGO=NULL) library("biomaRt"); m <- useMart("ENSEMBL_MART_PLANT", dataset="athaliana_eg_gene") goslimvec <- as.character(getBM(attributes=c("goslim_goa_accession"), mart=m)[,1]) BatchResultslim <- GOCluster_Report(catdb=catdb, setlist=DEGlist, method="slim", id_type="gene", myslimv=goslimvec, CLSZ=10, cutoff=0.01, gocats=c("MF", "BP", "CC"), recordSpecGO=NULL) @ \subsubsection{Plot batch GO term results} The \Robject{data.frame} generated by \Rfunction{GOCluster\_Report} can be plotted with the \Rfunction{goBarplot} function. Because of the variable size of the sample sets, it may not always be desirable to show the results from different DEG sets in the same bar plot. Plotting single sample sets is achieved by subsetting the input data frame as shown in the first line of the following example. <>= gos <- BatchResultslim[grep("M6-V6_up_down", BatchResultslim$CLID), ] gos <- BatchResultslim pdf("GOslimbarplotMF.pdf", height=8, width=10); goBarplot(gos, gocat="MF"); dev.off() goBarplot(gos, gocat="BP") goBarplot(gos, gocat="CC") @ \begin{figure}[H] \centering \includegraphics[width=20cm]{GOslimbarplotMF.pdf} \caption{GO Slim Barplot for MF Ontology.} \label{fig:GOMF} \end{figure} \section{Clustering and heat maps} The following example performs hierarchical clustering on the \Rfunarg{rlog} transformed expression matrix subsetted by the DEGs identified in the above differential expression analysis. It uses a Pearson correlation-based distance measure and complete linkage for cluster joining. <>= library(pheatmap) geneids <- unique(as.character(unlist(DEG_list[[1]]))) y <- assay(rlog(dds))[geneids, ] pdf("heatmap1.pdf") pheatmap(y, scale="row", clustering_distance_rows="correlation", clustering_distance_cols="correlation") dev.off() @ \begin{figure}[H] \centering \includegraphics[width=12cm]{heatmap1.pdf} \caption{Heat map with hierarchical clustering dendrograms of DEGs.} \label{fig:heatmap} \end{figure} \section{Version Information} <>= toLatex(sessionInfo()) @ \section{Funding} This project was supported by funds from the National Institutes of Health (NIH). \bibliography{bibtex} \end{document}