--- title: > "Finding differentially expressed unannotated genomic regions from RNA-seq data with srnadiff" author: - name: Matthias Zytnicki - name: Ignacio González date: "`r BiocStyle::doc_date()`" package: "`r BiocStyle::pkg_ver('srnadiff')`" bibliography: references.bib vignette: > %\VignetteIndexEntry{The srnadiff package} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} output: BiocStyle::html_document: toc_float: true --- ```{r setup, echo=FALSE, results="hide", warning=FALSE} suppressPackageStartupMessages({ library("srnadiff") library("BiocManager") library("BiocStyle") library("knitr") library("rmarkdown") library("grid") }) knitr::opts_chunk$set( error=FALSE, fig.height=5, fig.width=8, message=FALSE, warning=FALSE, tidy=FALSE ) ``` ```{r echo=FALSE} basedir <- system.file("extdata", package="srnadiff", mustWork=TRUE) sampleInfo <- read.csv(file.path(basedir, "dataInfo.csv")) bamFiles <- paste(file.path(basedir, sampleInfo$FileName), "bam", sep=".") gtfFile <- file.path(basedir, "Homo_sapiens.GRCh38.76.gtf.gz") annotReg <- readAnnotation(gtfFile, feature="gene", source="miRNA") ``` # Version info - **R version**: `r R.version.string` - **Bioconductor version**: `r BiocManager::version()` - **Package version**: `r packageVersion("srnadiff")` # Abstract `r Biocpkg("srnadiff")` is an *R* package that finds differently expressed regions from RNA-seq data at *base-resolution* level without relying on existing annotation. To do so, the package implements the *identify-then-annotate* methodology that builds on the idea of combining two pipelines approach: differential expressed regions detection and differential expression quantification. # Introduction There is no real method for finding differentially expressed short RNAs. The most used method focusses on miRNAs, and only uses a standard RNA-Seq pipe-line on these genes. However, annotated tRF, siRNAs, piRNA, etc. are thus out of the scope of these analyses. Several ad hoc method have been used, and this package implements a unifying method, finding differentially expressed genes or regions of any kind. The `r Biocpkg("srnadiff")` package implements two major methods to produce potential differentially expressed regions: the HMM and IR method. Briefly, these methods identify contiguous base-pairs in the genome that present differential expression signal, then these are regrouped into genomic intervals called differentially expressed regions (DERs). Once DERs are detected, the second step in a *sRNA-diff* approach is to quantify the signification of these. To do so, reads (including fractions of reads) that overlap each expressed region are counted to arrive at a count matrix with one row per region and one column per sample. Then, this count matrix is analyzed using the standard workflow of `r Biocpkg("DESeq2")` for differential expression of RNA-seq data, assigning a p-value to each candidate DER. The main functions for finds differently expressed regions are `srnadiffExp` and `srnadiff`. The first one creates an S4 class providing the infrastructure (slots) to store the input data, methods parameters, intermediate calculations and results of an *sRNA-diff* approach. The second one implement four methods to find candidate differentially expressed regions and quantify the statistic signification of the finded regions. This vignette explains the basics of using `r Biocpkg("srnadiff")` by showing an example, including advanced material for fine tuning some options. The vignette also includes description of the methods behind the package. ## Citing `r Biocpkg("srnadiff")` We hope that `r Biocpkg("srnadiff")` will be useful for your research. Please use the following information to cite `r Biocpkg("srnadiff")` and the overall approach when you publish results obtained using this package, as such citation is the main means by which the authors receive credit for their work. Thank you! ```{r eval=TRUE, echo=FALSE} x <- citation("srnadiff") ``` Zytnicki, M., and I. González. `r paste0("(", x$year, "). ")` "`r paste0(x$title, ". ")`" `r paste0(x$note, ".")` ## How to get help for `r Biocpkg("srnadiff")` Most questions about individual functions will hopefully be answered by the documentation. To get more information on any specific named function, for example `MIMFA`, you can bring up the documentation by typing at the *R*. ```{r eval=FALSE, message=FALSE, warning=FALSE} help("srnadiff") ``` or ```{r eval=FALSE, message=FALSE, warning=FALSE} ?srnadiff ``` The authors of `r Biocpkg("srnadiff")` always appreciate receiving reports of bugs in the package functions or in the documentation. The same goes for well-considered suggestions for improvements. If you've run into a question that isn't addressed by the documentation, or you've found a conflict between the documentation and what the software does, then there is an active community that can offer help. Send your questions or problems concerning `r Biocpkg("srnadiff")` to the Bioconductor support site at \url{https://support.bioconductor.org}. Please send requests for general assistance and advice to the support site, rather than to the individual authors. It is particularly critical that you provide a small reproducible example and your session information so package developers can track down the source of the error. Users posting to the support site for the first time will find it helpful to read the posting guide at [the Bioconductor help page](http://www.bioconductor.org/help/support/posting-guide). ## Quick start A typical *sRNA-diff* session can be divided into three steps: 1. **Data preparation:** In this first step, a convenient *R* object of class `srnadiffExp` is created containing all the information required for the two remaining steps. The user needs to provide a vector with the full paths to the BAM files, a `data.frame` with sample and experimental design information and optionally annotated regions as a `GRanges` object. 2. **Performing srnadiff:** Using the object created in the first step the user can perform `srnadiff` to find potential DERs and quantify the statistic signification of these. 3. **Visualization of the results:** The DERs obtained in the second step are visualized by plotting the coverage information surrounding genomic regions. A typical `r Biocpkg("srnadiff")` session might look like the following. Here we assume that `bamFiles` is a vector with the full paths to the BAM files and the sample and experimental design information are stored in a data frame `sampleInfo`. ```{r message=FALSE, warning=FALSE, include=FALSE, results="hide"} #-- Data preparation srnaExp <- srnadiffExp(bamFiles, sampleInfo) #-- Performing srnadiff srnaExp <- srnadiff(srnaExp) #-- Visualization of the results plotRegions(srnaExp, regions(srnaExp)[1]) ``` # Using `r Biocpkg("srnadiff")` ## Installation We assume that the user has the *R* program (see the [*R* project](http://www.r-project.org)) already installed. The `r Biocpkg("srnadiff")` package is available from the [Bioconductor repository](http://www.bioconductor.org). To be able to install the package one needs first to install the core `r CRANpkg("Bioconductor")` packages. If you have already installed `r CRANpkg("Bioconductor")` packages on your system then you can skip the two lines below. ```{r, bioconductor, eval=FALSE} if (!requireNamespace("BiocManager", quietly=TRUE)) install.packages("BiocManager") ``` Once the core `r CRANpkg("Bioconductor")` packages are installed, you can install the `r Biocpkg("srnadiff")` package by ```{r, install, eval=FALSE} BiocManager::install("srnadiff", version="3.8") ``` Load the `r Biocpkg("srnadiff")` package in your *R* session: ```{r, loadLibrary, eval=FALSE} library(srnadiff) ``` A list of all accessible vignettes and methods is available with the following command: ```{r, helpSearch, eval=FALSE} help.search("srnadiff") ``` ## Data overview To help demonstrate the functionality of `r Biocpkg("srnadiff")`, the package includes datasets of published by [@data]. Briefly, these data consist of three replicates of sRNA-Seq of SLK (human) cell lines, and three replicates of SLK cell lines infected with Kaposi's sarcoma associated herpesvirus. The analysis shows that several loci are repressed in the infected cell lines, including the 14q32 miRNA cluster. Raw data have been downloaded from the GEO data set [GSE62830](https://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE62830). Adapters were removed with [fastx\_clipper](http://hannonlab.cshl.edu/fastx_toolkit/commandline.html) and mapped with [@bowtie2] on the human genome version GRCh38. This data is restricted to a small locus on chr14. It uses the whole genome annotation (with coding genes, etc.) and extracts miRNAs. The file `dataInfo.csv` contains three columns for each BAM file: * the file name (`FileName`) * a human readable sample name (`SampleName`) * a condition, e.g. `WT` (`Condition`) ## Data preparation: the `srnadiffExp` object The first step in an *sRNA-diff* approach is to create a `srnadiffExp` object. `srnadiffExp` is an `S4` class providing the infrastructure to store the input data, methods parameters, intermediate calculations and results of a *sRNA-diff* approach. This object will be also the input of the visualization function. The object `srnadiffExp` will usually be represented in the code here as `srnaExp`. To build such an object the user needs the following: 1. **paths to the BAM files:** a vector with the full paths to the sample BAM files. 2. **sample information:** a `data.frame` with three columns labelled `FileName`, `SampleName` and `Condition`. The first column is the BAM file name (without extension), the second column the sample name, and the third column the condition to which sample belongs. Each row describes one sample. 3. **annotation:** optionally annotation information. A `r Biocpkg("GRanges")` object containing annotated regions. Here, we demonstrate how to construct a `srnadiffExp` object from [@data]. ```{r message=FALSE, warning=FALSE} ## Determine the path to data files basedir <- system.file("extdata", package="srnadiff", mustWork=TRUE) ## Vector with the full paths to the BAM files to use bamFiles <- paste(file.path(basedir, sampleInfo$FileName), "bam", sep=".") ## Reads sample information file and creates a data frame from it sampleInfo <- read.csv(file.path(basedir, "dataInfo.csv")) ## Vector with the full paths to the BAM files to use bamFiles <- paste(file.path(basedir, sampleInfo$FileName), "bam", sep = ".") ## Creates an srnadiffExp object srnaExp <- srnadiffExp(bamFiles, sampleInfo) ``` Optionally, if annotation information is available as a `r Biocpkg("GRanges")` object, `annotReg` say, then a `srnadiffExp` object can be created by ```{r, eval=FALSE} srnaExp <- srnadiffExp(bamFiles, sampleInfo, annotReg) ``` or by ```{r message=FALSE, warning=FALSE} srnaExp <- srnadiffExp(bamFiles, sampleInfo) annotReg(srnaExp) <- annotReg ``` A summary of the `srnaExp` object can be seen by typing the object name at the *R* prompt ```{r} srnaExp ``` For your conveniance and illustrative purposes, an example of an `srnadiffExp` object can be loaded with an only command, so the script boils down to: ```{r} srnaExp <- srnadiffExample() ``` The `srnadiffExp` object in this example was constructed by: ```{r message=FALSE, warning=FALSE} basedir <- system.file("extdata", package="srnadiff", mustWork=TRUE) sampleInfo <- read.csv(file.path(basedir, "dataInfo.csv")) gtfFile <- file.path(basedir, "Homo_sapiens.GRCh38.76.gtf.gz") annotReg <- readAnnotation(gtfFile, feature="gene", source="miRNA") bamFiles <- paste(file.path(basedir, sampleInfo$FileName), "bam", sep=".") srnaExp <- srnadiffExp(bamFiles, sampleInfo, annotReg) ``` ## Read annotation `r Biocpkg("srnadiff")` offers the `readAnnotation` function related to loading annotation data. This accepts two annotation format files: GTF and GFF formats. Specification of GTF/GFF format can be found at the [UCSC dedicated page](https://genome.ucsc.edu/FAQ/FAQformat.html). `readAnnotation` reads and parses content of GTF/GFF files and stores annotated genomic features (regions) in a `GRanges` object. This has three main arguments: the first argument indicates the path, URL or connection to the GTF/GFF annotation file. The second and third argument, `feature` and `source` respectively, are of type character string, these specify the feature and attribute type used to select rows in the GTF/GFF annotation which will be imported. `feature` and `source` can be `NULL`, in this case, no selection is performed and all content into the file is imported. ### Extraction of putative regions using an GTF annotation file This method simply provides the genomic regions corresponding to the annotation file that is optionally given by the user. It can be a set of known miRNAs, siRNAs, piRNAs, genes, or a combination thereof. ### Whole genome file annotation This GTF file can be found in the central repositories ([NCBI](https://www.ncbi.nlm.nih.gov/), [Ensembl](https://www.ensembl.org/index.html)) and contains all the annotation found in an organism (coding genes, tranposable element, etc.). The following function reads the annotation file and extracts the miRNAs. Annotation files may have different formats, but this command has been tested on several model organisms (including human) from Ensembl. ```{r message=FALSE, warning=FALSE} gtfFile <- file.path(basedir, "Homo_sapiens.GRCh38.76.gtf.gz") annotReg <- readAnnotation(gtfFile, feature="gene", source="miRNA") ``` ### Extraction of precursor miRNAs using a miRBase-formatted GFF file [miRBase](http://www.mirbase.org/) [@mirbase] is the central repository for miRNAs. If your organism is available, you can download their miRNA annotation in GFF3 format (check the "Browse" tab). The following code parses a GFF3 miRBase file, and extracts the precursor miRNAs. ```{r message=FALSE, warning=FALSE} gffFile <- file.path(basedir, "mirbase21_GRCh38.gff3") annotReg <- readAnnotation(gffFile, feature="miRNA_primary_transcript") ``` ### Extraction of mature miRNAs using a miRBase-formatted GFF file In the previous example, the reads will be counted per pre-miRNA, and the 5' and 3' arms, the miRNA and the miRNA\* will be merged in the same feature. If you want to separate the two, use: ```{r message=FALSE, warning=FALSE} gffFile <- file.path(basedir, "mirbase21_GRCh38.gff3") annotReg <- readAnnotation(gffFile, feature="miRNA") ``` ### Other format When the previous functions do not work, you can use your own parser with: ```{r message=FALSE, warning=FALSE} annotation <- readAnnotation(gtfFile, source="miRNA", feature="gene") ``` The `source` parameter keeps all the lines such that the second field matches the given parameter (e.g. `miRNA`). The `feature` parameter keeps all the lines such that the third field matches the given parameter (e.g. `gene`). The name of the feature will be given by the tag `name` (e.g. `gene_name`). `source`, `feature` and `name` can be `NULL`. In this case, no selection is performed on `source` or `feature`. If `name` is null, then a systematic name is given (`annotation_N`). ## Performing *sRNA-diff* The main function for performing an *sRNA-diff* analysis is `srnadiff`, this the wrapper for running several key functions from this package. `srnadiff` implement four methods to produce potential DERs: the *annotation*, *naive*, *hmm* and *IR* method (see bellow). Once potential DERs are detected, the second step in `srnadiff` is to quantify the statistic signification of these. `srnadiff` has three main arguments. The first argument is an instance of class `srnadiffExp`. The second argument is of type character vector, it specify the segmentation methods to use, one of `annotation`, `naive`, `hmm`, `IR` or combinations thereof. The default `all`, all methods are used. The third arguments is of type list, it contain named components for the methods parameters to use. If missing, default parameter values are supplied. Details about the methods parameters are further described in the manual page of the `parameters` function and in *Methods to produce differentially expressed regions* section. We then performs an *sRNA-diff* analysis on the input data contained in `srnaExp` by ```{r} srnaExp <- srnadiff(srnaExp) ``` `srnadiff` returns an object of class `srnadiffExp` again containing additional slots for: * `regions` * `parameters` * `countMatrix` # Working with the `srnadiffExp` object Once the `srnadiffExp` object is created the user can use the methods defined for this class to access the information encapsulated in the object. By example, the sample information is accessed by ```{r} sampleInfo(srnaExp) ``` For accessing the `chromosomeSize` slot ```{r} chromosomeSizes(srnaExp) ``` The list of parameters can be exported by the function `parameters` ```{r} parameters(srnaExp) ``` ## Extracting regions The regions, with corresponding information provided by `r Biocpkg("DESeq2")` (mean expression, fold-change, p-value, adjusted p-value, etc.), can be extracted with this command: ```{r message=FALSE, warning=FALSE} regions <- regions(srnaExp, pvalue=0.5) ``` where `pvalue` is the (adjusted) p-value threshold. The output in a `r Biocpkg("GenomicRanges")` object, and the information is accessible with the `mcols()` function. ## Data visualization An insightful way of looking at the results of `srnadiff` is to investigate how the coverage information surrounding finded regions are distributed on the genomic coordinates. `plotRegions` provides a flexible genomic visualization framework by displaying tracks in the sense of the `Gviz` package. Given a region (or regions), four separate tracks are represented: 1. `GenomeAxisTrack`, a horizontal axis with genomic coordinate tickmarks for reference location to the displayed genomic regions; 2. `GeneRegionTrack`, if the `annot` argument is passed, a track displaying all gene and/or sRNA annotation information in a particular region; 3. `AnnotationTrack`, regions are plotted as simple boxes if no strand information is available, or as arrows to indicate their direction; and 4. `DataTrack`, plot the sample coverages surrounding the genomic regions. The sample coverages can be plotted in various different forms as well as combinations thereof. Supported plotting types are: The sample coverages can be plotted in various different forms as well as combinations thereof. Supported plotting types are: `p`: simple dot plot; `l`: lines plot; `b`: combination of dot and lines plot; `a`: lines plot of the sample-groups average (i.e., mean) values; `confint`: confidence intervals for average values. The default visualization for results from `srnadiff` is a lines plot of the sample-groups average. ```{r, fig.height=3.5, fig.width=4, fig.align='center', out.width='450pt'} plotRegions(srnaExp, regions(srnaExp)[1]) ``` # Methods behind `r Biocpkg("srnadiff")` ## Pre-processing data As input, `srnadiffExp` expects BAM files as obtained, *e.g.*, from RNA-Seq or another high-throughput sequencing experiment. Reading and processing of BAM files uses current `r CRANpkg("Bioconductor")` infrastructure for processing sequencing reads: `r Biocpkg("RSamtools")`, `r Biocpkg("IRanges")` and `r Biocpkg("GenomicRanges")` libraries. At this stage BAM files are summarized into base-resolution coverage and stored in a run-length encoding format in order to enhance computational performance. Run-length encoding is a compact way to store an atomic vector as a pairs of vectors (value, length). It is based on the `rle` function from the base *R* package. As a second pre-processing step, `srnadiffExp` estimate the size factors (the effective library size) from the coverage data, such that count values in each sample coverage can be brought to a common scale by dividing by the corresponding size (normalization) factor. This step is also called normalization, its purpose is to render coverages (counts) from different samples, which may have been sequenced to different depths, comparable. Normalization factors are estimated using the *median ratio method* described by Equation 5 in [@deseq]. ## Methods to produce differentially expressed regions ### HMM method: `hmm` The first step in HMM method is quantifying the evidence for differential expression at the base-resolution level. To do this, `srnadiff` use the common approach in comparative analysis of transcriptomics data: test the null hypothesis that the logarithmic fold change between condition groups for a nucleotide expression is exactly zero. The next step in the HMM approach enforces a smoothness assumption over the state of nucleotides: differential expression does not randomly switch along the chromosome, rather, continuous regions of RNA are either "differentially expressed" or "not". This is captured with a hidden Markov model (HMM) with binary latent state corresponding to the true state of each nucleotide: differentially expressed or not differentially expressed. The observations of the HMM are then the empirical p-values arising from the differential expression analysis corresponding to each nucleotide position. Modelling p-values directly enabled us to define the emission of each state as follows: the differentially expressed state emits a p-value $< t$ with probability $p$, and the not differentially expressed state emits a p-value $\geqslant t$ with probability $1-p$, where $t$ is a real number between 0 and 1. The HMM approach normally needs emission, transition, and starting probabilities values. They can be tuned by the user according to the overall p-values from differential analysis. We then run the Viterbi algorithm [ref] in order to finding the most likely sequence of states from the HMM. This essentially segments the genome into regions, where a region is defined as a set of consecutive bases showing a common expression signature. A region of bases with differentially expressed state is referred as an expressed region and is given as output of the method. To run the HMM approach, `srnadiff` first form a large matrix, with rows corresponding to bases, columns corresponding to samples and entries are the coverage from a nucleotide of a particular sample. This count matrix is then analyzed as into feature-level counts using the feature-level RNA-seq differential expression analysis from `r Biocpkg("DESeq2")`. In practice, the p-value is not computed for every nucleotide. Nucleotides for which the sum of the coverage across all samples is less than a threshold are given a p-value of 1, because these poorly expressed bases are unlikely to provide a differentially expressed sRNA. The parameters for the HMM method are: `noDiffToDiff`: Initial transition probability from "no differentially expressed" state to "differentially expressed" state. `diffToNoDiff`: Initial transition probability from "differentially expressed" state to no "differentially expressed" state. `emission`: Is the probability to emit a p-value $