--- title: RNA-Seq Workflow author: "First/last name (first.last@ucr.edu)" date: "Last update: `r format(Sys.time(), '%d %B, %Y')`" output: BiocStyle::html_document: toc: true toc_depth: 3 fig_caption: yes 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")), warning=FALSE, message=FALSE) ``` Note: the most recent version of this tutorial can be found here and a short overview slide show [here](https://htmlpreview.github.io/?https://github.com/tgirke/systemPipeR/master/inst/extdata/slides/systemPipeRslides.html). # Introduction Users want to provide here background information about the design of their RNA-Seq project. # Sample definitions and environment settings ## 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 `data` subdirectory and all results will be written to a separate `results` directory, while the `systemPipeRNAseq.Rmd` script and the `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 [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. The following loads one of the available NGS workflow templates (here RNA-Seq) into the user's current working directory. At the moment, the package includes workflow templates for RNA-Seq, ChIP-Seq, VAR-Seq and Ribo-Seq. Templates for additional NGS applications will be provided in the future. ```{r genRna_workflow, eval=FALSE} library(systemPipeRdata) genWorkenvir(workflow="rnaseq", bam=TRUE) setwd("rnaseq") ``` ## Download latest version of this tutorial In case there is an updated version of this tutorial, download its `systemPipeRNAseq.Rmd` source. ```{r download_latest, eval=FALSE} download.file("https://raw.githubusercontent.com/tgirke/CSHL_RNAseq/gh-pages/_vignettes/06_RNAseq/systemPipeRNAseq.Rmd", "systemPipeRNAseq.Rmd") ``` Now open the R markdown script `systemPipeRNAseq.Rmd`in your R IDE (_e.g._ vim-r or RStudio) and run the workflow as outlined below. If you are on Windows please run the following command to move all input files into the right location. This problem will be fixed next week. ```{r windows_fix, eval=FALSE} file.copy(list.files("data/fastq/", "*", full.names=TRUE), "data") file.copy(list.files("data/annotation/", "*", full.names=TRUE), "data", recursive=TRUE) file.copy(list.files("results/bam/", "*", full.names=TRUE), "results") ``` ## Required packages and resources The `systemPipeR` package needs to be loaded to perform the analysis steps shown in this report [@Girke2014-oy]. ```{r load_systempiper, eval=TRUE} library(systemPipeR) ``` In the workflow environments generated by `genWorkenvir` all data inputs are stored in a `data/` directory and all analysis results will be written to a separate `results/` directory, while the `systemPipeRNAseq.Rmd` script and the `targets` file are expected to be located in the parent directory. The R session is expected to run from this parent directory. Additional parameter files are stored under `param/`. To work with real data, users want to organize their own data similarly and substitute all test data for their own data. To rerun an established workflow on new data, the initial `targets` file along with the corresponding FASTQ files are usually the only inputs the user needs to provide. If applicable users can load custom functions not provided by `systemPipeR`. Skip this step if this is not the case. ```{r source_helper_fcts, eval=FALSE} source("systemPipeRNAseq_Fct.R") ``` ## Experiment definition provided by `targets` file The `targets` file defines all FASTQ files and sample comparisons of the analysis workflow. ```{r load_targets, eval=TRUE} targetspath <- system.file("extdata", "targets.txt", package="systemPipeR") targets <- read.delim(targetspath, comment.char = "#")[,1:4] targets ``` # Read preprocessing ## Read quality filtering and trimming 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 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, _e.g._ running the NGS alignments using the trimmed FASTQ files. ```{r fastq_filter, eval=FALSE} args <- systemArgs(sysma="param/trim.param", mytargets="targets.txt") preprocessReads(args=args, Fct="trimLRPatterns(Rpattern='GCCCGGGTAA', subject=fq)", batchsize=100000, overwrite=TRUE, compress=TRUE) writeTargetsout(x=args, file="targets_trim.txt", overwrite=TRUE) ``` ## 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. The results are written to a PDF file named `fastqReport.pdf`. ```{r fastq_report, eval=FALSE} args <- systemArgs(sysma=NULL, mytargets="targets.txt") fqlist <- seeFastq(fastq=infile1(args), batchsize=100000, klength=8) png("./results/fastqReport.png", height=18, width=4*length(fqlist), units="in", res=72) seeFastqPlot(fqlist) dev.off() ``` ![](results/fastqReport.png)
Figure 1: FASTQ quality report for 18 samples

