%\VignetteIndexEntry{JunctionSeq Vignette} %\VignettePackage{JunctionSeq} %\VignetteEngine{knitr::knitr_notangle} \documentclass[12pt]{article} <>= library("knitr") opts_chunk$set(tidy=FALSE,dev="png",fig.show="hide", fig.width=6.5,fig.height=5.5,fig.keep="high", message=FALSE) @ <>= BiocStyle::latex() @ \usepackage[sort]{cite} \usepackage{framed} \hypersetup{bookmarksopen=false} \usepackage{amsmath,amssymb} %\usepackage[top=0.75in, bottom=0.75in, left=0.75in, right=0.5in]{geometry} %\usepackage{graphicx} %\usepackage[colorlinks=true, linkcolor=blue]{hyperref} %------------------------------------------------------------ % newcommands %------------------------------------------------------------ \newcommand{\myincfigRaw}[3]{% \begin{figure}[h] \centering \includegraphics[width=#2]{figure/#1} \caption{\label{figure/#1}#3} \label{#1-fig} \end{figure} } \newcommand{\myincfig}[3]{% \begin{figure}[h] \centering \includegraphics[width=#2]{figure/#1-1} \caption{\label{figure/#1}#3} \label{#1-fig} \end{figure} } \newcommand{\myincfigTwo}[3]{% \begin{figure}[h] \centering \includegraphics[width=#2]{figure/#1-1} \includegraphics[width=#2]{figure/#1-2} \caption{\label{figure/#1}#3} \label{#1-fig} \end{figure} } \newcommand{\myincfigThree}[3]{% \begin{figure}[h] \centering \includegraphics[width=#2]{figure/#1-1} \includegraphics[width=#2]{figure/#1-2} \includegraphics[width=#2]{figure/#1-3} \caption{\label{figure/#1}#3} \label{#1-fig} \end{figure} } \title{JunctionSeq Package User Manual} \author{Stephen Hartley\\ National Human Genome Research Institute\\ National Institutes of Health } \begin{document} \maketitle %VERSIONS: \begin{center} \textbf{JunctionSeq v}\Sexpr{packageVersion("JunctionSeq")} \end{center} \tableofcontents <>= options(width=80, signif=3, digits=3, prompt=" ", continue=" ") set.seed(0xdada) suppressWarnings({ library("JunctionSeq") library("JctSeqData") }) ## To create bitmap versions of plots with many dots, circumventing ## Sweave's fig=TRUE mechanism... ## (pdfs are too large) openBitmap = function(nm, rows=1, cols=1, height = 600, width = 800, cex = 1.2) { png(paste("QoRT-", nm, ".png", sep=""), width=width*cols, height=height*rows, pointsize=14) par(mfrow=c(rows, cols), cex=cex) } @ %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Overview} \label{sec:praeludium} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- The JunctionSeq R package offers a powerful tool for detecting, identifying, and characterizing differential usage of transcript exons and/or splice junctions in next-generation, high-throughput RNA-Seq experiments. \emph{"Differential usage"} is defined as a differential expression of a particular feature or sub-unit of a gene (such as an exon or spice junction) relative to the overall expression of the gene as a whole (which may or may not itself be differentially expressed). The term was originally used by Anders et.al.\cite{DEXSeq} for their tool, DEXSeq, which is designed to detect differential usage of exonic sub-regions in RNA-Seq data. Tests for differential usage can be used as a proxy for detecting alternative isoform regulation (AIR). Although many tools already exist to detect alternative isoform regulation, the results of these tools are almost universally resistant to interpretation. Alternative isoform regulation is a broad and diverse class of phenomena which may involve alternative splicing, alternative promoting, nucleosome occupancy, long non coding RNAs, alternative polyadenylation, or any number of these factors occurring simultaneously. Unlike simple gene-level differential expression, such regulation cannot be represented by a single p-value and fold-change. The results are often complex and counterintuitive, and effective interpretation is both crucial and non-trivial. JunctionSeq provides a powerful and comprehensive automated visualization toolset to assist in this interpretation process. This includes gene profile plots (see Section \ref{sec:subplot1}), as well as genome-wide browser tracks for use with IGV or the UCSC genome browser (see Section \ref{sec:browserFiles}). Together, these visualizations allow the user to quickly and intuitively interpret and understand the underlying regulatory phenomena that are taking place. This visualization toolset is the core advantage to JunctionSeq over similar methods like DEXSeq. In fact, using a specific set of parameters (outlined in Section \ref{sec:methods-dexreplicate}), you can replicate a DEXSeq analysis in JunctionSeq, and thus use the JunctionSeq visualization tools along with the standard DEXSeq methodology. Under the default parameterization, JunctionSeq also builds upon and expands the basic design put forth by DEXSeq, providing (among other things) the ability to test for both differential exon usage and differential splice junction usage. These two types of analyses are complementary: Exons represent a broad region on the transcripts and thus tend to have higher counts than the individual splice junctions (particularly for large exons) and as such tests for differential exon usage will often have higher power under ideal conditions. However, certain combinations of isoform differentials may not result in strong observed differentials in the individual exons, and up- or down-regulation in unannotated isoforms may not be detectable at all. The addition of splice junction usage analysis makes it possible to detect a broader array of isoform-regulation phenomena, as well as vastly improving the detection of differentials in unannotated isoforms. JunctionSeq is \emph{not} designed to detect changes in overall gene expression. Gene-level differential expression is best detected with tools designed specifically for that purpose such as DESeq2\cite{DESeq2} or edgeR\cite{edgeR}. Additional help and documentation is \href{http://hartleys.github.io/JunctionSeq/index.html}{available online}. There is also a \href{http://hartleys.github.io/JunctionSeq/doc/example-walkthrough.pdf}{comprehensive walkthrough} of the entire analysis pipeline, along with a full \href{ftp://nhgriftp.nhgri.nih.gov/pub/outgoing/mullikin/QoRTsExample/QoRTsPipelineWalkthrough.zip}{example dataset}. Note that the raw \href{ftp://nhgriftp.nhgri.nih.gov/pub/outgoing/mullikin/QoRTsExample/QoRTsPipelineWalkthroughData.zip}{example bam and fastq files} are in a separate download. You can cite JunctionSeq by citing the preprint on arxiv.org: Hartley SW, Mullikin JC. \href{http://arxiv.org/abs/1512.06038}{Detection and Visualization of Differential Exon and Splice Junction Usage in RNA-Seq Data with JunctionSeq.} \href{http://arxiv.org}{arXiv} preprint \href{http://arxiv.org/abs/1512.06038}{arXiv:1512.06038}. 2015 Dec 18. %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Requirements} \label{sec:reqs} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- JunctionSeq requires R v3.2.2 as well as a number of R packages. See the \href{http://hartleys.github.io/JunctionSeq/}{JunctionSeq website} for more detailed and updated information on JunctionSeq installation. JunctionSeq is now part of Bioconductor, and will be included in the next release. If you are on the Devel branch of R (currently v3.3), then you can install JunctionSeq as part of the Devel branch of Bioconductor, using the commands (in R): <>= if (!requireNamespace("BiocManager", quietly=TRUE)) install.packages("BiocManager") BiocManager::install("JunctionSeq") @ If you do not wish to work in the (often unstable) devel Bioconductor branch, you can install JunctionSeq now using a custom install script, available online. Use the commands (in R): <>= source("http://hartleys.github.io/JunctionSeq/install/JS.install.R"); JS.install(); @ Alternatively, JunctionSeq can be installed manually: <>= #CRAN package dependencies: install.packages("statmod"); install.packages("plotrix"); install.packages("stringr"); install.packages("Rcpp"); install.packages("RcppArmadillo"); install.packages("locfit"); install.packages("Hmisc"); #Bioconductor dependencies: if (!requireNamespace("BiocManager", quietly=TRUE)); install.packages("BiocManager"); BiocManager::install(); BiocManager::install("Biobase"); BiocManager::install("BiocGenerics"); BiocManager::install("BiocParallel"); BiocManager::install("GenomicRanges"); BiocManager::install("IRanges"); BiocManager::install("S4Vectors"); BiocManager::install("genefilter"); BiocManager::install("geneplotter"); BiocManager::install("SummarizedExperiment"); BiocManager::install("DESeq2"); install.packages("http://hartleys.github.io/JunctionSeq/install/JunctionSeq_LATEST.tar.gz", repos = NULL, type = "source"); @ \emph{Generating splice junction counts:} The easiest way to generate the counts for JunctionSeq is via the QoRTs\cite{QoRTs} software package. The QoRTs software package requires R version 3.0.2 or higher, as well as (64-bit) java 6 or higher. QoRTs can be found online \href{http://hartleys.github.io/QoRTs/index.html}{here}. \emph{Hardware:} Both the JunctionSeq and QoRTs\cite{QoRTs} software packages will generally require at least 2-10 gigabytes of RAM to run. In general at least 4gb is recommended when available for the QoRTs runs, and at least 10gb for the JunctionSeq runs. R multicore functionality usually involves duplicating the R environment, so memory usage may be higher if multicore options are used. \emph{Annotation:} JunctionSeq requires a transcript annotation in the form of a gtf file. If you are using a annotation guided aligner (which is STRONGLY recommended) it is likely you already have a transcript gtf file for your reference genome. We recommend you use the same annotation gtf for alignment, QC, and downstream analysis. We have found the Ensembl "Gene Sets" gtf\footnote{Which can be acquired from the \href{http://www.ensembl.org}{Ensembl website}} suitable for these purposes. However, any format that adheres to the gtf file specification\footnote{See the gtf file specification \href{http://genome.ucsc.edu/FAQ/FAQformat.html}{here}} will work. \emph{Dataset:} JunctionSeq requires aligned RNA-Seq data. Data can be paired-end or single-end, unstranded or stranded. For paired-end data it is strongly recommended (but not explicitly required) that the SAM/BAM files be sorted either by name or position (it does not matter which). %---------- \subsection{Alignment} \label{sec:align} %---------- QoRTs\cite{QoRTs}, which is used to generate read counts, is designed to run on paired-end or single-end next-gen RNA-Seq data. The data must first be aligned (or "mapped") to a reference genome. RNA-Star\cite{STAR}, GSNAP\cite{GSNAP}, and TopHat2\cite{TopHat2} are all popular and effective aligners for use with RNA-Seq data. The use of short-read or unspliced aligners such as BowTie, ELAND, BWA, or Novoalign is NOT recommended. %----------------------------------------------------------- \subsection{Recommendations} \label{sec:recs} %----------------------------------------------------------- Using barcoding, it is possible to build a combined library of multiple distinct samples which can be run together on the sequencing machine and then demultiplexed afterward. In general, it is recommended that samples for a particular study be multiplexed and merged into "balanced" combined libraries, each containing equal numbers of each biological condition. If necessary, these combined libraries can be run across multiple sequencer lanes or runs to achieve the desired read depth on each sample. This reduces "batch effects", reducing the chances of false discoveries being driven by sequencer artifacts or biases. It is also recommended that the data be thoroughly examined and checked for artifacts, biases, or errors using the QoRTs\cite{QoRTs} quality control package. %---------- \section{Example Dataset} \label{sec:exdata} %---------- To allow users to test JunctionSeq and experiment with its functionality, an example dataset is available online \href{http://hartleys.github.io/JunctionSeq/install/JctSeqData_LATEST.tar.gz}{here}. It can be installed with the command: <>= install.packages("http://hartleys.github.io/JunctionSeq/install/JctSeqData_LATEST.tar.gz", repos=FALSE, type = "source"); @ A more complete walkthrough using this same dataset is available \href{ftp://nhgriftp.nhgri.nih.gov/pub/outgoing/mullikin/QoRTsExample/QoRTsPipelineWalkthrough.zip}{here}. The bam and fastq files are available \href{ftp://nhgriftp.nhgri.nih.gov/pub/outgoing/mullikin/QoRTsExample/QoRTsPipelineWalkthroughData.zip}{here}. The example dataset was taken from rat pineal glands. Sequence data from six samples (aka "biological replicates") are included, three harvested during the day, three at night. To reduce the file sizes to a more managable level, this dataset used only 3 out of the 6 sequencing lanes, and only the reads aligning to chromosome 14 were included. This yielded roughly 750,000 reads per sample. The example dataset, including aligned reads, QC data, example scripts, splice junction counts, and JunctionSeq results, is available online (see the JunctionSeq github page). Information on the original dataset from which the example data was derived is described elsewhere\cite{ratPineal}. THIS EXAMPLE DATASET IS INTENDED FOR TESTING PURPOSES ONLY! The original complete dataset can be downloaded from \href{http://www.ncbi.nlm.nih.gov/geo/}{GEO} (the National Center for Biotechnology Information Gene Expression Omnibus), GEO series accession number GSE63309. The test dataset has been cut down to reduce file sizes and processing time, and a number of artificial "edge cases" were introduced for testing purposes. For example: in the gene annotation, one gene has an artificial transcript that lies on the opposite strand from the other transcripts, to ensure that JunctionSeq deals with that (unlikely) possibility in a sensible way. The results from this test analysis are not appropriate for anything other than testing! Splice junction counts and annotation files generated from this example dataset are included in the \verb|JcnSeqExData| R package, available online (see the JunctionSeq github page), which is what will be used by this vignette. The annotation files can be accessed with the commands: <>= decoder.file <- system.file("extdata/annoFiles/decoder.bySample.txt", package="JctSeqData", mustWork=TRUE); decoder <- read.table(decoder.file, header=TRUE, stringsAsFactors=FALSE); gff.file <- system.file( "extdata/cts/withNovel.forJunctionSeq.gff.gz", package="JctSeqData", mustWork=TRUE); @ The count files can be accessed with the commands: <>= countFiles.noNovel <- system.file(paste0("extdata/cts/", decoder$sample.ID, "/QC.spliceJunctionAndExonCounts.forJunctionSeq.txt.gz"), package="JctSeqData", mustWork=TRUE); countFiles <- system.file(paste0("extdata/cts/", decoder$sample.ID, "/QC.spliceJunctionAndExonCounts.withNovel.forJunctionSeq.txt.gz"), package="JctSeqData", mustWork=TRUE); @ There is also a smaller subset dataset, intended for rapid testing purposes. This is the dataset we will use in this vignette, to reduce build times to within the Bioconductor limitations. The annotation gtf file can be accessed with the command: <>= gff.file <- system.file( "extdata/tiny/withNovel.forJunctionSeq.gff.gz", package="JctSeqData", mustWork=TRUE); @ The count files can be accessed with the commands: <>= countFiles <- system.file(paste0("extdata/tiny/", decoder$sample.ID, "/QC.spliceJunctionAndExonCounts.withNovel.forJunctionSeq.txt.gz"), package="JctSeqData", mustWork=TRUE); @ %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Preparations} \label{sec:prep} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- Once alignment and quality control has been completed on the study dataset, the splice-junction and gene counts must be generated via QoRTs\cite{QoRTs}. To reduce batch effects, the RNA samples used for the example dataset were barcoded and merged together into a single combined library. This combined library was run on three HiSeq 2000 sequencer lanes. Thus, after demultiplexing each sample consisted of three "technical replicates". JunctionSeq is only designed to compare biological replicates, so QoRTs includes functions for generating counts for each technical replicate and then combining the counts across the technical replicates from each biological sample. %---------- \subsection{Generating raw counts via QoRTs} \label{sec:qortsrun} %---------- To generate read counts, you must run QoRTs\cite{QoRTs} on each aligned bam file. QoRTs includes a basic function that calculates a variety of QC metrics along with gene-level and splice-junction-level counts. All these functions can be performed in a single step and a single pass through the input alignment file, greatly simplifying the analysis pipeline. For example, to run QoRTs on the first read-group of replicate \verb|SHAM1_RG1| from the example dataset: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar QC \ --stranded \ inputData/bamFiles/SHAM1_RG1.bam \ inputData/annoFiles/anno.gtf.gz \ rawCts/SHAM1_RG1/ \end{verbatim} \end{framed} Note that the \verb|--stranded| option is required because this example dataset is strand-specific. Also note that QoRTs uses the original gtf annotation file, NOT the flattened gff file produced in section \ref{sec:flatAnno}. More information on this command and on the available options can be found online \href{http://hartleys.github.io/QoRTs/jarHtml/index.html}{here}. For more information about the quality control metrics provided by QoRTs and how to visualize, organize, and view them, see the QoRTs github page and documentation, available online \href{http://hartleys.github.io/QoRTs/index.html}{here}. Certain datasets may require additional options. For example: by default, QoRTs assumes the "MAPQ" rules used by the RNA-STAR aligner and the TopHat 1 aligner. If you use other aligners (TopHat2, GSNAP, etc.) you may need to specify the "--minMAPQ" parameter. See the \href{http://hartleys.github.io/QoRTs/jarHtml/QC.html}{QoRTs QC command documentation}. If Quality Control is being done seperately by other software packages or collaborators, the JunctionSeq counts can be generated without the rest of the QC data by setting the \verb|--runFunctions| option: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar QC \ --stranded \ --runFunctions writeKnownSplices,writeNovelSplices,writeSpliceExon \ inputData/bamFiles/SHAM1_RG1.bam \ inputData/annoFiles/anno.gtf.gz \ rawCts/SHAM1_RG1/ \end{verbatim} \end{framed} This will take much less time to run, as it does not generate the full battery of quality control metrics. The above command will generate a count file \verb|rawCts/SHAM1_RG1/QC.spliceJunctionAndExonCounts.forJunctionSeq.txt.gz|. This file contains both gene-level coverage counts, as well as coverage counts for transcript subunits such as exons and splice junction loci. It will also generate two additional secondary count files (one with all known splice counts, the other with all novel splice counts) which may be needed if you wish to add novel splice junctions to your analysis (see Section \ref{sec:novelCounts}). %---------- \subsection{Merging Counts from Technical Replicates (If Needed)} \label{sec:mergeCounts} %---------- QoRTs\cite{QoRTs} includes functions for merging all count data from various technical replicates. If your dataset does not include technical replicates, or if technical replicates have already been merged prior to the count-generation step then this step is unnecessary. The example dataset has three such technical replicates per sample, which were aligned separately and counted separately. For the purposes of quality control QoRTs was run separately on each of these bam files (making it easier to discern any lane or run specific artifacts that might have occurred). It is then necessary to combine the read counts from each of these bam files. QoRTs includes an automated utility for performing this merge. For the example dataset, the command would be: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar \ mergeAllCounts \ rawCts/ \ annoFiles/decoder.byUID.txt \ cts/ \end{verbatim} \end{framed} The "\verb|rawCts/|" and "\verb|cts/|" directories are the relative paths to the input and output data, respectively. The "decoder" file should be a tab-delimited text file with column titles in the first row. One of the columns must be titled "sample.ID", and one must be labelled "unique.ID". The unique.ID must be unique and refers to the specific technical replicate, the sample.ID column indicates which biological sample each technical replicate belongs to. %---------- \subsection{(Option 1) Including Only Annotated Splice Junction Loci} \label{sec:flatAnno} %---------- If you wish to only test annotated splice junctions, then a simple flat annotation file can be generated for use by JunctionSeq. This file parses the input gtf annotation and assigns unique identifiers to each feature (exon or splice junction) belonging to each gene. These identifiers will match the identifiers listed in the count files produced by QoRTs in the count-generation step (see Section \ref{sec:qortsrun}). \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar makeFlatGff \ --stranded \ annoFiles/anno.gtf.gz \ annoFiles/JunctionSeq.flat.gff.gz \end{verbatim} \end{framed} Note it is \emph{vitally important} that the flat gff file and the read-counts be run using the same "strandedness" options (as described in Section \ref{sec:qortsrun}). If the counts are generated in stranded mode, the gff file must also be generated in stranded mode. %---------- \subsection{(Option 2) Including Novel Splice Junction Loci} \label{sec:novelCounts} %---------- One of the core advantages of JunctionSeq over similar tools such as DEXSeq is the ability to test for differential usage of previously-undiscovered transcripts via novel splice junctions. Most advanced aligners have the ability to align read-pairs to both known and unknown splice junctions. QoRTs\cite{QoRTs} can produce counts for these novel splice junctions, and JunctionSeq can test them for differential usage. Because many splice junctions that appear in the mapped file may only have a tiny number of mapped reads, it is generally desirable to filter out low-coverage splice junctions. Many such splice junctions may simply be errors, and in any case their coverage counts will be too low to detect differential effects. In order to properly filter by read depth, size factors are needed. These can be generated by JunctionSeq (see Section \ref{sec:sf}), or generated from gene-level read counts using DESeq2, edgeR, or QoRTs. See the QoRTs vignette for more information on how to generate these size factors. Once size factors have been generated, novel splice junctions can be selected and counted using the command: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar \ mergeNovelSplices \ --minCount 6 \ --stranded \ cts/ \ annoFiles/decoder.bySample.txt \ annoFiles/anno.gtf.gz \ cts/ \end{verbatim} \end{framed} Note: The decoder file must have a column titled "sample.ID". Optionally, it can also include size factors for each sample (generated by external tools) in a column labelled 'size.factor'. All other columns will be ignored. If no size factors are specified, then size factors will automatically be calculated using the DESeq "geometric" normalization method, based on the gene-level counts for all non-aggregate (ie. non-overlapping) genes. The absence of the overlapping genes may cause these size factors to differ slightly from the DESeq-calculated size-factors. Note: In older versions of QoRts (v1.0.20 and below) the size.factor column was required. If you encounter errors, upgrade to the newest version of QoRTs. This utility finds all splice junctions that fall inside the bounds of any known gene. It then filters this set of splice junctions, selecting only the junction loci with mean normalized read-pair counts of greater than the assigned threshold (set to 6 read-pairs in the example above). It then gives each splice junction that passes this filter a unique identifier. This utility produces two sets of output files. First it writes a .gff file containing the unique identifiers for each annotated and novel-and-passed-filter splice locus. Secondly, for each sample it produces a count file that includes these additional splice junctions (along with the original counts as well). \pagebreak %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Differential Usage Analysis via JunctionSeq} \label{sec:workflow} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %---------- \subsection{Simple analysis pipeline} \label{sec:DJUtest} %---------- JunctionSeq includes a single function that loads the count data and performs a full analysis automatically. This function internally calls a number of sub-functions and returns a JunctionSeqCountSet object that holds the analysis results, dispersions, parameter estimates, and size factor data. % <>= jscs <- runJunctionSeqAnalyses(sample.files = countFiles, sample.names = decoder$sample.ID, condition=factor(decoder$group.ID), flat.gff.file = gff.file, nCores = 1, analysis.type = "junctionsAndExons" ); @ The default analysis type, which is explicitly set in the above command, performs a "hybrid" analysis that tests for differential usage of both exons \emph{and} splice junctions simultaniously. The methods used to test for differential splice junction usage are extremely similar to the methods used for differential exon usage other than the fact that the counts are derived from splice sites rather than exonic regions. These methods are based on the methods used by the DEXSeq package \cite{DEXSeq, DESeq2}. The advantage to our method is that reads (or read-pairs) are never counted more than once in any given statistical test. This is not true in DEXSeq, as the long paired-end reads produced by current sequencers may span several exons and splice junctions and thus be counted several times. Our method also produces estimates that are more intuitive to interpret: our fold change estimates are defined as the fold difference between two conditions across a transcript sub-feature relative to the fold change between the same two conditions for the gene as a whole. The DEXSeq method instead calculates a fold change relative to the sum of the other sub-features, which may change depending on the analysis inclusion parameters. \noindent The above function has a number of optional parameters, which may be relevant to some analyses: \begin{description} \small \item[\texttt{analysis.type}]: By default JunctionSeq simultaniously tests both splice junction loci and exonic regions for differential usage (a "hybrid" analysis). This parameter can be used to limit analyses specifically to either splice junction loci or exonic regions. \item[\texttt{nCores}]: The number of cores to use (note: this will only work when package BiocParallel is used on a Linux-based machine. Unfortunately R cannot run multiple threads on windows at this time). \item[\texttt{meanCountTestableThreshold}]: The minimum normalized-read-count threshold used for filtering out low-coverage features. \item[\texttt{test.formula0}]: The model formula used for the null hypothesis in the ANODEV analysis. \item[\texttt{test.formula1}]: The model formula used for the alternative hypothesis in the ANODEV analysis. \item[\texttt{effect.formula}]: The model formula used for estimating the effect size and parameter estimates. \item[\texttt{geneLevel.formula}]: The model formula used for estimating the gene-level expression. \item[\texttt{gene.names}]: A text file or data.frame that contains two columns: first the geneID then the gene name. Many annotations do not use common gene symbols as the primary identifier for genes. However, full ensembl ID's (for example) are not suited to easy recognition or conversation. Thus, you can use this additional parameter to tell JunctionSeq to use the common gene names in tables and plots (while still using the geneIDs for identification purposes). \end{description} A full description of all these options and how they are used can be accessed using the command: <>= help(runJunctionSeqAnalyses); @ %---------- \subsubsection{Exporting size factors (optional)} \label{sec:sf} %---------- As part of the analysis pipeline, JunctionSeq produces size factors which can be used to "normalize" all samples' read counts to a common scale. These can be accessed using the command: <>= writeSizeFactors(jscs, file = "sizeFactors.txt"); @ These size factors are sometimes needed for various purposes, such as the generation of summary wiggle tracks (see \ref{sec:browserFiles}). %---------- \subsection{Advanced Analysis Pipeline} \label{sec:DJUtestAdv} %---------- Some advanced users may need to deviate from the standard analysis pipeline. They may want to use different size factors, apply multiple different models without having to reload the data from file each time, or use other advanced features offered by JunctionSeq. First you must create a "design" data frame: <>= design <- data.frame(condition = factor(decoder$group.ID)); @ Note: the experimental condition variable MUST be named "condition". Just to demonstrate one advanced feature: We can include the gene names file. This text file contains two columns: first the geneID then the gene name. Many annotations do not use common gene symbols as the primary identifier for genes. However, full ensembl ID's (for example) are not suited to easy recognition or conversation. Thus, you can use this additional parameter to tell JunctionSeq to use the common gene names in tables and plots (while still using the geneIDs for identification purposes). <>= geneID.to.symbol.file <- system.file( "extdata/annoFiles/ensid.2.symbol.txt", package="JctSeqData", mustWork=TRUE); @ Next, the data must be loaded into a JunctionSeqCountSet: <>= jscs = readJunctionSeqCounts(countfiles = countFiles, samplenames = decoder$sample.ID, design = design, flat.gff.file = gff.file, gene.names = geneID.to.symbol.file ); @ Next, size factors must be created and loaded into the dataset: <>= #Generate the size factors and load them into the JunctionSeqCountSet: jscs <- estimateJunctionSeqSizeFactors(jscs); @ Next, we generate test-specific dispersion estimates: <>= jscs <- estimateJunctionSeqDispersions(jscs, nCores = 1); @ Next, we fit these observed dispersions to a regression to create fitted dispersions: <>= jscs <- fitJunctionSeqDispersionFunction(jscs); @ Next, we perform the hypothesis tests for differential splice junction usage: <>= jscs <- testForDiffUsage(jscs, nCores = 1); @ Finally, we calculate effect sizes and parameter estimates: <>= jscs <- estimateEffectSizes( jscs, nCores = 1); @ All these steps simply duplicates the behavior of the runJunctionSeqAnalyses function. However, far more options are available when run in this way. Full documentation of the various options for these commands: <>= help(readJunctionSeqCounts); help(estimateJunctionSeqSizeFactors); help(estimateJunctionSeqDispersions); help(fitJunctionSeqDispersionFunction); help(testForDiffUsage); help(estimateEffectSizes); @ %---------- \subsection{Extracting test results} \label{sec:getres} %---------- Once the differential splice junction usage analysis has been run, the results, including model fits, fold changes, p-values, and coverage estimates can all be written to file using the command: <>= writeCompleteResults(jscs, outfile.prefix="./test", save.jscs = TRUE ); @ This produces a series of output files. The main results files are allGenes.results.txt.gz and sigGenes.results.txt.gz. The former includes rows for all genes, whereas the latter includes only genes that have one or more statistically significant differentially used splice junction locus. The columns for both are: \begin{itemize} \item \emph{featureID}: The unique ID of the exon or splice junction locus. \item \emph{geneID}: The unique ID of the gene. \item \emph{countbinID}: The sub-ID of the exon or splice junction locus. \item \emph{testable}: Whether the locus has sufficient coverage to test. \item \emph{dispBeforeSharing}: The locus-specific dispersion estimate. \item \emph{dispFitted}: The the fitted dispersion. \item \emph{dispersion}: The final dispersion used for the hypothesis tests. \item \emph{pvalue}: The raw p-value of the test for differential usage. \item \emph{padjust}: The adjusted p-value, adjusted using the "FDR" method. \item \emph{chr, start, end, strand}: The (1-based) positon of the exonic region or splice junction. \item \emph{transcripts}: The list of known transcripts that contain this splice junction. Novel splice junctions will be listed as "UNKNOWN\_TX". \item \emph{featureType}: The type of the feature (exonic region, novel splice junction or known splice junction) \item \emph{baseMean}: The base mean normalized coverage counts for the locus across all conditions. \item \emph{HtestCoef(A/B)}: The interaction coefficient from the alternate hypothesis model fit used in the hypothesis tests. This is generally not used for anything but may be useful for certain testing diagnostics. \item \emph{log2FC(A/B)}: The estimated log2 fold change found using the effect-size model fit. This is calculated using a different model than the model used for the hypothesis tests. \item \emph{log2FCvst(A/B)}: The estimated log2 fold change between the vst-transformed estimates. This can be superior to simple fold-change estimates for the purposes of ranking effects by their effect size. Simple fold-change estimates can be deceptively high when expression levels are low. \item \emph{expr\_A}: The estimate for the mean normalized counts of the feature in group A. \item \emph{expr\_B}: The estimate for the mean normalized counts of the feature in group B. \item \emph{geneWisePadj}: The gene-level q-value, for the hypothesis that one or more features belonging to this gene are differentially used. This value will be the same for all features belonging to the same gene. In many studies this may actually be what you want to report. \end{itemize} Note: If the biological condition has more than 2 categories then there will be multiple columns for the HtestCoef and log2FC. Each group will be compared with the reference group. If the supplied condition variable is supplied as a \verb|factor| then the first "level" will be used as the reference group. There will also be a file labelled "sigGenes.genewiseResults.txt.gz", which includes summary information for all statistically-significant genes: \begin{itemize} \item \emph{geneID}: The unique ID of the gene. \item \emph{chr, start, end, strand}: The (1-based) positon of the gene. \item \emph{baseMean}: The base mean normalized coverage counts for the gene across all conditions. \item \emph{geneWisePadj}: The gene-level q-value, for the hypothesis that one or more features belonging to this gene are differentially used. This value will be the same for all features belonging to the same gene. \item \emph{mostSigID}: The sub-feature ID for the most significant exon or splice-junction belonging to the gene. \item \emph{mostSigPadjust}: The adjusted p-value for the most significant exon or splice-junction belonging to the gene. \item \emph{numExons}: The number of (known) non-overlapping exonic regions belonging to the gene. \item \emph{numKnown}: The number of known splice junctions belonging to the gene. \item \emph{numNovel}: The number of novel splice junctions belonging to the gene. \item \emph{exonsSig}: The number of statistically significant non-overlapping exonic regions belonging to the gene. \item \emph{knownSig}: The number of statistically significant known splice junctions belonging to the gene. \item \emph{novelSig}: The number of statistically significant novel splice junctions belonging to the gene. \item \emph{numFeatures}: The columns numExons, numKnown, and numNovel, separated by slashes. \item \emph{numSig}: The columns exonsSig, knownSig, and novelSig, separated by slashes. \end{itemize} The writeCompleteResults function also writes a file: allGenes.expression.data.txt.gz. This file contains the raw counts, normalized counts, expression-level estimates by condition value (as normalized read counts), and relative expression estimates by condition value. These are the same values that are plotted in Section \ref{sec:allPlot}. Finally, this function also produces splice junction track files suitable for use with the UCSC genome browser or the IGV genome browser. These files are described in detail in Section \ref{sec:junctionTracks}. If the save.jscs parameter is set to TRUE, then it will also save a binary representation of the JunctionSeqCountSet object. \pagebreak %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Visualization and Interpretation} \label{sec:vis} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- The interpretation of the results is almost as important as the actual analysis itself. Differential splice junction usage is an observed phenomenon, not an actual defined biological process. It can be caused by a number of underlying regulatory processes including alternative start sites, alternative end sites, transcript truncation, mRNA destabilization, alternative splicing, or mRNA editing. Depending on the quality of the transcript annotation, apparent differential splice junction usage can even be caused by simple gene-level differential expression, if (for example) a known gene and an unannotated gene overlap with one another but are independently regulated. Many of these processes can be difficult to discern and differentiate. Therefore, it is imperative that the results generated by JunctionSeq be made as easy to interpret as possible. JunctionSeq, especially when used in conjunction with the QoRTs\cite{QoRTs} software package, contains a number of tools for visualizing the results and the patterns of expression observed in the dataset. A comprehensive battery of plots, along with a set of html files for easy navigation, can be generated with the command: <>= buildAllPlots(jscs=jscs, outfile.prefix = "./plots/", use.plotting.device = "png", FDR.threshold = 0.01 ); @ This produces a number of sub-directories, and includes analysis-wide summary plots along with plots for every gene with one or more statistically significant exons or splice junctions. Alternatively, plots for a manually-specified gene list can be generated using the \texttt{gene.list} parameter. Other parameters include: { \small \begin{itemize} \item \texttt{colorRed.FDR.threshold}: The adjusted-p-value threshold to use to determine when features should be colored as significant. \item \texttt{plot.gene.level.expression}: If \texttt{FALSE}, do not include the gene-level expression plot on the right side. \item \texttt{sequencing.type}: "paired-end" or "single-end". This option is purely cosmetic, and simply determines whether the y-axes will be labelled as "read-pairs" or "reads". \item \texttt{par.cex, anno.cex.text, anno.cex.axis, anno.cex.main}: The character expansions for various purposes. \texttt{par.cex} determines the cex value passed to \texttt{par}, which expands all other text, margins, etc. The other three parameters determine the size of text for the various labels. \item \texttt{plot.lwd, axes.lwd, anno.lwd, gene.lwd}: Line width variables for various lines. \item \texttt{base.plot.height, base.plot.width, base.plot.units}: The "base" width and height for all plots. Plots with transcripts drawn will be expanded vertically to fit the transcripts, plots with a large number of exons or junctions will be likewise be expanded horizontally. This functionality can be deactivated by setting \texttt{autoscale.height.to.fit.TX.annotation} to FALSE and \texttt{autoscale.width.to.fit.bins} to \texttt{Inf}. \item \texttt{writeHTMLresults}: If \texttt{TRUE}, create a navigable HTML page (default: \texttt{TRUE}). \item \texttt{html.compare.results.list}: If you have multiple studies that you want to directly compare, you can use this parameter to allow cross-linking between analyses from the html navigation pages. The results from each analysis must be placed in the same parent directory, and this parameter must be set to a named list of subdirectory names, corresponding to the \texttt{outfile.prefix} for each run of \texttt{buildAllPlots}. In general you should also set the \texttt{gene.list} parameter, or else cross-links will fail if a gene is significant in one analysis but not the other. This parameter can also be used to cross-link different runs of \texttt{buildAllPlots} on the same dataset, to easily compare different parameter settings. \end{itemize} } \pagebreak %---------- \subsection{Summary Plots} \label{sec:summaryPlots} %---------- JunctionSeq provides functions for two basic summary plots, used to display experiment-wide results. The first is the dispersion plot, which displays the dispersion estimates (y-axis) as a function of the base mean normalized counts (x-axis). The test-specific dispersions are displayed as blue density shading, and the "fitted" dispersion is displayed as a red line. JunctionSeq can also produce an "MA" plot, which displays the fold change("M") on the y-axis as a function of the overall mean normalized counts ("A") on the x-axis. These plots can be generated with the command: <>= plotDispEsts(jscs); plotMA(jscs, FDR.threshold=0.05); @ \myincfigTwo{summaryplots}{0.45\textwidth}{Summary Plots} \pagebreak %---------- \subsection{Gene Profile plots} \label{sec:allPlot} %---------- Generally it is desired to print plots for all genes that have one or more splice junctions with statistically significant differential splice junction usage. By default, JunctionSeq uses an FDR-adjusted p-value threshold of 0.01. For each selected gene, the \texttt{buildAllPlots} function will yield 6 plots: { \small \begin{itemize} \item \emph{geneID-expr.png}: Estimates of average coverage count over each feature, for each value of the biological condition. See Section \ref{sec:subplot1}. \item \emph{geneID-expr-TX.png}: as above, but with the transcript annotation displayed. \item \emph{geneID-normCts.png}: Normalized read counts for each sample, colored by biological condition. See Section \ref{sec:subplot3}. \item \emph{geneID-normCts-TX.png}: as above, but with the transcript annotation displayed. \item \emph{geneID-rExpr.png}: Expression levels relative to the overall gene-level expression. See Section \ref{sec:subplot2}. \item \emph{geneID-rExpr-TX.png}: as above, but with the transcript annotation displayed. \end{itemize} } Note: If the \verb|writeHTMLresults| option is \verb|TRUE|, JunctionSeq will also generate an html index file for easy browsing of these output files. By default, if this option is not deactivated, the 'geneID' will be replaced with a number to reduce the length of the filename (some web servers will have issues if image file names are too long). This behavior can be turned off by setting the \texttt{minimalImageFilenames} option to \texttt{FALSE}. \pagebreak %---------- \subsubsection{Coverage/Expression Profile Plots} \label{sec:subplot1} %---------- %ENSRNOG00000009281 Figure \ref{../plots/expr/2-expr-fig} displays the model estimates of the mean normalized exonic-region/splice-junction coverage counts for each biological condition (in this example: DAY vs NIGHT). Note that these values are not equal to the simple mean normalized read counts across all samples. Rather, these estimates are derived from the GLM parameter estimates (via linear contrasts). \myincfigRaw{../plots/expr/2-expr}{0.85\textwidth}{\footnotesize Feature coverage by biological condition. The plot includes 3 major frames. The top-left frame graphs the expression levels for each tested splice junction (aka the model estimates for the normalized read-pair counts). On the top-right is a box containing the estimates for the overall expression of the gene as a whole (plotted on a distinct y-axis scale). Junctions or exonic regions that are statisically significant are marked with vertical pink lines, and the significant p-values are displayed along the top of the plot. The bottom frame is a line drawing of the known exons for the given gene, as well as all splice junction. Known junctions are drawn with solid lines, novel junctions are dashed, and "untestable" junctions are drawn in grey. Additionally, statistically significant loci are colored pink. Between the two frames a set of lines connect the gene drawing to the expression plots.} There are a few things to note about this plot: \begin{itemize} \item The y-axis is log-transformed, except for the area between 0 and 1 which is plotted on a simple linear scale. \item In the line drawing at the bottom, the exons and introns are not drawn to a common scale. The exons are enlarged to improve readability. Each exon or intron is rescaled proportionately to the square-root of their width, then the exons are scaled up to take up 30 percent of the horizontal space. This can be adjusted via the "exon.rescale.factor" parameter (the default is 0.3), or turned off entirely by setting this parameter to -1. The proportionality function can also be changed using the \texttt{exonRescaleFunction} parameter. For mammalian genomes at least, we have found square-root-proportional function to be a good trade-off between making sure that features are distinct while still being able to visually identify which features are larger or smaller. \item Note that the rightmost junction (J005) is marked as significantly differentially used, even though it is NOT differentially expressed. This is because JunctionSeq does not test for simple differential expression. It tests for differential splice junction coverage relative to gene-wide expression. Therefore, if (as in this case) the gene as a whole is strongly differentially expressed, then a splice junction that is NOT differentially expressed is the one that is being differentially used. \item Similarly the second junction from the right (N008), which is novel, is differentially used in the opposite way: while the gene itself is somewhat differentially expressed (at a fold change of roughly 3-4x overall), this particular junction has a \emph{massive} differential, far beyond that found in the gene as a whole (at 45x fold change). Thus, it is being differentially used. \end{itemize} \pagebreak Figure \ref{../plots/exprTX/2-expr-TX-fig} displays the same information found in Figure \ref{../plots/expr/2-expr-fig}, except with all known transcripts plotted beneath the main plot. \myincfigRaw{../plots/exprTX/2-expr-TX}{0.65\textwidth}{\footnotesize Feature coverage by biological condition, with annotated transcripts displayed. This plot is identical to the previous, except all annotated transcripts are displayed below the standard plot.} \pagebreak %---------- \subsubsection{Normalized Count Plots} \label{sec:subplot3} %---------- \myincfigRaw{../plots/normCounts/2-normCts}{0.85\textwidth}{\footnotesize Normalized counts for each sample} Figure \ref{../plots/normCounts/2-normCts-fig} displays the normalized coverage counts for each sample over each splice junction or exonic region. This will be identical to the counts displayed in Section \ref{sec:subplot4} except that each read count will be normalized by the sample size factors so that the samples can be compared directly. \pagebreak %---------- \subsubsection{Relative Expression Plots} \label{sec:subplot2} %---------- \myincfigRaw{../plots/rExpr/2-rExpr}{0.85\textwidth}{\footnotesize Relative splice junction coverage by biological condition.} Figure \ref{../plots/rExpr/2-rExpr-fig} displays the "relative" coverage for each splice junction or exonic region, \emph{relative} to gene-wide expression (which is itself plotted in the right panel). JunctionSeq is designed to detect differential expression even in the presence of gene-wide, multi-transcript differential expression. However it can be difficult to visually assess differential usage on differentially exprssed genes. This plot displays the coverage relative to the gene-wide expression. These estimates are derived from the GLM parameter estimates (via linear contrasts). The reported "fold changes" are ratios of these relative expression estimates. \pagebreak %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \subsection{Generating Genome Browser Tracks} \label{sec:browserFiles} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- Once potential genes of interest have been identified via JunctionSeq, it can be helpful to examine these genes manually on a genome browser (such as the UCSC genome browser or the IGV browser). This can assist in the identification of potential sources of artifacts or errors (such as repetitive regions or the presence of unannotated overlapping features) that may underlie false discoveries. To this end, the QoRTs\cite{QoRTs} and JunctionSeq software packages include a number of tools designed to assist in generating simple and powerful browser tracks designed to aid in the interpretation of the data and results. \myincfigRaw{../figs/final}{0.95\textwidth}{ \footnotesize An example of the browser tracks that can be generated using QoRTs and JunctionSeq. The top two "MultiWig" tracks display the mean normalized coverage for 100-base-pair windows across the genome, by biological group (DAY or NIGHT). The top track displays the coverage across the forward (genomic) strand, and the second track displays the coverage across the reverse strand. In both tracks the mean normalized coverage depth across the three NIGHT samples is displayed in blue and the mean normalized coverage depth across the three DAY samples is displayed in red. The overlap is colored black. The third track displays the mean normalized coverage across all testable splice junction loci for each biological condition. Each junction is labelled with the condition ID (DAY or NIGHT), the splice junction ID (J for annotated, N for novel), followed by the mean normalized coverage across that junction and biological condition, in parentheses. Once again, DAY samples are displayed in red and NIGHT samples are displayed in blue. The final bottom track displays the splice junctions that exhibit statistically significant differential usage. Each junction is labelled with the splice junction ID and the p-value. These images were produced by the \href{https://genome.ucsc.edu/}{UCSC genome browser}, and a browser session containing these tracks in this configuration is available online \href{https://genome.ucsc.edu/cgi-bin/hgTracks?hgS_doOtherUser=submit&hgS_otherUserName=stephen.hartley&hgS_otherUserSessionName=JunctionSeqExample_Full}{here}} Advanced tracks like those displayed in Figure \ref{../figs/final-fig} can be used to visualize the data and can aid in determining the form of regulatory activity that underlies any apparent differential splice junction usage. The wiggle files needed to produce such tracks can be generated via QoRTs and JunctionSeq. Configuring the multi-colored "MultiWig" tracks in the UCSC browser require the use of \href{https://genome.ucsc.edu/goldenPath/help/hgTrackHubHelp.html}{track hubs}, the configuration of which is beyond the scope of this manual. More information on track hubs can be found \href{https://genome.ucsc.edu/goldenPath/help/hgTrackHubHelp.html}{on the UCSC browser documentation}. %---------- \subsubsection{Wiggle Tracks} \label{sec:wiggleTracks} %---------- Both IGV and the UCSC browser can display "wiggle" tracks, which can be used to display coverage depth across the genome. QoRTs\cite{QoRTs} includes functions for generating these wiggle files for each sample/replicate, merging across technical replicates, and computing mean normalized coverages across multiple samples for each biological condition. \myincfigRaw{../figs/SHAM1_RG1-wiggle}{0.9\textwidth}{\footnotesize Two "wiggle" tracks displaying the forward- and reverse-strand coverage for replicate SHAM1\_RG1. These tracks display the read-pair mean coverage depth for each 100-base-pair window across the whole genome. The reverse strand is displayed as negative values. These tracks have been loaded into the \href{https://genome.ucsc.edu/}{UCSC genome browser}, and a browser session containing these tracks in this configuration is available online \href{https://genome.ucsc.edu/cgi-bin/hgTracks?hgS_doOtherUser=submit&hgS_otherUserName=stephen.hartley&hgS_otherUserSessionName=JunctionSeqExample_SHAM1_RG1wiggle}{here}} Figure \ref{../figs/SHAM1_RG1-wiggle-fig} shows an example pair of "wiggle" files produced by QoRTs for replicate SHAM1\_RG1. QoRTs includes two ways to generate such wiggle files from a sam/bam file: The first way is to create these files at the same time as the read counts. To do this, simply add the "--chromSizes" parameter like so: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar \ --stranded \ --chromSizes inputData/annoFiles/chrom.sizes \ inputData/bamFiles/SHAM1_RG1.chr14.bam \ inputData/annoFiles/rn4.anno.chr14.gtf.gz \ outputData/qortsData/SHAM1_RG1/ \end{verbatim} \end{framed} By default this will cause QoRTs to generate the wiggle file(s) for this sample. Note that if the --runFunctions parameter is being included, you must also include in the function list the function "makeWiggles". Note that if the data is stranded (as in the example dataset), then two wiggle files will be generated, one for each strand. Alternatively, the wiggle file(s) can be generated manually using the command: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar \ bamToWiggle \ --stranded \ --negativeReverseStrand \ --includeTrackDefLine \ inputData/bamFiles/SHAM1_RG1.chr14.bam \ SHAM1_RG1 \ inputData/annoFiles/rn4.chr14.chrom.sizes \ outputData/qortsData/SHAM1_RG1/QC.wiggle \end{verbatim} \end{framed} If this step is performed prior to merging technical replicates (see Section \ref{sec:mergeCounts}), and if the standard file-name conventions are followed (as displayed in the examples above), then the technical-replicate wiggle files will automatically be merged along with the other count information in by the QoRTs technical replicate merge utility. \pagebreak %---------- \subsubsection{Merging Wiggle Tracks} \label{sec:wiggleMerge} %---------- \myincfigRaw{../figs/DAY-wiggle}{0.9\textwidth}{\footnotesize Two "wiggle" tracks displaying the forward- and reverse-strand coverage for the three "DAY" samples. These tracks display the mean normalized read-pair coverage depth for each 100-base-pair window across the whole genome. The reverse strand is displayed as negative values. These tracks have been loaded into the \href{https://genome.ucsc.edu/}{UCSC genome browser}, and a browser session containing these tracks in this configuration is available online \href{https://genome.ucsc.edu/cgi-bin/hgTracks?hgS_doOtherUser=submit&hgS_otherUserName=stephen.hartley&hgS_otherUserSessionName=JunctionSeqExample_DAYwiggle}{here}.} QoRTs\cite{QoRTs} can generate wiggle files containing mean normalized coverage counts across a group of samples, as shown in Figure \ref{../figs/DAY-wiggle-fig}. These can be generated using the command: \begin{framed} \small \begin{verbatim} java -jar /path/to/jarfile/QoRTs.jar \ mergeWig \ --calcMean \ --trackTitle DAY_FWD \ --infilePrefix outputData/countTables/ \ --infileSuffix /QC.wiggle.fwd.wig.gz \ --sizeFactorFile sizeFactors.GEO.txt \ --sampleList SHAM1,SHAM2,SHAM3 \ outputData/DAY.fwd.wig.gz \end{verbatim} \end{framed} The "--sampleList" parameter can also accept data from standard input ("-"), or a text file (which must end ".txt") containing a list of sample ID's, one on each line. \pagebreak %---------- \subsubsection{Splice Junction Tracks} \label{sec:junctionTracks} %---------- \myincfigRaw{../figs/junctions}{0.9\textwidth}{\footnotesize The first track displays the mean-normalized coverage for each splice junction for DAY (blue) and NIGHT (red) conditions. Each splice junction is marked with the condition ID (DAY or NIGHT), followed by the splice junction ID (J for annotated, N for novel), followed by the mean normalized read count in parentheses. The second track displays the splice junctions that display statistically-significant differential usage, along with the junction ID and the adjusted p-value in parentheses (rounded to 4 digits). This track has been loaded into the \href{https://genome.ucsc.edu/}{UCSC genome browser}, and a browser session containing these tracks in this configuration is available online \href{https://genome.ucsc.edu/cgi-bin/hgTracks?hgS_doOtherUser=submit&hgS_otherUserName=stephen.hartley&hgS_otherUserSessionName=JunctionSeqExample_BEDs}{here}.} Splice Junction Tracks are generated automatically by the \verb|writeCompleteResults| function, assuming the \verb|write.bedTracks| parameter has not been set to FALSE. By default, five bed tracks are generated: \begin{itemize} \small \item \emph{allGenes.junctionCoverage.bed.gz}: A track that lists the mean normalized read (or read-pair) coverage for each splice junction and for each biological condition. These are not equal to the actual mean normalized counts, rather these are derived from the parameter estimates. \item \emph{sigGenes.junctionCoverage.bed.gz}: Identical to the previous, but only for genes with statistically significant features. \item \emph{allGenes.exonCoverage.bed.gz}: Similar to the Junction Coverage plots, but displaying coverage across exonic regions. \item \emph{sigGenes.exonCoverage.bed.gz}: Identical to the previous, but only for genes with statistically significant features. \item \emph{sigGenes.pvalues.bed.gz}: This track displays all statistically-significant loci (both exonic regions and splice junctions), with the adjusted p-value listed in parentheses. \end{itemize} These five files can be generated via the \verb|writeCompleteResults| command (see Section \ref{sec:getres}) Individual-sample or individual-replicate junction coverage tracks can also be generated via QoRTs. See the QoRTs\cite{QoRTs} documentation available \href{http://hartleys.github.io/QoRTs/index.html}{online}. \pagebreak %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \subsection{Additional Plotting Options} \label{sec:nonDefaultPlots} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %---------- \subsubsection{Raw Count Plots} \label{sec:subplot4} %---------- Raw expression plots can be added using the command: <>= buildAllPlots(jscs=jscs, outfile.prefix = "./plots2/", use.plotting.device = "png", FDR.threshold = 0.01, expr.plot = FALSE, normCounts.plot = FALSE, rExpr.plot = FALSE, rawCounts.plot = TRUE ); @ \myincfigRaw{../plots2/rawCounts/2-rawCts}{0.5\textwidth}{\footnotesize Raw counts for each sample} Figure \ref{../plots2/rawCounts/2-rawCts-fig} displays the raw (un-normalized) coverage counts for each sample over each splice junction or exonic region. This is equal to the number of reads (or read-pairs, for paired-end data) that bridge each junction. Note that these counts are not normalized, and are generally not directly comparable. Note that by default this plot is NOT generated by \verb|buildAllPlots|. Generation of these plots must be turned on by adding the parameter: "\verb|rawCounts.plot=TRUE|". \pagebreak %---------- \subsubsection{Additional Optional Parameters} \label{sec:subplot6} %---------- By setting the \texttt{plot.exon.results} and \texttt{plot.junction.results} parameters to \texttt{FALSE}, we can exclude exons or junctions from the plots, respectively. There are numerous other options that can be used to generate plots either individually or in groups. \noindent Just a few examples: <>= #Make a battery of exon-only plots for one gene only: buildAllPlotsForGene(jscs=jscs, geneID = "ENSRNOG00000009281", outfile.prefix = "./plots/", use.plotting.device = "png", colorRed.FDR.threshold = 0.01, #Limit plotting to exons only: plot.junction.results = FALSE, #Change the fill color of significant exons: colorList = list(SIG.FEATURE.FILL.COLOR = "red"), ); #Make a set of Junction-Only Plots for a specific list of interesting genes: buildAllPlots(jscs=jscs, gene.list = c("ENSRNOG00000009281"), outfile.prefix = "./plotsJCT/", use.plotting.device = "png", FDR.threshold = 0.01, #Do not graph exonic regions: plot.exon.results = FALSE, ); @ <>= #NOTE: This is a version of the above script with additional options to # reduce execution time: buildAllPlotsForGene(jscs=jscs, geneID = "ENSRNOG00000009281", outfile.prefix = "./plots/", use.plotting.device = "png", colorRed.FDR.threshold = 0.01, #Limit plotting to exons only: plot.junction.results = FALSE, #Change the fill color of significant exons: colorList = list(SIG.FEATURE.FILL.COLOR = "red"), #NOTE: The following options reduce the generation of # extra plots not actually needed for this manual: expr.plot = TRUE, normCounts.plot = FALSE, rExpr.plot = FALSE, rawCounts.plot = FALSE, without.TX = FALSE ); #Junction-Only Plots: buildAllPlots(jscs=jscs, gene.list = c("ENSRNOG00000009281"), outfile.prefix = "./plotsJCT/", use.plotting.device = "png", FDR.threshold = 0.01, plot.exon.results = FALSE, #NOTE: The following options reduce the generation of # extra plots not actually needed for this manual: expr.plot = TRUE, normCounts.plot = FALSE, rExpr.plot = FALSE, rawCounts.plot = FALSE, ma.plot = FALSE, variance.plot = FALSE, without.TX = FALSE, writeHTMLresults = FALSE ); @ See Figure \ref{../plots/ENSRNOG00000009281-expr-TX-fig} and Figure \ref{../plotsJCT/exprTX/ENSRNOG00000009281-expr-TX-fig}. For more information, use the help documentation: <>= help(buildAllPlots); help(buildAllPlotsForGene); help(plotJunctionSeqResultsForGene); @ \pagebreak \myincfigRaw{../plots/ENSRNOG00000009281-expr-TX}{0.85\textwidth}{\footnotesize Exon expression plots. Note the altered fill color of the significant exon.} \myincfigRaw{../plotsJCT/exprTX/ENSRNOG00000009281-expr-TX}{0.85\textwidth}{\footnotesize Splice Junction expression plots.} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %\section{Interpretation of the results} \label{sec:interpretation} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %MORE INFO HERE. \pagebreak %------------------------------------------------------------------- %------------------------------------------------------------------- \section{Statistical Methodology} \label{sec:methods-main} %------------------------------------------------------------------- %------------------------------------------------------------------- The statistical methods are based on those described in (S. Anders et al., 2012), as implemented in the DEXSeq Bioconductor package. However these methods have been expanded, adapted, and altered in a number of ways in order to accurately and efficiently test for differential junction usage. %------------------------------------------------------------------- \subsection{Preliminary Definitions} \label{sec:methods-prelimdefs} %------------------------------------------------------------------- The methods used to generate exonic segment counting bins is idential to the methods used in DEXSeq. Briefly: for each gene the set of all transcripts are combined and exons are broken up into mutually-exclusive sub-segments. Using the transcript annotation we generate a list of "features" for each gene. These features include all the mutually-exclusive exonic segments, along with all the splice junctions (both annotated splice junctions and novel splice junctions that pass the coverage thresholds). For each sample \(i\) and feature \(j\) on gene \(g\) we define the counts: \begin{subequations} \begin{align} k_{ji}^{\text{Feature}} &= \text{\# reads/pairs covering feature }j\text{ in sample }i \label{eq:kJunct} \\ &\text{and} \nonumber \\ k_{gi}^{\text{Gene}} &= \text{\# reads/pairs covering gene }g\text{ in sample }i \label{eq:kGene} \end{align} \end{subequations} Reads are counted towards an exonic segment feature if they overlap with any portion of that exonic segment. Reads are counted towards a splice junction feature if they map across the splice junction. Paired-end reads in which the two mates lie on opposite sides of a junction are NOT counted, as the specific splice junction between them cannot be reliably identified. The count \(k_{gi}^{\text{Gene}}\) is calculated using the same methods already in general use for gene-level differential expression analysis (using the "union" rule). Briefly: any reads or read-pairs that cover any part of any of the exons of any one unique gene are counted towards that gene. Purely-intronic reads are ignored. It should be noted that the QoRTs script that produces the exonic segments differs slightly in its output from the DEXSeq script \verb|prepare_annotation.py|. This is because the DEXSeq script uses unordered data structures and thus the order in which elements appear is not defined. Additionally, for unstranded data QoRTs version operates in a slightly-different "unstranded" mode, in which genes that overlap across opposite strands are counted as overlapping. %------------------------------------------------------------------- \subsection{Model Framework} \label{sec:methods-modelFramework} %------------------------------------------------------------------- Each feature is fitted to a separate statistical model. For a given feature j located on gene g, we define two counting bins: covering the feature \(j\), and covering the gene \(g\), but NOT covering the feature \(j\). Thus we define \(y_{1i}\) and \(y_{0i}\) for each sample \(i \in \{1,2,\dots,n\}\): \begin{subequations} \begin{align} y_{1i} &= k_{ji}^{\text{Feature}} \\ &\text{and} \nonumber \\ y_{0i} &= k_{gi}^{\text{Gene}} - k_{ji}^{\text{Feature}} \end{align} \end{subequations} Thus, \(y_{1i}\) is simply equal to the number of reads covering the feature in sample \(i\), and \(y_{0i}\) is equal to the number of reads covering the gene but NOT covering the feature. Note that while JunctionSeq generally uses methods similar to those used by DEXSeq, this framework differs from that used by DEXSeq on exon counting bins. In the framework used by DEXSeq, the feature count (i.e. \(k_{ji}^{\text{Feature}}\)) is compared with the sum of all other feature counts for the given gene. This means that some reads may be counted more than once if they span multiple features. When reads are relatively short (as was typical when DEXSeq was designed) this effect is minimal, but it becomes more problematic as reads become longer. This problem is also exacerbated in genes with a large number of features in close succession. Under our framework, no read-pair is ever counted more than once in a given model. As in DEXSeq, we assume that the count \(y_{bi}\) is a realization of a negative-binomial random variable \(Y_{bi}\): \begin{equation} Y_{bi} \sim NegBin(\text{mean} = s_{i}\mu_{bi}, \text{dispersion} = \alpha_{j}) \end{equation} Where \(\alpha_{j}\) is the dispersion parameter for the current splice junction \(j\), \(s_{i}\) is the normalization size factor for each sample \(i\), and \(\mu_{bi}\) is the mean for sample \(i\) and counting-bin \(b\). Size factors \(s_{i}\) are estimated using the "geometric" normalization method, which is the default method used by DESeq, DESeq2, DEXSeq, and CuffDiff. Unlike with DEXSeq, these size factors are calculated using the gene-level counts \(k_{gi}^{\text{Gene}}\), and thus does not implicitly assign excess weight to genes with a large number of features. %------------------------------------------------------------------- \subsection{Dispersion Estimation} \label{sec:methods-disp} %------------------------------------------------------------------- In many high-throughput sequencing experiments there are too few replicates to directly estimate the locus-specific dispersion term \(\alpha_j\) for each splice junction \(j\). This problem is well-characterized, and a number of different solutions have been proposed, the vast majority of which involve sharing information between loci across the genome. By default JunctionSeq uses the exact same method for dispersion estimation used by the more recent versions of the DEXSeq package (v1.14.0 and above) and by the DESeq2 package. This method is described at length in the DESeq2 methods paper \cite{DESeq2}. Briefly: "feature-specific" estimates of dispersion are generated via a Cox-Reid adjusted-log-likelihood-based method, and these dispersions are fitted to a parametric trend. These "fitted" and "feature-specific" estimates of dispersion are combined, and we use the maximum \emph{a posteriori} estimate as the final dispersion for each feature. It is currently unclear whether splice junction loci and exonic segment loci can be reasonably expected to follow the same trend for dispersion. By default JunctionSeq fits the two types of features separately. For most datasets the difference appears to be negligible. JunctionSeq offers a number of alternatives to this methodology, and can calculate dispersions using the methods of older versios of DEXSeq (v1.8.0 and earlier). See Section \ref{sec:methods-alts} for more information. %------------------------------------------------------------------- \subsection{Hypothesis Testing} \label{sec:methods-htest} %------------------------------------------------------------------- Other than our use of overall-gene counts and our use of splice junction counts, the methods used in the hypothesis test are identical to those used by the more recent versions of DEXSeq (v1.12.0 and higher). "Differential usage" is not in and of itself a biological process. Rather, it is a \textbf{symptom} of transcript-specific differential regulation, which can take many forms. Put simply: we are attempting to test whether the fold-change for the biological condition across feature \(j\) is the same as the fold-change for the biological condition across the gene \(g\) as a whole. This can occur both with and without overall gene-level differential expression. Some features marked as "differentially used" might show a strong differential and be part of a gene that is not, as a whole, differentially expressed. Alternatively, a feature might be marked as "differentially used" if it displays flat expression but is part of a gene that otherwise shows extreme differentials. In both cases, the feature does not show the same pattern of expression relative to the gene as a whole. Thus, some form of transcript-specific differential regulation must be taking place. In statistical terms: we are attempting to detect "interaction" between the counting-bin variable B and the experimental-condition variable C. Thus, two models are fitted to the mean \(\mu_{bi}\): \begin{subequations} \label{eq:hypoModels} \begin{align} H_{0}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{S}_{i} \label{eq:hypoModels0} \\ H_{1}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{S}_{i} + \beta^{CB}_{\rho_{i}b} \label{eq:hypoModels1} \end{align} \end{subequations} Where \(\rho_{j}\) is the biological condition (eg case/control status) of sample \(j\). Note that the bin-condition interaction term (\(\beta^{CB}_{\rho_{i}b}\)) is included, but the condition main-effect term (\(\beta^{C}_{\rho{i}}\)) is absent. This term can be omitted is because JunctionSeq is not designed to detect or assess gene-level differential expression. Thus there are two components that can be treated as "noise": variation in junction-level expression and variation in gene-level expression. As proposed by Anders et. al., we use a main-effects term for the sample ID (\(\beta^{S}_{i}\)), which subsumes the condition main-effect term. This subsumes both the differential trends and the random variations (noise) in the gene-level expression, improving the power towards detecting differential interaction between the count-bin term and the experimental-condition term. These models can easily be extended to include confounding variables: For confounding variable \(\tau\), define the value of \(\tau\) for each sample \(i\) as \(\tau_{i}\). Then we can define our null and alternative hypotheses: \begin{subequations} \label{eq:hypoModelsConfounders} \begin{align} H_{0}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{S}_{i} + \beta^{TB}_{\tau_{i}b} \label{eq:hypoModelsConfounders0} \\ H_{1}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{S}_{i} + \beta^{TB}_{\tau_{i}b} + \beta^{CB}_{\rho_{i}b} \label{eq:hypoModelsConfounders1} \end{align} \end{subequations} %------------------------------------------------------------------- \subsection{Parameter Estimation} \label{sec:methods-est} %------------------------------------------------------------------- While the described statistical model is robust, efficient, and powerful, it cannot be used to effectively estimate the size of the differential effect or to produce informative parameter estimates. For the purposes of estimating expression levels and the strength of differential effects, we create an entirely separate set of generalized linear models. For the purposes of hypothesis testing this model would be less powerful than the model used in section \ref{sec:methods-htest}, but has the advantage of producing more intuitive and interpretable parameter estimates and fitted values. As before, we generate one set of generalized linear models per splice junction locus: \begin{align} H_{E}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{C}_{\rho{i}} + \beta^{CB}_{\rho_{i}b} \end{align} These models appear similar to the models used for parameter estimation used by DEXSeq, however they differ considerably in that JunctionSeq models each feature separately, rather than attempting to include all features for a given gene in a single model. This "big model" method was originally used for both hypothesis tests and parameter estimation in DEXSeq, but later versions only use this framework for parameter estimation, as the "big model" was found to be less efficient and less consistent. This "big model" methodology also has a number of other weaknesses, including over-weighting of loci with a large number of alternative exonic regions. Reads covering such regions will be "double counted" several times over as they cover multiple regions. For all intents and purposes, DEXSeq counts such double counts as independent observations. JunctionSeq does not count any read or read-pair more than once in any model (unless alternative methodologies are selected, see Section \ref{sec:methods-alts}). Using linear contrasts, the parameter estimates \(\hat{\beta}\), \(\hat{\beta}^{B}_{b}\), \(\hat{\beta}^{C}_{\rho{i}}\), and \(\hat{\beta}^{CB}_{\rho_{i}b}\) can be used to calculate estimates of the effect size (fold change), as well as the mean normalized coverage over the sfeature for each condition. Additionally, the "relative" expression levels can be calculated for each condition, which indicates the expression of the splice junction, normalized relative to the overall gene-wide expression (which may be differentially expressed). This model can be extended to include confounding variables in a manner similar to how the hypothesis test models can be extended (as in Equation \ref{eq:hypoModelsConfounders}). \begin{align} H_{E}: \qquad \log(\mu_{bi}) &= \beta + \beta^{B}_{b} + \beta^{C}_{\rho{i}} + \beta^{TB}_{\tau_{i}b} + \beta^{CB}_{\rho_{i}b} \end{align} %------------------------------------------------------------------- \section{For Advanced Users: Optional Alternative Methodologies} \label{sec:methods-alts} %------------------------------------------------------------------- Since development began on JunctionSeq, a number of new methods have been developed for modelling expression in RNA-Seq data. Many of these more-advanced methods have been integrated into JunctionSeq. However, some advanced users with strong statistical backgrounds may prefer the older, more proven methods. Thus, the legacy methods are still optionally available to JunctionSeq users. In addition, JunctionSeq includes a number of small methodological changes relative to the DEXSeq software upon which it was based. Any or all of these changes can be reverted as well, if desired. All of these optional alternatives can be selected using the "method" parameters passed to \texttt{runJunctionSeqAnalyses} and its sub-functions (see Section \ref{sec:DJUtestAdv}). These optional parameters are: \noindent Additionally, there are a number of options that alter the fundamental methodology to be used. In general the defaults are sufficient, but advanced users may prefer alternative methods. \begin{itemize} \item \texttt{method.GLM} The default is "advanced" or, equivalently, "DESeq2-style". This uses the dispersion estimation methodology used by DESeq2 and DEXSeq v1.12.0 or higher to generate the initial (feature-specific) dispersion estimates. The alternative method is "simpleML" or, equivalently, "DEXSeq-v1.8.0-style". This uses a simpler maximum-likelihood-based method used by the original DESeq and by DEXSeq v1.8.0 or less. \item \texttt{method.dispFit} Determines the method used to generated "fitted" dispersion estimates. One of "parametric" (the default), "local", or "mean". See the DESeq2 documentation for more information. \item \texttt{method.dispFinal} Determines the method used to arrive at a "final" dispersion estimate. The default, "shrink" uses the maximum a posteriori estimate, combining information from both the fitted and feature-specific dispersion estimates. This is the method used by DESeq2 and DEXSeq v1.12.0 and above. \item \texttt{method.sizeFactors} Determines the method used to calculate normalization size factors. By default JunctionSeq uses gene-level expression. As an alternative, feature-level counts can be used as they are in DEXSeq. In practice the difference is almost always negligible. \item \texttt{method.countVectors} Determines the type of count vectors to be used in the model framework. By default JunctionSeq compares the counts for a specific feature against the counts across the rest of the gene minus the counts for the specific feature. Alternatively, the sum of all other features on the gene can be used, like in DEXSeq. The advantage to the default JunctionSeq behavior is that no read or read-pair is ever counted more than once in any model. Under DEXSeq, some reads may cover many exonic segments and thus be counted repeatedly. \item{ \texttt{method.expressionEstimation} Determines the methodology used to generate feature expression estimates and relative fold changes. By default each feature is modeled separately. Under the default count-vector method, this means that the resultant relative fold changes will be a measure of the relative fold change between the feature and the gene as a whole. Alternatively, the "feature-vs-otherFeatures" method builds a large, complex model containing all features belonging to the gene. The coefficients for each feature are then "balanced" using linear contrasts weighted by the inverse of their variance. In general we have found this method to produce very similar results but less efficiently and less consistently. Additionally, this alternative method "multi-counts" reads that cover more than one feature. This can result in over-weighting of exonic regions with a large number of annotated variations in a small genomic area, as each individual read or read-pair may be counted many times in the model. Under the default option, no read or read-pair is ever counted more than once in a given model. } \item \texttt{fitDispersionsForExonsAndJunctionsSeparately} When running a junctionsAndExons-type analysis in which both exons and splice junctions are being tested simultaniously, this parameter determines whether a single fitted dispersion model should be fitted for both exons and splice junctions, or if separate fitted dispersions should be calculated for each. By default the dispersion fits are run separately. \end{itemize} %------------------------------------------------------------------- \subsection{Replicating DEXSeq Analysis via JunctionSeq} \label{sec:methods-dexreplicate} %------------------------------------------------------------------- Using the alternative methods options described in Section \ref{sec:methods-alts}, it is possible to run a standard DEXSeq analysis in JunctionSeq, replicating precisely the output that would be produced by DEXSeq on the given dataset. <>= jscs.DEX <- runJunctionSeqAnalyses(sample.files = countFiles, sample.names = decoder$sample.ID, condition = factor(decoder$group.ID), flat.gff.file = gff.file, nCores = 1, analysis.type = "exonsOnly", method.countVectors = "sumOfAllBinsForGene", method.sizeFactors = "byCountbins", method.expressionEstimation = "feature-vs-otherFeatures", meanCountTestableThreshold = "auto", use.multigene.aggregates = TRUE, method.cooksFilter = TRUE, optimizeFilteringForAlpha = 0.1 ); @ Output files and plotting can be generated using the normal syntax. <>= writeCompleteResults(jscs.DEX, outfile.prefix="./testDEX", save.jscs = TRUE ); buildAllPlots(jscs=jscs.DEX, outfile.prefix = "./plotsDEX/", use.plotting.device = "png", FDR.threshold = 0.01 ); @ Note that these will be formatted and organized in the JunctionSeq style, including the numerous improvements to the plots, along with the browser tracks designed for use with IGC or the UCSC genome browser. Using the above method, we can precisely reproduce the exact analyses that would be generated by the following commands in DEXSeq v1.12.1: <>= library("DEXSeq"); countFiles.dexseq <- system.file(paste0("extdata/cts/", decoder$sample.ID, "/QC.exonCounts.formatted.for.DEXSeq.txt.gz"), package="JctSeqData", mustWork=TRUE); gff.dexseq <- system.file("extdata/annoFiles/DEXSeq.flat.gff.gz", package="JctSeqData", mustWork=TRUE); dxd <- DEXSeqDataSetFromHTSeq( countFiles.dexseq, sampleData = design, design = ~sample + exon + condition:exon, flattenedfile = gff.dexseq ); dxd <- estimateSizeFactors(dxd); dxd <- estimateDispersions(dxd); dxd <- testForDEU( dxd); dxd <- estimateExonFoldChanges(dxd, fitExpToVar = "condition"); dxr <- results(dxd); write.table(dxr,file="dxr.out.txt"); @ These options allow standard DEXSeq methods to be applied while still providing access to the improved output plots and tables generated by JunctionSeq. Also note: the above analyses use the DEXSeq-formatted counts generated by QoRTs\cite{QoRTs}. The python scripts provided with DEXSeq produce very similar counts, but with a few very minor differences. Firstly: some elements will appear in a different order, as DEXSeq uses unsorted collections and thus the ordering cannot be externally reproduced. This includes the ordering of gene names in "aggregate genes", and the ordering of transcripts in the flattened gff file. Additionally, for un-stranded data DEXSeq does not "flatten" genes that appear in overlapping regions on opposite strands. This only affects a very small number of genes, but can produce unusual results under these rare circumstances, and thus has not been reproduced. For un-stranded data QoRTs aggregates any genes that overlap on opposite strands. %------------------------------------------------------------------- \subsection{Advanced generalized linear modelling} \label{sec:methods-covars} %------------------------------------------------------------------- Since JunctionSeq uses generalized linear modelling, the experimental condition variable need not have only two values. Any categorical variable will do. For example: <>= threeLevelVariable <- c("GroupA","GroupA", "GroupB","GroupB", "GroupC","GroupC"); @ We can then carry out analysis normally: <>= jscs <- runJunctionSeqAnalyses(sample.files = countFiles, sample.names = decoder$sample.ID, condition=factor(threeLevelVariable), flat.gff.file = gff.file, nCores = 1, analysis.type = "junctionsAndExons" ); @ Similarly, we could theoretically add additional covariates to our analysis. Note that it is vitally important that every class have replicates, or else JunctionSeq will be unable to accurately assess the biological variance. As a result, any analysis involving covariates should have at least 8 samples. <>= #Artificially adding additional samples by using two of the samples twice: # (Note: this is purely for use as an example. Never do this.) countFiles.8 <- c(countFiles, countFiles[3],countFiles[6]); #Make up some sample names, conditions, and covariates for these samples: decoder.8 <- data.frame( sample.names = factor(paste0("SAMP",1:8)), condition = factor(rep(c("CASE","CTRL"),each=4)), smokeStatus = factor(rep(c("Smoker","Nonsmoker"),4)) ) print(decoder.8); @ To include covariates in the analysis, you will need to modify the various model formulae used by JunctionSeq: <>= jscs <- runJunctionSeqAnalyses(sample.files = countFiles.8, sample.names = decoder.8$sample.names, condition= decoder.8$condition, use.covars = decoder.8[,"smokeStatus",drop=F], flat.gff.file = gff.file, nCores = 1, analysis.type = "junctionsAndExons", test.formula0 = ~ sample + countbin + smokeStatus : countbin, test.formula1 = ~ sample + countbin + smokeStatus : countbin + condition : countbin, effect.formula = ~ condition + smokeStatus + countbin + smokeStatus : countbin + condition : countbin, geneLevel.formula = ~ smokeStatus + condition ); @ Warning: multivariate analysis using GLMs is an advanced task, and JunctionSeq modelling is already fairly complex. Multivariate GLMs are intended for use by advanced users only. It's quite easy to misunderstand what the models are doing, to do something wrong, or to misinterpret the results. %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \bibliography{JunctionSeq} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Session Information} \label{sec:sessioninfo} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- The session information records the versions of all the packages used in the generation of the present document. <>= sessionInfo() @ %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- \section{Legal} \label{sec:legal} %-------------------------------------------------------------------------------- %-------------------------------------------------------------------------------- This software package is licensed under the GNU-GPL v3. A full copy of the GPL v3 can be found in at inst/doc/gpl.v3.txt This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program. If not, see \url{http://www.gnu.org/licenses/}. Portions of this software (and this vignette) are "United States Government Work" under the terms of the United States Copyright Act. It was written as part of the authors' official duties for the United States Government and thus those portions cannot be copyrighted. Those portions of this software are freely available to the public for use without a copyright notice. Restrictions cannot be placed on its present or future use. Although all reasonable efforts have been taken to ensure the accuracy and reliability of the software and data, the National Human Genome Research Institute (NHGRI) and the U.S. Government does not and cannot warrant the performance or results that may be obtained by using this software or data. NHGRI and the U.S. Government disclaims all warranties as to performance, merchantability or fitness for any particular purpose. \end{document}