%\VignetteIndexEntry{INSPEcT} %\VignettePackage{INSPEcT} %\VignetteDepends{INSPEcT} %\VignetteEngine{knitr::knitr} \documentclass[11pt]{article} \usepackage[utf8]{inputenc} \usepackage{mathtools} <>= BiocStyle::latex() @ \newcommand{\Rmethod}[1]{{\Rfunction{#1}}} \title{INSPEcT - INference of Synthesis, Processing and dEgradation rates in Time-course analysis} \author{Stefano de Pretis} \begin{document} \maketitle \tableofcontents % \section{Introduction} \Rpackage{INSPEcT} provides an R/Bioconductor compliant solution for the study of dynamic transcriptional regulatory processes. Based on RNA- and 4sU-seq data, which can be jointly analyzed thanks to a computational normalization routine, \Rpackage{INSPEcT} determines mRNA synthesis, degradation and pre-mRNA processing rates over time for each gene, genome-wide. Finally, the \Rpackage{INSPEcT} modeling framework allows the identification of gene-level transcriptional regulatory mechanisms, determining which combination of synthesis, degradation and processing rates is most likely responsible for the observed mRNA level over time. \Rpackage{INSPEcT} is based on the estimation of total mRNA levels, pre-mRNA levels (from RNA-seq), synthesis rates and processing rates (from 4sU-seq), and degradation rates from the combined analysis of these two data types. 4sU-seq is a recent experimental technique developed to measure the concentration of nascent mRNA and for the genome-wide inference of gene-level synthesis rates. During a short pulse (typically few minutes), cells medium is complemented with 4sU, a naturally occurring modified uridine that is incorporated within growing mRNA chains with minimal impact on cell viability. The chains which have incorporated the uridine variant (the newly synthesized ones) can be isolated from the total RNA population by biotinylation and purification with streptavidin-coated magnetic beads, followed by sequencing. The the main set of steps in the \Rpackage{INSPEcT} workflow is as follows: \begin{itemize} \item Exonic and intronic RPKM for both RNA- and 4sU-seq are determined for each gene. Exonic and intronic RNA-seq RPKM allow to quantify the total mRNA and pre-mRNA, respectively (\Rmethod{makeExonsGtfFromDb}, \Rmethod{makeIntronsGtfFromDb}, \Rfunction{makeRPKMs}). \item Normalized synthesis, processing and degradation rates are obtained by integrating RNA- and 4sU-seq data (\Rmethod{newINSPEcT}). \item Rates, total mRNA and pre-mRNA concentrations are modeled for each gene to assess which of the rates, if any, determined changes in mRNA levels (\Rmethod{modelRates}). \item Simulated data that recapitulate rate distributions, their variation over time and their pair-wise correlations are created and used to evaluate the performance of the method (\Rmethod{makeSimModel}, \Rmethod{makeSimDataset}, \Rmethod{rocCurve}). \end{itemize} Within this vignette, a complete \Rpackage{INSPEcT} analysis is presented (Figure \ref{fig:inspect_workflow}). For details regarding \Rpackage{INSPEcT} extended methods description, refer to de Pretis S. et al., Bioinformatics (2015). \begin{figure} \centering \includegraphics[scale=0.75]{pipeline.png} \caption{Representation of the INSPEcT workflow} \label{fig:inspect_workflow} \end{figure} % \section{Quantification of Exon and Intron features} <>= library(INSPEcT) @ The INSPEcT framework includes function to quantify Exon and Intron features, both in 4sU- and RNA-seq experiments. \Rmethod{makeRPKMs} function builds an annotation for exons and introns, either at the level of transcripts or genes, and count reads falling on the two different annotations in each provided file in BAM or SAM format. This function prioritize the exon annotation, meaning that reads which fall on exon are not counted for introns, in case of overlap. Canonical RPKMs are then evaluated in order to provide expression values. RPKMs, read counts, and the annotation is returned as output. <>= require(TxDb.Mmusculus.UCSC.mm9.knownGene) txdb <- TxDb.Mmusculus.UCSC.mm9.knownGene paths_4su <- system.file('extdata', '4sURNA_0h.bam', package="INSPEcT") paths_total <- system.file('extdata', 'totalRNA_0h.bam', package="INSPEcT") makeRPKMsOut <- makeRPKMs(txdb, paths_4su, paths_total) rpkms <- makeRPKMsOut$rpkms @ In case intronic and exonic RPKMs have been computed using alternative methods, they can be provided for the following analyses. % \section{Time-course analysis} % \subsection{Estimation of rates} \Rpackage{INSPEcT} is based on a simple model of differential equations that describes the process of synthesis and processing of pre-mRNA and the degradation of mature mRNA. Equations model the synthesis of new pre-mRNAs which then decay into mature mRNAs, which in turn exponentially degrade and are removed from the system. The model is based on two main assumptions that are widely used in the description of mRNA life cycle: pre-mRNAs are not degraded, and translocation of mRNAs from nucleus to cytoplasm occurs immediately after maturation, or at a rate considerably faster than the rate of degradation (Rabani M. et al., Nature Biotechnology, 2011; Sun M. et. al, Genome Research, 2012). The model lacks any spatial assumption, like segregation of mRNAs into cellular compartments that could impact the degradation rate, but this is a consequence of the non-spatial nature of the data:\\ \begin{equation} \nonumber \left\{ \begin{array}{l l} \dot{P} = a(t) - c(t) \, P \\ \dot{T} = a(t) - b(t) \, (T - P) \end{array} \right. \end{equation} where $T$ is total RNA, $P$ is pre-mature RNA, $a(t)$ is the synthesis rate, $b(t)$ is the degradation rate and $c(t)$ is the processing rate. After having quantified data from RNA-seq ($R$) and 4sU-seq (labeled, $L$) libraries into intronic and exonic RPKMs (\Rfunction{makeRPKMs} function), the \Rmethod{newINSPEcT} method is used to estimate synthesis, processing and degradation rates by solving the above system of differential equations applied at every time point ($t$) to both the total and labeled fractions. When applied to the labeled RNA fraction, the system can be solved and integrated between $t-t_L$ and $t$, assuming that no labeled molecules existed before the labeling pulse ($t_L$). At each time point, \Rpackage{INSPEcT} solves a system of four equations: \begin{equation} \nonumber \left\{ \begin{array}{l l} \dot{P}_{R_t} = a_t - c_t \, P_{R_t} \\ \dot{T}_{R_t} = a_t - b_t \, (T_{R_t} - P_{R_t})\\ P_{L_t} = \frac{a_t}{c_t} - ( 1 - e^{c_t \, t_L} ) \\ T_{L_t} = a_t \, t_L \end{array} \right. \end{equation} with three unknowns $a_t$, $b_t$ and $c_t$ which are respectively the synthesis, degradation and processing rates at time $t$. $P_{R_t}$ is equal to the pre-mRNA level (intronic RNA-seq RPKM), $T_{R_t}$ is equal to the total mRNA level (exonic RNA-seq RPKM), and $P_{L_t}$ is equal to the pre-mRNA level as quantified in the labeled fraction (intronic 4sU-seq RPKM). Finally, $T_{L_t}$ is equal to the total mRNA level as quantified in the labeled fraction (exonic 4sU-seq RPKM). $\dot{P}_{R_t}$, $\dot{T}_{R_t}$ are estimated from the interpolation of $P_R(t)$ and $T_R(t)$. The overdetermination of the system is used to calculate a time-point specific scaling factor between RNA- and 4sU-seq RPKMs that can be visualized using the \Rmethod{sfPlot} method.\\ In order to create an object of class \Rclass{INSPEcT} and to calculate the first estimates of rates and concentrations for each gene, the specific time points of the time-course, the 4sU labeling time and the RPKMs corresponding to labeled and total fraction of the intronic and exonic regions of each gene have to be provided (here time is in $hours$). All results are stored in an object of class \Rclass{INSPEcT} and rates can be accessed with the \Rmethod{ratesFirstGuess} method. This sample code calculates rates and concentrations on a sample set of 500 genes: <>= data('rpkms', package='INSPEcT') tpts <- c(0, 1/6, 1/3, 1/2, 1, 2, 4, 8, 16) tL <- 1/6 mycerIds <- newINSPEcT(tpts, tL, rpkms$foursu_exons, rpkms$total_exons, rpkms$foursu_introns, rpkms$total_introns, BPPARAM=SerialParam()) head(ratesFirstGuess(mycerIds, 'synthesis')) @ In case of a long 4sU labeling time (longer than 10-15 minutes), it could be useful to activate the \Rcode{degDuringPulse} option, as shown below. This option estimates all the rates of the RNA life-cycle without assuming that no degradation of the newly synthesized transcripts occurs during the pulse. The longer the labeling time is, the weaker this assumption gets. This option, however, involves solving a more complicated system of differential equations and for this reason it is not recommended for short labeling times. <>= mycerIdsDdp <- newINSPEcT(tpts, tL, rpkms$foursu_exons, rpkms$total_exons, rpkms$foursu_introns, rpkms$total_introns, BPPARAM=SerialParam(), degDuringPulse=TRUE) head(ratesFirstGuess(mycerIdsDdp, 'synthesis')) @ It is possible to subset the \Rpackage{INSPEcT} object and focus on a specific set of genes. For the sake of speeding up the downstream analysis, we are going to focus on the first 10 genes of the obtained \Rpackage{INSPEcT} object. We can now display the total and pre-mRNA concentrations together with the synthesis, degradation and processing pre-modeling rates they originated from using the \Rmethod{inHeatmap} method (Figure \ref{fig:pre_model_heatmap}). \begin{center} <>= mycerIds10 <- mycerIds[1:10] inHeatmap(mycerIds10, clustering=FALSE) @ \end{center} % \subsection{Modeling of rates to determine transcriptional regulatory mechanism} Once a prior estimate is obtained for synthesis, processing and degradation rates over time for each gene, \Rpackage{INSPEcT} tests different models of transcriptional regulation to identify the most likely combination of rates explaining the observed changes in gene expression (\Rmethod{modelRates} method). To this purpose, a parametric function is fit to each rate over time, through minimization of residual sum of squares. Once the parametric functionalization for synthesis, degradation and processing rates are obtained, it is possible to test how those parametric functions recapitulate the experimental data they originated from after an additional minimization step. To identify the most likely mechanism of transcriptional regulation, \Rpackage{INSPEcT} tests the possibility that each rate is constant during the time course by building models that alternatively set as constant one, two or all the three rates. Due to the fact the the initial parameters for each rate function are initialized randomly, the \Rcode{seed} argument in \Rmethod{modelRates} can be set to obtain reproducible results. <>= mycerIds10 <- modelRates(mycerIds10, seed=1) @ Following this modeling procedure, new rates are computed and they can be accessed through the \Rmethod{viewModelRates} method and visualized with the \Rmethod{inHeatmap} method (Figure \ref{fig:post_model_heatmap}). \begin{center} <>= data('mycerIds10', package='INSPEcT') head(viewModelRates(mycerIds10, 'synthesis')) inHeatmap(mycerIds10, type='model', clustering=FALSE) @ \end{center} The \Rmethod{geneClass} method can be used to recapitulate the transcriptional regulatory mechanism assigned to each modeled gene. In particular, each gene is assigned to a class named after the set of varying rates, if any ("0" denotes a gene whose rates are constant over time, "a" denotes a gene whose synthesis changes over time, "b" denotes a gene whose degradation changes over time, "c" denotes a gene whose processing changes over time). <>= geneClass(mycerIds10) @ The \Rmethod{plotGene} method can be used to investigate profiles of mRNA concentrations and rates for a given gene. Estimated synthesis, degradation and processing rates, pre-mRNA and total mRNA concentrations are displayed with solid thin lines, while their variances are in dashed lines and the modeled rates and concentrations are in thick solid lines. This example shows a gene of class "ab", indicating that its levels are controlled by both synthesis and degradation. In this case, the variation of multiple rates determines opposite trends for total and pre-mRNA concentrations over time: despite a decrease in the synthesis rate over time, the levels of total mRNA increase mostly due to a decrease in the degradation rate (Figure \ref{fig:gene_plot}). \begin{center} <>= plotGene(mycerIds10, 2, fix.yaxis=TRUE) @ \end{center} For each model, the chi-squared statistic that measures the goodness of the fit is calculated. In order to evaluate how good the models are able to recapitulate the experimental data, the chi-squared test p-values of the model that better represents the transcriptional scenario for each gene can be visualized as a histogram (Figure \ref{fig:goodness_of_fit}). Eventually, models with a p-value of the chi-squared test higher than a selected threshold can be discarded. \begin{center} <>= chisq <- chisqmodel(mycerIds10) hist(log10(chisq), main='', xlab='log10 chi-squared p-value') discard <- which(chisq>.05) featureNames(mycerIds10)[discard] mycerIds10new <- mycerIds10[-discard] @ \end{center} % \subsection{Evaluation of performance via simulated data} \Rpackage{INSPEcT} provides functionality to build a synthetic dataset for which the transcriptional scenario is known for each gene. Simulated data can be used to evaluate the performance of \Rpackage{INSPEcT} in classifying each rate as constant or variable over time, and to estimate the number of time points and replicates necessary to achieve a given performance. The method \Rmethod{makeSimModel} takes as arguments an \Rpackage{INSPEcT} object and the number of genes that have to be sampled. The \Rpackage{INSPEcT} object is used to sample absolute values of the rates, fold changes, correlations between absolute values and fold changes and variance of the noise to be added to each feature. Optionally, the user can provide a new set of time points where the synthetic dataset will be sampled, and the ratio between the constant rates, the rates modeled with an impulse model and the rates modeled with a sigmoid function. By default the ratio is .5 constant, .3 impulse, .2 sigmoid. <>= simRates <- makeSimModel(mycerIds, 1000, newTpts=NULL, seed=1) @ The \Rmethod{makeSimModel} method generates an object of class \Rclass{INSPEcT\_model} which can be used to generate an object of class \Rclass{INSPEcT} using the \Rmethod{makeSimDataset} method. The \Rmethod{makeSimDataset} method takes as arguments the timepoints at which the dataset has to be simulated and the number of replicates that need to be simulated. The object created by this method can be modeled via \Rmethod{modelRates} as any other object of class \Rclass{INSPEcT}. <>= simData1rep <- makeSimDataset(simRates, tpts, 1, seed=1) simData1rep <- modelRates(simData1rep, seed=1) newTpts <- c(0, 1/6, 1/3, 1/2, 1, 1.5, 2, 4, 8, 12, 16, 24) simData3rep <- makeSimDataset(simRates, newTpts, 3, seed=1) simData3rep <- modelRates(simData3rep, seed=1) @ It is possible to compare now the performance of the modeling, by comparing the \Robject{simRates} object, which contains the ground truth of rates, to the \Robject{simData1rep} or \Robject{simData3rep} objects, which contain the predictions made by \Rpackage{INSPEcT} on datasets that have been simulated with one replicate of 9 time points or three replicates of 12 time points. Modeled data have been previously computed and stored within the package for computational time reasons and are not evaluated directly within this vignette. The evaluation of the performance is done using a ROC curve analysis and measured with the area under the curve (AUC) (Figure \ref{fig:simulated_data_run}). \begin{center} <>= data('simRates', package='INSPEcT') data('simData1rep', package='INSPEcT') data('simData3rep', package='INSPEcT') par(mfrow=c(1,2)) rocCurve(simRates, simData1rep); title("1 replicate - 9 time points", line=3) rocCurve(simRates, simData3rep); title("3 replicates - 12 time points", line=3) @ \end{center} The method \Rmethod{rocThresholds} can be used to assess the sensitivity and specificity that is achieved thanks to the given thresholds for the chi-squared test and for the Brown's test. If thresholds are not provided, default values are used (Figure \ref{fig:simulated_data_run2}). \begin{center} <>= rocThresholds(simRates, simData3rep, bTsh=c(.01,.01,.05), cTsh=.1) @ \end{center} % \subsection{Parameter settings} If desired, different parameters can be set for both the modeling and the testing part. Regarding the modeling part, we might want to exclude testing sigmoid functions (all evaluated smooth function will be impulse model functions, \Rcode{useSigmoidFun=FALSE} option), increase the number of different initializations that are performed for each gene (\Rcode{nInit} option), or increase the maximum number of steps in the rates optimization process (\Rcode{nIter} option). All these choices could increase the performance of the method, but also the needed computational time. Nevertheless, the use of sigmoid functions can reduce over-fitting problems and is highly recommended when the number of data points within the time-course is lower than 7. The impact of these options can be evaluated using a synthetic dataset. <>= modelingParams(mycerIds10)$useSigmoidFun <- FALSE modelingParams(mycerIds10)$nInit <- 20 modelingParams(mycerIds10)$nIter <- 1000 @ Alternatively we might want to change the thresholds for chi-squared and log-likelihood ratio tests, or define the specific set of models to be compared with log-likelihood ratio test while assessing if a given rate is variable or not. In this example, we are changing the thresholds of both the chi-squared test and the Brown's method for combining p-values. Regarding the processing rates, only the models in which all rates are constant ("0") will be compared to the one in which only processing varies ("c") to assess the variability of the rates using log-likelihood ratio test. <>= thresholds(mycerIds10)$chisquare <- .1 thresholds(mycerIds10)$brown <- c(alpha=.01, beta=.01, gamma=.05) llrtests(mycerIds10)$processing <- list(c('0','c')) @ To have a sense of all parameters that can be set, type: <>= ## modeling modelingParams(mycerIds10) ## model selection and testing framework modelSelection(mycerIds10) thresholds(mycerIds10) llrtests(mycerIds10) @ % \section{Steady-state analysis} Synthesis, processing and degradation rates are the determinants of the levels of pre-mRNAs and mature mRNAs. At steady-state the ratio between synthesis and processing rates determine the pre-mRNA levels, while the ratio between synthesis and degradation determine the mature mRNA levels. Different conditions might express different levels of mature mRNA or pre-mRNA for each gene and this difference could arise from the different usage of one, two or all three rates. In order to address the problem, rates from the two different experimental conditions are compared at the gene level by a simple t-test between rate values. Rate variances are estimated using the properties of variance using the approximation that synthesis rate ($a$), degradation rate ($b$) and processing rates ($c$) are computed as: \begin{equation} \nonumber \begin{array}{l l} a = \, \frac{T_{L}}{t_L \cdot s_f} \\ b = \, \frac{a}{T_{T} - P_{T}} \\ c = \, \frac{c}{P_{T}} \end{array} \end{equation} where $T_{L}$ is the 4sU-seq-exonic RPKM of the gene under consideration, $T_{T}$ is the RNA-seq-exonic RPKM, $P_{T}$ is the RNA-seq-intronic RPKM, $t_L$ is the labeling time and $s_f$ is a scaling factor that INSPEcT calculated to linearly scale the 4sU library. Considered that we know the variance of $T_{L}$, $T_{T}$, and $P_{T}$, we can also infer the variance of $a$, $b$ and $c$. In order to test log scaled means and variances, we transform the variance using the formula: \begin{equation} \nonumber \begin{array}{l l} Var[f(X)] \approx \, (f'(E[x]))^2 \cdot Var(X) \end{array} \end{equation} In order to exemplify the steady state analysis task, we generate a second set of simulated data with 3 replicates from the object \Rclass{INSPEcT\_model} simRates and compare it with the one which was previously generated for the time course analysis. <>= newTpts <- c(0, 1/6) simData3rep_2 <- makeSimDataset(simRates, newTpts, 3, seed=2) @ We have now two datasets with 3 replicates. The comparison between the two datasets is performed by the \Rmethod{compareSteady} method. This method take as input two objects of class INSPEcT and check that each of have been profiled with replicates, at least for the steady state condition (identified as the first temporal condition). <>= diffrates <- compareSteady(simData3rep, simData3rep_2) @ Results of the comparison can be seen at a glance by typing the name of the variable which contains the results: <>= diffrates @ The complete datasets regarding synthesis, processing and degradation can be accessed via the accessor methods: <>= head(synthesis(diffrates)) head(processing(diffrates)) head(degradation(diffrates)) @ A method for plotting results for each rate has been implemented ((Figure \ref{fig:steady5})): <>= plotMA(diffrates, rate='synthesis', alpha=.5) @ Enjoy! % \section{About this document} <>= sessionInfo() @ \end{document}