# Alignments ## Read mapping with `Bowtie2/Tophat2` The NGS reads of this project will be aligned against the reference genome sequence using `Bowtie2/TopHat2` [@Kim2013-vg; @Langmead2012-bs]. The parameter settings of the aligner are defined in the `tophat.param` file. ```{r tophat_alignment1, eval=FALSE} args <- systemArgs(sysma="param/tophat.param", mytargets="targets_trim.txt") sysargs(args)[1] # Command-line parameters for first FASTQ file ``` Submission of alignment jobs to compute cluster, here using 72 CPU cores (18 `qsub` processes each with 4 CPU cores). ```{r tophat_alignment2, eval=FALSE, warning=FALSE, message=FALSE} 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) ``` ## Alignment with `Rsubread` The following example shows how one can use within the `systemPipeR` environment the R-based aligner `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="param/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=2) for(i in seq(along=outfile1(args))) asBam(file=outfile1(args)[i], destination=gsub(".sam", "", outfile1(args)[i]), overwrite=TRUE, indexDestination=TRUE) unlink(outfile1(args)); unlink(paste0(outfile1(args),".indel")) ``` ## Read mapping with `HISAT2` ```{r hisat2_alignment, eval=FALSE} args <- systemArgs(sysma="param/hisat2.param", mytargets="targets.txt") # args <- systemArgs(sysma="param/hisat2.param", mytargets="targets_trim.txt") sysargs(args)[1] # Command-line parameters for first FASTQ file moduleload(modules(args)) system("hisat2-build ./data/tair10.fasta ./data/tair10.fasta") runCommandline(args=args) ``` Check whether all BAM files have been created ```{r check_files_exist, eval=FALSE} file.exists(outpaths(args)) ``` ## 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. ```{r align_stats, eval=FALSE} args <- systemArgs(sysma="param/hisat2.param", mytargets="targets.txt") read_statsDF <- alignStats(args=args) write.table(read_statsDF, "results/alignStats.xls", row.names=FALSE, quote=FALSE, sep="\t") ``` ```{r align_stats_view, eval=TRUE} read.table(system.file("extdata", "alignStats.xls", package="systemPipeR"), header=TRUE)[1:4,] ``` ## Create symbolic links for viewing BAM files in IGV The `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 `urlfile` in the `results` directory. ```{r bam_urls, eval=FALSE} symLink2bam(sysargs=args, htmldir=c("~/.html/", "projects/tests/"), urlbase="http://biocluster.ucr.edu/~tgirke/", urlfile="./results/IGVurl.txt") ``` # Read distribution across genomic features The `genFeatures` function generates a variety of feature types from `TxDb` objects using utilities provided by the `GenomicFeatures` package. ## Obtain feature types The first step is the generation of the feature type ranges based on annotations provided by a GFF file that can be transformed into a `TxDb` object. This includes ranges for mRNAs, exons, introns, UTRs, CDSs, miRNAs, rRNAs, tRNAs, promoter and intergenic regions. In addition, any number of custom annotations can be included in this routine. ```{r genFeatures, eval=FALSE} library(GenomicFeatures) txdb <- makeTxDbFromGFF(file="data/tair10.gff", format="gff3", organism="Arabidopsis") feat <- genFeatures(txdb, featuretype="all", reduce_ranges=TRUE, upstream=1000, downstream=0, verbose=TRUE) ``` ## Count reads and plot results The `featuretypeCounts` function counts how many reads in short read alignment files (BAM format) overlap with entire annotation categories. This utility is useful for analyzing the distribution of the read mappings across feature types, _e.g._ coding versus non-coding genes. By default the read counts are reported for the sense and antisense strand of each feature type separately. To minimize memory consumption, the BAM files are processed in a stream using utilities from the `Rsamtools` and `GenomicAlignment` packages. The counts can be reported for each read length separately or as a single value for reads of any length. Subsequently, the counting results can be plotted with the associated `plotfeaturetypeCounts` function. The following generates and plots feature counts for any read length. ```{r featuretypeCounts, eval=FALSE} library(ggplot2); library(grid) fc <- featuretypeCounts(bfl=BamFileList(outpaths(args), yieldSize=50000), grl=feat, singleEnd=TRUE, readlength=NULL, type="data.frame") p <- plotfeaturetypeCounts(x=fc, graphicsfile="results/featureCounts.png", graphicsformat="png", scales="fixed", anyreadlength=TRUE, scale_length_val=NULL) ``` ![](results/featureCounts.png)
Figure 2: Read distribution plot across annotation features for any read length.

# Read quantification per annotation range ## Read counting with `summarizeOverlaps` in parallel mode using multiple cores Reads overlapping with annotation ranges of interest are counted for each sample using the `summarizeOverlaps` function [@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 *reads per kp per million mapped reads* (RPKM). The raw read count table (`countDFeByg.xls`) and the correspoding RPKM table (`rpkmDFeByg.xls`) are written to separate files in the directory of this project. Parallelization is achieved with the `BiocParallel` package, here using 8 CPU cores. ```{r read_counting1, eval=FALSE} library("GenomicFeatures"); library(BiocParallel) txdb <- makeTxDbFromGFF(file="data/tair10.gff", format="gff", dataSource="TAIR", organism="Arabidopsis thaliana") saveDb(txdb, file="./data/tair10.sqlite") txdb <- loadDb("./data/tair10.sqlite") (align <- readGAlignments(outpaths(args)[1])) # Demonstrates how to read bam file into R eByg <- exonsBy(txdb, by=c("gene")) bfl <- BamFileList(outpaths(args), yieldSize=50000, index=character()) multicoreParam <- MulticoreParam(workers=2); 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 ```{r view_counts, eval=FALSE} read.delim("results/countDFeByg.xls", row.names=1, check.names=FALSE)[1:4,1:5] ``` Sample of data slice of RPKM table ```{r view_rpkm, eval=FALSE} 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 `edgeR` or `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, *e.g.* manually comparing the expression levels among different genes or features. ## Sample-wise correlation analysis The following computes the sample-wise Spearman correlation coefficients from the `rlog` transformed expression values generated with the `DESeq2` package. After transformation to a distance matrix, hierarchical clustering is performed with the `hclust` function and the result is plotted as a dendrogram (also see file `sample_tree.pdf`). ```{r sample_tree, eval=FALSE} 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)) png("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() ``` ![](results/sample_tree.png)
Figure 2: Correlation dendrogram of samples

# Analysis of differentially expressed genes with edgeR The analysis of differentially expressed genes (DEGs) is performed with the glm method of the `edgeR` package [@Robinson2010-uk]. The sample comparisons used by this analysis are defined in the header lines of the `targets.txt` file starting with ``. ## Run `edgeR` ```{r run_edger, eval=FALSE} library(edgeR) countDF <- read.delim("results/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 gene descriptions ```{r custom_annot, eval=FALSE} library("biomaRt") m <- useMart("plants_mart", dataset="athaliana_eg_gene", host="plants.ensembl.org") desc <- getBM(attributes=c("tair_locus", "description"), mart=m) 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) ``` ## Plot DEG results Filter and plot DEG results for up and down regulated genes. The definition of *up* and *down* is given in the corresponding help file. To open it, type `?filterDEGs` in the R console. ```{r filter_degs, eval=FALSE} edgeDF <- read.delim("results/edgeRglm_allcomp.xls", row.names=1, check.names=FALSE) png("./results/DEGcounts.png", height=10, width=10, units="in", res=72) DEG_list <- filterDEGs(degDF=edgeDF, filter=c(Fold=2, FDR=20)) dev.off() write.table(DEG_list$Summary, "./results/DEGcounts.xls", quote=FALSE, sep="\t", row.names=FALSE) ``` ![](results/DEGcounts.png)
Figure 3: Up and down regulated DEGs with FDR of 1%

## Venn diagrams of DEG sets The `overLapper` function can compute Venn intersects for large numbers of sample sets (up to 20 or more) and plots 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). ```{r venn_diagram, eval=FALSE} vennsetup <- overLapper(DEG_list$Up[6:9], type="vennsets") vennsetdown <- overLapper(DEG_list$Down[6:9], type="vennsets") pdf("results/vennplot.png") vennPlot(list(vennsetup, vennsetdown), mymain="", mysub="", colmode=2, ccol=c("blue", "red")) dev.off() ``` ![](results/vennplot.png)
Figure 4: Venn Diagram for 4 Up and Down DEG Sets

# GO term enrichment analysis of DEGs ## Obtain gene-to-GO mappings The following shows how to obtain gene-to-GO mappings from `biomaRt` (here for *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 `*.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 `load` function as shown in the next subsection. ```{r get_go_annot, eval=FALSE} library("biomaRt") listMarts() # To choose BioMart database listMarts(host="plants.ensembl.org") m <- useMart("plants_mart", host="plants.ensembl.org") listDatasets(m) m <- useMart("plants_mart", dataset="athaliana_eg_gene", host="plants.ensembl.org") 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") ``` ## 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 `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 function. When `method=all`, it returns all GO terms passing the p-value cutoff specified under the `cutoff` arguments. When `method=slim`, it returns only the GO terms specified under the `myslimv` argument. The given example shows how a GO slim vector for a specific organism can be obtained from BioMart. ```{r go_enrich, eval=FALSE} library("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("plants_mart", dataset="athaliana_eg_gene", host="plants.ensembl.org") 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) ``` ## Plot batch GO term results The `data.frame` generated by `GOCluster` can be plotted with the `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. ```{r go_plot, eval=FALSE} gos <- BatchResultslim[grep("M6-V6_up_down", BatchResultslim$CLID), ] gos <- BatchResultslim png("./results/GOslimbarplotMF.png", height=12, width=12, units="in", res=72) goBarplot(gos, gocat="MF") dev.off() goBarplot(gos, gocat="BP") goBarplot(gos, gocat="CC") ``` ![](results/GOslimbarplotMF.png)
Figure 5: GO Slim Barplot for MF Ontology

# Clustering and heat maps The following example performs hierarchical clustering on the `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. ```{r heatmap, eval=FALSE} library(pheatmap) geneids <- unique(as.character(unlist(DEG_list[[1]]))) y <- assay(rlog(dds))[geneids, ] png("results/heatmap1.png") pheatmap(y, scale="row", clustering_distance_rows="correlation", clustering_distance_cols="correlation") dev.off() ``` ![](results/heatmap1.png)
Figure 6: Heat Map with Hierarchical Clustering Dendrograms of DEGs

# Render report in HTML and PDF format ```{r render_report, eval=FALSE} rmarkdown::render("systemPipeRNAseq.Rmd", c("BiocStyle::html_document", "BiocStyle::pdf_document")) ``` # Version Information ```{r sessionInfo} sessionInfo() ``` # Funding This project was supported by funds from the National Institutes of Health (NIH). # References