Title: | Modeling RNA synthesis, processing and degradation with RNA-seq data |
---|---|
Description: | INSPEcT (INference of Synthesis, Processing and dEgradation rates from Transcriptomic data) RNA-seq data in time-course experiments or steady-state conditions, with or without the support of nascent RNA data. |
Authors: | Stefano de Pretis |
Maintainer: | Stefano de Pretis <[email protected]>, Mattia Furlan <[email protected]> |
License: | GPL-2 |
Version: | 1.37.0 |
Built: | 2024-11-18 03:46:11 UTC |
Source: | https://github.com/bioc/INSPEcT |
This method is used to retrieve AIC values for all models tested for all genes.
## S4 method for signature 'INSPEcT_model' AIC(object, ..., k = 2) ## S4 method for signature 'INSPEcT' AIC(object, ..., k = 2)
## S4 method for signature 'INSPEcT_model' AIC(object, ..., k = 2) ## S4 method for signature 'INSPEcT' AIC(object, ..., k = 2)
object |
An object of class INSPEcT or INSPEcT_model |
... |
Additional arguments for the generic |
k |
Additional parameter for the generic |
A matrix of AIC values
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) AIC(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) AIC(nascentInspObj10)
A list containing mature and nascent counts for exons and introns, three replicates and 11 time points: 0,1/6,1/3,1/2,1,1.5,2,4,8,12,16 hours.
A list of 4 matrices 500 x 33
This method is used to calculate all the p-values relative to the variability of synthesis, processing and degradation rates. For object modeled with nascent RNA or when non-functional modeling was used, the variability is calculated using the confidence intervals. For objects modeled without nascent RNA, model selection is performed by comparing the likelihood of different (nested) models.
calculateRatePvals( object, modelSelection = c("aic", "llr", "hib"), preferPValue = TRUE, padj = TRUE, p_goodness_of_fit = 0.1, p_variability = rep(0.05, 3), limitModelComplexity = FALSE ) ## S4 method for signature 'INSPEcT' calculateRatePvals( object, modelSelection = c("aic", "llr", "hib"), preferPValue = TRUE, padj = TRUE, p_goodness_of_fit = 0.1, p_variability = rep(0.05, 3), limitModelComplexity = FALSE )
calculateRatePvals( object, modelSelection = c("aic", "llr", "hib"), preferPValue = TRUE, padj = TRUE, p_goodness_of_fit = 0.1, p_variability = rep(0.05, 3), limitModelComplexity = FALSE ) ## S4 method for signature 'INSPEcT' calculateRatePvals( object, modelSelection = c("aic", "llr", "hib"), preferPValue = TRUE, padj = TRUE, p_goodness_of_fit = 0.1, p_variability = rep(0.05, 3), limitModelComplexity = FALSE )
object |
An object of class INSPEcT or INSPEcT_model |
modelSelection |
'aic' compares nested models closest to the one with lowest AIC, 'llr' compares all nested models, 'hib' is a mix between the previous two. (default 'aic') |
preferPValue |
a logical, if TRUE (default) limit the search for best models among the ones with succeded the goodness of fit test. |
padj |
a logical, if TRUE (default) correct the p-values for multiple testing |
p_goodness_of_fit |
a numeric, the threshold for the goodness-of-fit test (default = .1) |
p_variability |
a numeric, a vector with the thresholds for the p-value of the variability test (one threshold for each rate, default = rep(.05, 3)) |
limitModelComplexity |
a logical that limits the complexity of the function used to describe dynamics to the length of the time-course (default = FALSE) |
ratePvals retrieve a single p-value for each rate and gene associated to its variability (null hypothesis = the rate is not changing between the conditions)
A matrix containing p-values calculated for each rate
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) # Set the chi-squared threshold at .2 for nascentInspObj10 object nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_goodness_of_fit=.2)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) # Set the chi-squared threshold at .2 for nascentInspObj10 object nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_goodness_of_fit=.2)
This method is used to retrieve the chi-squared test results for the models that have been selected to better represent the behavior of each gene.
chisqmodel(object, gc = NULL, tpts = NULL, ...) ## S4 method for signature 'INSPEcT' chisqmodel(object, gc = NULL, tpts = NULL, ...)
chisqmodel(object, gc = NULL, tpts = NULL, ...) ## S4 method for signature 'INSPEcT' chisqmodel(object, gc = NULL, tpts = NULL, ...)
object |
An object of class INSPEcT or INSPEcT_model |
gc |
Additional arguments for the generic |
tpts |
Additional arguments for the generic |
... |
Additional arguments for the generic |
A vector of chi-squared test results
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) chisqmodel(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) chisqmodel(nascentInspObj10)
This method is used to retrieve all the chi-squared test results for all models tested for all genes.
chisqtest(object, ...) ## S4 method for signature 'INSPEcT_model' chisqtest(object, ...) ## S4 method for signature 'INSPEcT' chisqtest(object, ...)
chisqtest(object, ...) ## S4 method for signature 'INSPEcT_model' chisqtest(object, ...) ## S4 method for signature 'INSPEcT' chisqtest(object, ...)
object |
An object of class INSPEcT or INSPEcT_model |
... |
Additional arguments for the generic |
A matrix of chi-squared test results for all the tested models
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) chisqtest(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) chisqtest(nascentInspObj10)
This method combines the information coming from different Objects of INSPEcT class. Requirements for two or more object to be combined together are:
they must be either modeled or either not modeled
they must have the same time points
they must have the same modeling parameters
## S4 method for signature 'INSPEcT,INSPEcT' combine(x, y, ...)
## S4 method for signature 'INSPEcT,INSPEcT' combine(x, y, ...)
x |
An object of class INSPEcT |
y |
An object of class INSPEcT |
... |
Additional objects of class INSPEcT |
In case the same gene is contained in more than one object that the user tries to combine, the information from one object will be used and a warning will be reported
An Object of class INSPEcT
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_2genes <- nascentInspObj10[1:2] nascentInspObj10_5genes <- nascentInspObj10[6:10] nascentInspObj10_7genes <- combine(nascentInspObj10_2genes, nascentInspObj10_5genes)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_2genes <- nascentInspObj10[1:2] nascentInspObj10_5genes <- nascentInspObj10[6:10] nascentInspObj10_7genes <- combine(nascentInspObj10_2genes, nascentInspObj10_5genes)
This method compares two object of class INSPEcT in order to identify differential usage of synthesis, processing or degradation rates in two different steady-state conditions. The two INSPEcT objects must have been profiled with replicates in order to provide a statistical significance to the differences between their rates.
compareSteady(inspectIds, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' compareSteady(inspectIds, BPPARAM = SerialParam())
compareSteady(inspectIds, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' compareSteady(inspectIds, BPPARAM = SerialParam())
inspectIds |
An object of calss INSPEcT with two conditions |
BPPARAM |
Configuration for BiocParallel parallelization. By default is set to SerialParam() |
An object of class INSPEcT_diffsteady which contains both the absolute
quantification of the rates as well as the comparison with the statistical significance
associated for each gene and rate. (See INSPEcT_diffsteady-class
)
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT( tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) }
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT( tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) }
This function compare exons and introns expression level matrices, from two up to an arbitrary number of samples, in order to identify genes which are oddly regluated, compared to an expected standard behaviour, from the post transcriptional point of view.
compareSteadyNoNascent( inspectIds, expressionThreshold = 0.25, log2FCThreshold = 2, trivialAngle = NaN, returnNormScores = FALSE, referenceCondition = "median" ) ## S4 method for signature 'INSPEcT_steadyNoNascent' compareSteadyNoNascent( inspectIds, expressionThreshold = 0.25, log2FCThreshold = 2, trivialAngle = NaN, returnNormScores = FALSE, referenceCondition = "median" )
compareSteadyNoNascent( inspectIds, expressionThreshold = 0.25, log2FCThreshold = 2, trivialAngle = NaN, returnNormScores = FALSE, referenceCondition = "median" ) ## S4 method for signature 'INSPEcT_steadyNoNascent' compareSteadyNoNascent( inspectIds, expressionThreshold = 0.25, log2FCThreshold = 2, trivialAngle = NaN, returnNormScores = FALSE, referenceCondition = "median" )
inspectIds |
An object of class INSPEcT_steadyNoNascent |
expressionThreshold |
A parameter which sets how many log2 fold changes of distance from the median behaviour are imputable to noise. |
log2FCThreshold |
A parameter which sets the log2 fold change distance from the median behaviour that is imputable to noise. |
trivialAngle |
A numeric between 0 and 90 to define the standard behavior, if NaN (default) it is computed internally from the data. |
returnNormScores |
A logical, if TRUE returned the deviations from the standard behavior normalized by the sd. |
referenceCondition |
The label of the condition to use as reference, if NaN (default) the medians are used. |
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj <- newINSPEcT(tpts=conditions,matureExpressions=matExp_DESeq2) matureInspObj<-compareSteadyNoNascent(inspectIds=matureInspObj ,expressionThreshold=0.25 ,log2FCThreshold=.5) regGenes <- PTreg(matureInspObj) head(regGenes) table(regGenes)
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj <- newINSPEcT(tpts=conditions,matureExpressions=matExp_DESeq2) matureInspObj<-compareSteadyNoNascent(inspectIds=matureInspObj ,expressionThreshold=0.25 ,log2FCThreshold=.5) regGenes <- PTreg(matureInspObj) head(regGenes) table(regGenes)
This function is used to compute the confidence intervals for a given set of modeled genes in the NoNascent scenario.
computeConfidenceIntervals(object, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' computeConfidenceIntervals(object, BPPARAM = SerialParam())
computeConfidenceIntervals(object, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' computeConfidenceIntervals(object, BPPARAM = SerialParam())
object |
An object of class INSPEcT_model |
BPPARAM |
Parallelization parameters for bplapply. By default SerialParam() |
An object of class INSPEcT.
This method is used to retrieve the convergence of the models that have been selected to better represent the behavior of each gene. 0 is converged, 1 not converged, 10 degenerated
convergence(object) ## S4 method for signature 'INSPEcT' convergence(object)
convergence(object) ## S4 method for signature 'INSPEcT' convergence(object)
object |
An object of class INSPEcT or INSPEcT_model |
A vector of numeric
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) convergence(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) convergence(nascentInspObj10)
This function plot the rates of a simulated dataset against the modeled ones and compute their correlations.
correlationPlot(object, object2, plot = TRUE) ## S4 method for signature 'INSPEcT_model,INSPEcT' correlationPlot(object, object2, plot = TRUE)
correlationPlot(object, object2, plot = TRUE) ## S4 method for signature 'INSPEcT_model,INSPEcT' correlationPlot(object, object2, plot = TRUE)
object |
An object of class INSPEcT_model with simulated rates. |
object2 |
An object of class INSPEcT. |
plot |
A logical indicating whether to draw or not the plot. (default=TRUE) |
An list with the correlation values.
A method to obtain the dimension of the object of class INSPEcT reported as a vector containing of the genes and the number of time points
## S4 method for signature 'INSPEcT' dim(x)
## S4 method for signature 'INSPEcT' dim(x)
x |
An object of class INSPEcT |
A numeric that indicates the number of genes within the object and the number of time points contained the object
Operators acting on INSPEcT, INSPEcT_model or INSPEcT_diffsteady objects to extract parts. INSPEcT_model objects can be subsetted only by gene. INSPEcT objects can be subsetted either by gene id or time point. In case of subsetting an INSPEcT object by time point, the model should be empty.
## S4 method for signature 'INSPEcT_model,ANY,ANY,ANY' x[i] ## S4 method for signature 'INSPEcT,ANY,ANY,ANY' x[i, j] ## S4 method for signature 'INSPEcT_diffsteady,ANY,ANY,ANY' x[i, j]
## S4 method for signature 'INSPEcT_model,ANY,ANY,ANY' x[i] ## S4 method for signature 'INSPEcT,ANY,ANY,ANY' x[i, j] ## S4 method for signature 'INSPEcT_diffsteady,ANY,ANY,ANY' x[i, j]
x |
An object of class INSPEcT or INSPEcT_model |
i |
A numeric, a vector of logicals or a vector of names indicating the features to be extracted |
j |
A numeric, a vector of logicals indicating the time points to be extracted |
An Object of class INSPEcT
removeModel
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_5genes <- nascentInspObj10[1:5] ## Not run: ## This will turn out into an error: nascentInspObj10_5genes_5tpts <- nascentInspObj10[1:5, 1:5] ## End(Not run) ## Before subsetting time points, the model should be removed: nascentInspObj10_5genes_5tpts <- removeModel(nascentInspObj10)[1:5, 1:5]
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_5genes <- nascentInspObj10[1:5] ## Not run: ## This will turn out into an error: nascentInspObj10_5genes_5tpts <- nascentInspObj10[1:5, 1:5] ## End(Not run) ## Before subsetting time points, the model should be removed: nascentInspObj10_5genes_5tpts <- removeModel(nascentInspObj10)[1:5, 1:5]
A method to visualize gene names associated with the object of class INSPEcT
## S4 method for signature 'INSPEcT' featureNames(object) ## S4 replacement method for signature 'INSPEcT' featureNames(object) <- value
## S4 method for signature 'INSPEcT' featureNames(object) ## S4 replacement method for signature 'INSPEcT' featureNames(object) <- value
object |
An object of class INSPEcT |
value |
A character that will replace the current feature names |
A character that contains gene names associated with the object of class INSPEcT
Contains two variables: "exWdths" and "intWdths" containing the lenght of the exons and introns, respectively, relative to the genes in "allcounts"
numeric vector of length 500
This method returns a factor that summarise the gene class (transcriptional regulatory mechanism) that
INSPEcT has assigned to each gene. The variability of each rate is indicated with a letter, 's' for
synthesis, 'p' for processing and 'd' for degradation. In case more than one rate is variable, the
letters associated to each variable rate are merged, for example 'sd' stands for a gene where synthesis
and degradation cotributed to transcriptional changes. 'no-reg' is associated to genes with no
change in transcription. The classification depends on the thresholds of the goodness-of-fit and
rate variability tests that can be changed via the method calculateRatePvals
.
geneClass(object, ...) ## S4 method for signature 'INSPEcT' geneClass(object, ...) ## S4 method for signature 'INSPEcT_model' geneClass(object, ...) ## S4 method for signature 'INSPEcT_diffsteady' geneClass(object, ...)
geneClass(object, ...) ## S4 method for signature 'INSPEcT' geneClass(object, ...) ## S4 method for signature 'INSPEcT_model' geneClass(object, ...) ## S4 method for signature 'INSPEcT_diffsteady' geneClass(object, ...)
object |
An object of class INSPEcT or INSPEcT_model |
... |
specify the threshold for rate variability 'bTsh' in case of 'INSPEcT_diffsteady' objects (default = .1) |
A character containing the regulatory class for each gene
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) geneClass(nascentInspObj10) # see the classification with another threshold for rate variability nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_variability=rep(1,3)) geneClass(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) geneClass(nascentInspObj10) # see the classification with another threshold for rate variability nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_variability=rep(1,3)) geneClass(nascentInspObj10)
A method to see as an heatmap the logRatios of synthesis, degradation and processing rates and pre-RNA and total RNA concentration of a population of genes, either at the level of etimated or modeled rates.
inHeatmap( object, type = "pre-model", breaks = seq(-1, 1, length.out = 51), palette = colorRampPalette(c("green", "black", "firebrick3")), plot_matureRNA = FALSE, absoluteExpression = TRUE, show_rowLabels = TRUE, clustering = TRUE, clustIdx = 3:5 ) ## S4 method for signature 'INSPEcT' inHeatmap( object, type = "pre-model", breaks = seq(-1, 1, length.out = 51), palette = colorRampPalette(c("green", "black", "firebrick3")), plot_matureRNA = FALSE, absoluteExpression = TRUE, show_rowLabels = TRUE, clustering = TRUE, clustIdx = 3:5 )
inHeatmap( object, type = "pre-model", breaks = seq(-1, 1, length.out = 51), palette = colorRampPalette(c("green", "black", "firebrick3")), plot_matureRNA = FALSE, absoluteExpression = TRUE, show_rowLabels = TRUE, clustering = TRUE, clustIdx = 3:5 ) ## S4 method for signature 'INSPEcT' inHeatmap( object, type = "pre-model", breaks = seq(-1, 1, length.out = 51), palette = colorRampPalette(c("green", "black", "firebrick3")), plot_matureRNA = FALSE, absoluteExpression = TRUE, show_rowLabels = TRUE, clustering = TRUE, clustIdx = 3:5 )
object |
An object of class INSPEcT |
type |
Eiher "pre-model" or "model" to switch between pre-modeled or modeled features |
breaks |
A vector of breaks for the heatmap |
palette |
A color generating function, output of colorRampPalette |
plot_matureRNA |
A logical. If set to TRUE, matrue-RNA is displayed instead of total-RNA (default: FALSE) |
absoluteExpression |
A logical. If set to FALSE, the plot representing the intensity of expression is omitted. (default=TRUE) |
show_rowLabels |
A logical defining whether rownames are reported or not. (default=TRUE) |
clustering |
A logical. If set to FALSE, it displys genes the order they are, with no clustering (default: TRUE) |
clustIdx |
A numeric. Indicates which of the features are used for the clustering. 0=absoluteExpression; 1=total-RNA/mature-RNA; 2=preMRNA; 3=synthesis; 4=degradation; 5=processing (default=3:5, meaning that synthesis, degradation and processing are used for the clustering) |
A list of matrices containing the logRatios for total RNA levels, pre-RNA levels, synthesis rates, degradation rates and processing rates. Matrices are ordered according to the clustering.
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) inHeatmap(nascentInspObj10, 'pre-model') inHeatmap(nascentInspObj10, 'model')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) inHeatmap(nascentInspObj10, 'pre-model') inHeatmap(nascentInspObj10, 'model')
INSPEcT (INference of Synthesis, Processing and dEgradation rates from Transcriptomic data), is a package that analyse RNA-seq data in order to evaluate synthesis, processing and degradation rates and asses via modeling the rates that determines changes in RNA levels.
To see how the typical workflow of INSPEcT works, type: vignette('INSPEcT')
INSPEcT implements two main classes (INSPEcT
and
INSPEcT_diffsteady
) and their corresponding methods.
To have a detailed description of how the two classes are structured
and which methods apply on, type:
?'INSPEcT-class'
?'INSPEcT_diffsteady-class'
To obtain the citation, type: citation('INSPEcT')
INSPEcT_diffsteady is a class able to store the results of the comparisons between two steady states.
An object of class INSPEcT_diffsteady is created with the method "compareSteady" applied on
two "INSPEcT" objects (see compareSteady
).
synthesis(object) processing(object) degradation(object) ## S4 method for signature 'INSPEcT_diffsteady' show(object) ## S4 method for signature 'INSPEcT_diffsteady' synthesis(object) ## S4 method for signature 'INSPEcT_diffsteady' processing(object) ## S4 method for signature 'INSPEcT_diffsteady' degradation(object) ## S4 method for signature 'INSPEcT_diffsteady' featureNames(object)
synthesis(object) processing(object) degradation(object) ## S4 method for signature 'INSPEcT_diffsteady' show(object) ## S4 method for signature 'INSPEcT_diffsteady' synthesis(object) ## S4 method for signature 'INSPEcT_diffsteady' processing(object) ## S4 method for signature 'INSPEcT_diffsteady' degradation(object) ## S4 method for signature 'INSPEcT_diffsteady' featureNames(object)
object |
An object of class INSPEcT_model |
Methods associated to the class INSPEcT_diffsteady are:
synthesis: Accessor to the synthesis rates and their comparisons.
degradation: Accessor to the degradation rates and their comparisons.
processing: Accessor to the processing rates and their comparisons.
plotMA: visualization fuction for rates comparisons, see plotMA
Method show for objects of class INSPEcT_model returns the number of the genes that have been modeled
synthesis
A data.frame which contains both input data and comparisons results regarding synthesis rates
degradation
A data.frame which contains both input data and comparisons results regarding degradation rates
processing
A data.frame which contains both input data and comparisons results regarding processing rates
modeling_res
A data.frame which contains modeling results
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT(tpts=conditions,labeling_time=tL ,nascentExpressions=nasExp_DESeq2,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) head(synthesis(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { head(processing(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { head(degradation(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { featureNames(diffrates) }
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT(tpts=conditions,labeling_time=tL ,nascentExpressions=nasExp_DESeq2,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) head(synthesis(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { head(processing(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { head(degradation(diffrates)) } if( Sys.info()["sysname"] != "Windows" ) { featureNames(diffrates) }
INSPEcT_model is a class able to store all the results of the modeling of synthesis, processing and degradation rates made via the method modelRates
(slot ratesSpecs).
It also stores the criteria (slot parameter) to choose between the many models tested for each gene the one that better describes the data and the results.
The slot simple is a flag that distinguish wheter the model contains the information of the introns or not.
In case not, the flag simple is set to TRUE.
Also the method makeSimModel
of class INSPEcT-class
creates an object of class INSPEcT_model.
This object will be used by makeSimDataset
to generate a complete simulated data-set, whose classification performance can be tested.
## S4 method for signature 'INSPEcT_model' show(object)
## S4 method for signature 'INSPEcT_model' show(object)
object |
An object of class INSPEcT_model |
Methods that apply to INSPEcT_model class are
[
AIC
chisqtest
correlationPlot
geneClass
logLik
makeModelRates
makeSimDataset
modelSelection
rocCurve
rocThresholds
show
Method show for objects of class INSPEcT_model returns the number of the genes that have been modeled
params
A list that defines thresholds and how to perform log likelihood ratio tests
ratesSpecs
A list containing the modeling output
simple
A logical that indicates whether the mode of INSPEcT is simple (no pre-mRNA and degradation rates) or not.
INSPEcT_steadyNoNascent is a class able to store data and arguments that are necessary to
make the analysis concerning premature and mature expressions in different samples. In
particular, the ratio between mature and premature can be calculated, which reflects
the ratio between the rates of processing and degradation in individaul genes
(see PTratio
), and the analysis of post-transcriptionally regualted genes
can be run to identify genes that in specfic samples show a trand which cannot be attributed
to transcriptional regulation alone (see PTreg
).
## S4 method for signature 'INSPEcT_steadyNoNascent,ANY,ANY,ANY' x[i, j] ## S4 method for signature 'INSPEcT_steadyNoNascent' show(object)
## S4 method for signature 'INSPEcT_steadyNoNascent,ANY,ANY,ANY' x[i, j] ## S4 method for signature 'INSPEcT_steadyNoNascent' show(object)
x |
An object of class INSPEcT_steadyNoNascent |
i |
A numeric, a vector of logicals indicating the rows to be extracted |
j |
A numeric, a vector of logicals indicating the columns to be extracted |
object |
An object of class INSPEcT_steadyNoNascent |
Method show for objects of class INSPEcT_steadyNoNascent
sampleNames
Vector with the names of the samples (columns of the dataset)
geneNames
Vector with the names of the genes (rows of the dataset)
premature
Matrix containing the expressions of the premature RNAs (row=genes, columns=samples)
mature
Matrix containing the expressions of the emature RNAs (row=genes, columns=samples)
prematureVar
Matrix containing the expressions variances of the premature RNAs (row=genes, columns=samples)
matureVar
Matrix containing the expressions variances of the emature RNAs (row=genes, columns=samples)
trivialAngle
Numeric that indicates the angle (slope) of the linear model between mature and premature expressions
log2FCThreshold
Numeric that describes the threshold of the variation to be considered significant
expressionThreshold
Numeric that describes the threshold of the expression to consider the gene expressed
referenceCondition
A sample identifier that set the reference for the post-transcriptional regulation analysis, if NULL the median of all samples is used
ptreg
Matrix containing the post-transcriptioanl regulation state of each gene in the different samples (row=genes, columns=samples)
INSPEcT is a class able to store all the estimated rates and concentrations (slot ratesFirstGuess),
the modeled rates and concentrations (slot modelRates) and the model themselves (slot model).
Within the class INSPEcT other information regarding the experimental design are stored, such as the time
points where experiments were performed (slot tpts) and, if provided, the nascent RNA collecting time (slot tL)
and the normalization scale fators used for nascent (labeledSF) RNA-seq libraries.
A list of parameters that will be used during the
modeling process is stored within the slot params and can be accessed by modelingParams
.
A new instance of the class INSPEcT can be generated by the constructor function newINSPEcT
.
## S4 method for signature 'INSPEcT' show(object)
## S4 method for signature 'INSPEcT' show(object)
object |
An object of class INSPEcT |
Methods that apply to INSPEcT class are
AIC
[
calculateDelta
calculateRatePvals
calculateTau
chisqmodel
chisqtest
combine
compareSteady
compareSteadyNoNascent
computeConfidenceIntervals
correlationPlot
dim
featureNames
geneClass
inHeatmap
labeledSF
logLik
makeModelRates
makeOscillatorySimModel
makeSimModel
modelRates
modelRatesNF
modelSelection
modelingParams
nGenes
nTpts
plotGene
processingDelay
ratePvals
ratesFirstGuess
ratesFirstGuessVar
removeModel
rocCurve
rocThresholds
setConfidenceIntervals
show
split
tpts
viewConfidenceIntervals
viewModelRates
Method show for objects of class INSPEcT displays the main features of the slots ratesFirstGuess, model and modelRates
params
A list of parameters of the modeling part
ratesFirstGuess
An object of class ExpressionSet that contains all the rates and concentrations guessed from the first part of INSPEcT analysis (before modeling)
ratesFirstGuessVar
An object of class ExpressionSet that contains the variances related to rates and concentrations guessed from the first part of INSPEcT analysis (before modeling)
confidenceIntervals
An object of class ExpressionSet that contains the confidence intervals.
model
An object of class INSPEcT_model that contains the output of the mdoeling.
modelRates
An object of class ExpressionSet that contains all modeled the rates and concentrations.
ratePvals
A matrix containing the p-value relative to the variability of synthesis, processing and degradation for each gene.
tpts
A numeric vector of the time-points.
labeledSF
A numeric vector of the scaling factor used for inter time-point normalization of Nascent-seq libraries.
tL
A numeric containing the length of the Nascent pulse.
NoNascent
A logical indicating if the nascent RNA was included into the analysis.
NF
A logical indicating if the modeling approach is Non-Functional
degDuringPulse
A logical indicating if degradation of RNA during the 4sU pulse was considered.
version
A character indicating the version of INSPEcT that created the object
Two shiny apps are encoded into the package inspect: - runProcessingRateDelay: plots single genes as well as genome wide plots associated to the processing induced delay loading the data from an INSPEcT object. - runINSPEcTGUI: is a way to visualize and interact with the RNAdynamics at the level of a single gene, either loading the data from ad INSPEcT object or from scratch.
runProcessingRateDelay() runINSPEcTGUI()
runProcessingRateDelay() runINSPEcTGUI()
Function to run the whole INSPEcT differential rate analysis procedure with a single line. The function save the output analysis to file that can be later loaded in the R environment or in the INSPEcT-GUI.
inspectFromBAM( txdb, annotation_table, labeling_time = NULL, strandSpecific = 0, isPairedEnd = FALSE, estimateRatesWith = "der", useSigmoidFun = TRUE, file = NULL )
inspectFromBAM( txdb, annotation_table, labeling_time = NULL, strandSpecific = 0, isPairedEnd = FALSE, estimateRatesWith = "der", useSigmoidFun = TRUE, file = NULL )
txdb |
A TranscriptDB object for the selected organism |
annotation_table |
Paths and experimental design associated to bam files. They could be provided directly as a 'data.frame', or as a path to the file containing the information. Possible file formats are csv' (comma-separated-values), 'tsv' (comma-separated-values), or 'xls' (Excel). In case 'annotation_table' has 2 colums named 'condition' and 'total', INSPEcT- analysis is run. In case 'annotation_table' has 3 colums named 'condition', 'total' and 'nascent', INSPEcT+ analysis is run. 'condition' is a colums indicating the experimental condition, a character vector (containing, for example, 'WT' or 'KD') in case of steady-state experiments, or numerical values indicating the time from the unperturbed condition in case of time-course analysis. 'total' and 'nascent' contains the path to totalRNA and nascentRNA BAM files, respectively. |
labeling_time |
A numeric indicating the time of labeling exposure to the modified nucleotide. To be indicated only in case of INSPEcT+ analysis. |
strandSpecific |
A numeric indicating the strandness of the BAM files, 0 for non strand-specific, 1 for stranded, 2 for reversely-stranded. 0 by default. |
isPairedEnd |
A logical indicating if paired-end sequencing have been performed. FALSE by default. |
estimateRatesWith |
Either "int" or "der". With "int" the degradation and processing rates are estimated integrating the system between one time point and the following. With "der" degradation and processing rates are estimated using the derivative of total and pre mRNA. (default is "der") |
useSigmoidFun |
A logical, whether to choose between sigmoid and impulse function to fit rates and concentrations. In case not, always impulse function is used. (default is TRUE) |
file |
A character indicating where the output of the analysis will be stored. If not provided the file name will be created automaticcally and saved on the current folder. |
Function to run the whole INSPEcT differential rate analysis procedure with a single line. The function save the output analysis to file that can be later loaded in the R environment or in the INSPEcT-GUI.
inspectFromPCR( totalRNA_table, nascentRNA_table = NULL, labeling_time = NULL, estimateRatesWith = "der", useSigmoidFun = TRUE, file = NULL )
inspectFromPCR( totalRNA_table, nascentRNA_table = NULL, labeling_time = NULL, estimateRatesWith = "der", useSigmoidFun = TRUE, file = NULL )
totalRNA_table |
Exonic quantification, intronic quantification and experimental design associated to totalRNA of a single gene quantified by PCR. They could be provided directly as a 'data.frame', or as a path to the file containing the information. Possible file formats are csv' (comma-separated-values), 'tsv' (comma-separated-values), or 'xls' (Excel). 'totalRNA_table' must have 3 colums named 'condition', 'total_exonic' and 'total_intronic'. 'condition' is a column indicating the experimental condition, a character vector (containing, for example, 'WT' or 'KD') in case of steady-state experiments, or numerical values indicating the time from the unperturbed condition in case of time-course analysis. 'total_exonic' and 'total_intronic' contains abundance of gene measured in its exonic and intronic regions, respectively, in the total RNA fraction. |
nascentRNA_table |
similar to 'totalRNA_table' but referred to nascent RNA fraction. In this case, colums names must be 'condition', 'nascent_exonic' and 'nascent_intronic'. In case this infromation is not provided, INSPEcT- analysis is run. If otherwise this information is present, INSPEcT+ analysis is run. |
labeling_time |
A numeric indicating the time of labeling exposure to the modified nucleotide. To be indicated only in case of INSPEcT+ analysis. |
estimateRatesWith |
Either "int" or "der". With "int" the degradation and processing rates are estimated integrating the system between one time point and the following. With "der" degradation and processing rates are estimated using the derivative of total and pre mRNA. (default is "der") |
useSigmoidFun |
A logical, whether to choose between sigmoid and impulse function to fit rates and concentrations. In case not, always impulse function is used. (default is TRUE) |
file |
A character indicating where the output of the analysis will be stored. If not provided the file name will be created automaticcally and saved on the current folder. |
if( Sys.info()["sysname"] != "Windows" ) { totalAnnTabPCR <- system.file(package = 'INSPEcT', 'totalAnnTabPCR.csv') nascentAnnTabPCR <- system.file(package = 'INSPEcT', 'nascentAnnTabPCR.csv') inspectFromPCR(totalAnnTabPCR, nascentAnnTabPCR, labeling_time=1/6) }
if( Sys.info()["sysname"] != "Windows" ) { totalAnnTabPCR <- system.file(package = 'INSPEcT', 'totalAnnTabPCR.csv') nascentAnnTabPCR <- system.file(package = 'INSPEcT', 'nascentAnnTabPCR.csv') inspectFromPCR(totalAnnTabPCR, nascentAnnTabPCR, labeling_time=1/6) }
Accessor to obtain the labeledSF slot associated with the object of class INSPEcT
labeledSF(object) ## S4 method for signature 'INSPEcT' labeledSF(object)
labeledSF(object) ## S4 method for signature 'INSPEcT' labeledSF(object)
object |
An object of class INSPEcT |
A numeric that indicates the scaling factors applied between time points of the data coming from Nascent-seq library (applies directly to synthesis rates and indirectly to degradation rates)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) labeledSF(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) labeledSF(nascentInspObj10)
Contains two variables: "nascentLS" and "totalLS" containing the sequencing depth of nascent and total libraries respectively, relative to the experiments in "allcounts"
numeric vector of length 33
This method is used to retrieve all the log likelihood ratio test results for all pairs tested for all genes.
logLik(object, ...) ## S4 method for signature 'INSPEcT_model' logLik(object, ...) ## S4 method for signature 'INSPEcT' logLik(object, ...)
logLik(object, ...) ## S4 method for signature 'INSPEcT_model' logLik(object, ...) ## S4 method for signature 'INSPEcT' logLik(object, ...)
object |
An object of class INSPEcT or INSPEcT_model |
... |
Additional arguments for the generic |
A matrix of log likelihood test results for all the tested model comparisons
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) logLik(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) logLik(nascentInspObj10)
This function is used to evaluate rates and concentrations after modeling of the
rates has been run with modelRates
. The modeled rates are in functional
form and can be evaluated at any time points.
This method can be used to regenerate the rates assiciated to the modeling, in case some testing parameters has changed.
makeModelRates(object, ...) ## S4 method for signature 'INSPEcT_model' makeModelRates(object, ...) ## S4 method for signature 'INSPEcT' makeModelRates(object, ...)
makeModelRates(object, ...) ## S4 method for signature 'INSPEcT_model' makeModelRates(object, ...) ## S4 method for signature 'INSPEcT' makeModelRates(object, ...)
object |
An object of class INSPEcT_model |
... |
additional arguments tpts : A vector of time points where rates and concentrations have to be evaluated |
An object of class ExpressionSet containing the modeled rates and concentrations
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewModelRates(nascentInspObj10, 'degradation') ## force every degradation rate to be accepted as variable (makeModelRates is called internally) nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_variability = c(.05,.05,1)) viewModelRates(nascentInspObj10, 'degradation')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewModelRates(nascentInspObj10, 'degradation') ## force every degradation rate to be accepted as variable (makeModelRates is called internally) nascentInspObj10 <- calculateRatePvals(nascentInspObj10, p_variability = c(.05,.05,1)) viewModelRates(nascentInspObj10, 'degradation')
This method allow the creation of synthesis, degradation and processing rates that generate
an oscillatory expression with a period of 24 hours. Two modes are available: one where
oscillations arise just by oscillations in the synthesis of the genes (oscillatoryk3=FALSE, default) and
another one where both synthesis and degradation rates oscillates (oscillatoryk3=TRUE). In this latter case,
the oscillations of the two rates can be coupled by a cetrain delay (parametrer k3delay).
After the creation of the synthetic rates, a dataset with noise and contamination added can be made by makeSimDataset
.
makeOscillatorySimModel( object, nGenes, oscillatoryk3 = FALSE, k3delay = NULL, na.rm = TRUE, seed = NULL ) ## S4 method for signature 'INSPEcT' makeOscillatorySimModel( object, nGenes, oscillatoryk3 = FALSE, k3delay = NULL, na.rm = TRUE, seed = NULL )
makeOscillatorySimModel( object, nGenes, oscillatoryk3 = FALSE, k3delay = NULL, na.rm = TRUE, seed = NULL ) ## S4 method for signature 'INSPEcT' makeOscillatorySimModel( object, nGenes, oscillatoryk3 = FALSE, k3delay = NULL, na.rm = TRUE, seed = NULL )
object |
An object of class INSPEcT |
nGenes |
A numeric with the number of synthtic genes to be created |
oscillatoryk3 |
A logical that enables also degradation rate to oscillate |
k3delay |
A numeric that set the delay between synthesis and degradation oscillations. When NULL, no coupling between the two oscillations is set. |
na.rm |
A logical that set whether missing values in the real dataset should be removed |
seed |
A numeric to obtain reproducible results |
An object of class INSPEcT_model with synthetic rates
nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeOscillatorySimModel(nascentInspObj, 1000, seed=1) table(geneClass(simRates))
nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeOscillatorySimModel(nascentInspObj, 1000, seed=1) table(geneClass(simRates))
This method generates rates and concentrations where noise is added according to the desired number of
replicates that the user set as an arguments from the INSPEcT_model object that has been created by the
method of the class INSPEcT makeSimModel
. Rates and concentrations can be generated at the
time-points of interest within the original time window. This method generates an INSPEcT object that can
be modeled and the performance of the modeling can be tested directly aginst the INSPEcT_model object
created by makeSimModel
.
makeSimDataset( object, tpts, nRep, NoNascent = FALSE, seed = NULL, b = 0.3, tL = 1/6, noise_sd = 4 ) ## S4 method for signature 'INSPEcT_model' makeSimDataset( object, tpts, nRep, NoNascent = FALSE, seed = NULL, b = 0.3, tL = 1/6, noise_sd = 4 )
makeSimDataset( object, tpts, nRep, NoNascent = FALSE, seed = NULL, b = 0.3, tL = 1/6, noise_sd = 4 ) ## S4 method for signature 'INSPEcT_model' makeSimDataset( object, tpts, nRep, NoNascent = FALSE, seed = NULL, b = 0.3, tL = 1/6, noise_sd = 4 )
object |
An object of class INSPEcT_model, usually the output of |
tpts |
A numeric vector of time points where rates and concentrations have to be evaluated |
nRep |
Number of replicates to simulate |
NoNascent |
A logical which, if true, makes the output of the method suitable for an analysis without Nascent. (default=FALSE) |
seed |
A numeric to obtain reproducible results. When NULL (default) no seed is set. |
b |
A numeric which represents the probability of contamination of the unlabeled sample due to the labled one |
tL |
A numeric which represents the labeling time for an ideal nascent RNA profiling, it is required for the contamination analysis. (default=1/6) |
noise_sd |
A numeric which represents the noise standard deviation. (default=4) |
An object of the class ExpressionSet containing rates and concentrations
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) tpts <- tpts(nascentInspObj) nascentSim2replicates <- makeSimDataset(object=simRates,tpts=tpts,nRep=3,NoNascent=FALSE,seed=1) }
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) tpts <- tpts(nascentInspObj) nascentSim2replicates <- makeSimDataset(object=simRates,tpts=tpts,nRep=3,NoNascent=FALSE,seed=1) }
This method allow the creation of synthesis, degradation and processing rates for a certain number of genes.
The rates are created according to the distributions of the real data-set which is given as an input of the
method. Different proportions of constant varying rates can be set and a new vector of time points can be
provided. This method has to be used before the makeSimDataset
method.
makeSimModel( object, nGenes, newTpts = NULL, probs = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), na.rm = TRUE, seed = NULL ) ## S4 method for signature 'INSPEcT' makeSimModel( object, nGenes, probs = rbind(synthesis = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), processing = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), degradation = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2)), na.rm = TRUE, seed = NULL )
makeSimModel( object, nGenes, newTpts = NULL, probs = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), na.rm = TRUE, seed = NULL ) ## S4 method for signature 'INSPEcT' makeSimModel( object, nGenes, probs = rbind(synthesis = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), processing = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2), degradation = c(constant = 0.5, sigmoid = 0.3, impulse = 0.2)), na.rm = TRUE, seed = NULL )
object |
An object of class INSPEcT |
nGenes |
A numeric with the number of synthtic genes to be created |
newTpts |
A numeric verctor with time points of the synthtic dataset, if NULL the time points of the real dataset will be used |
probs |
A numeric matrix wich describes the probability of each rate (rows) to be constant, shaped like a sigmoid or like an impulse model (columns) |
na.rm |
A logical that set whether missing values in the real dataset should be removed |
seed |
A numeric to obtain reproducible results |
The method makeSimModel
generates an object of class INSPEcT_model that stores the parametric functions to genrate clean rates of a time-course. To any of the rates also a noise variance is associate but not used yet. In a typical workflow the output of makeSimModel
is the input of the method makeSimDataset
, that build the noisy rates and concentrations, given a specified number of replicates.
An object of class INSPEcT_model with synthetic rates
nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) table(geneClass(simRates))
nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) table(geneClass(simRates))
Extract mature RNA expressions
mature(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' mature(object)
mature(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' mature(object)
object |
An object of class INSPEcT_steadyNoNascent |
A matrix containing mature RNA expressions
Extract mature RNA expressions variances
matureVar(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' matureVar(object)
matureVar(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' matureVar(object)
object |
An object of class INSPEcT_steadyNoNascent |
A matrix containing mature RNA expressions variances
A method to get the parameters used for modeling rates and
concentrations by the method modelRates
modelingParams(object) ## S4 method for signature 'INSPEcT' modelingParams(object)
modelingParams(object) ## S4 method for signature 'INSPEcT' modelingParams(object)
object |
An object of class INSPEcT |
List of parameters and their values
estimateRatesWith Either "int" or "der". With "int" the degradation and processing rates are estimated integrating the system between one time point and the following. With "der" degradation and processing rates are estimated using the derivative of total and pre mRNA.
useSigmoidFun A logical, whether to choose between sigmoid and impulse function to fit rates and concentrations. In case not, always impulse function is used.
testOnSmooth A logical, wheter models should be tested on smoothed pre-mRNA, total mRNA and eventually synthesis rates or not.
nInit number of optimization to find the best functional representation of each rate
nIter number of max iteration during optimization
Dmin lower bondary for degradation rates in the NoNascent mode
Dmax upper bondary for degradation rates in the NoNascent mode
seed A numeric, indicatindg the seed set for reproducible results.
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) modelingParams(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) modelingParams(nascentInspObj10)
Launch the modeling process with parameters set with modelingParams
This method model the synthesis, degradation and processing rates after their estimation by the constructor function
newINSPEcT
. Estimated rates are not guaranteed to optimally describes provided input data yet.
To this purpose, modeled rates can be generated and genes can be assigned to a transcriptional regulatory mechanism.
Modeled rates can be accessed via the method viewModelRates
and gene classification according
to the regulatory mechanism can be accessed by geneClass
. The modeling options used for the
modeling can be later accessed by the user via modelingParams
. After modeling, model selection is run
by the method calculateRatePvals
with default parameters.
modelRates( object, estimateRatesWith = c("der", "int"), useSigmoidFun = TRUE, nInit = 10, nIter = 300, Dmin = 1e-06, Dmax = 10, seed = NULL, BPPARAM = SerialParam() ) ## S4 method for signature 'INSPEcT' modelRates( object, estimateRatesWith = c("der", "int"), useSigmoidFun = TRUE, nInit = 10, nIter = 300, Dmin = 1e-06, Dmax = 10, seed = NULL, BPPARAM = SerialParam() )
modelRates( object, estimateRatesWith = c("der", "int"), useSigmoidFun = TRUE, nInit = 10, nIter = 300, Dmin = 1e-06, Dmax = 10, seed = NULL, BPPARAM = SerialParam() ) ## S4 method for signature 'INSPEcT' modelRates( object, estimateRatesWith = c("der", "int"), useSigmoidFun = TRUE, nInit = 10, nIter = 300, Dmin = 1e-06, Dmax = 10, seed = NULL, BPPARAM = SerialParam() )
object |
An object of class INSPEcT |
estimateRatesWith |
Either "int" or "der". With "int" the degradation and processing rates are estimated integrating the system between one time point and the following. With "der" degradation and processing rates are estimated using the derivative of total and pre mRNA. (default is "der") |
useSigmoidFun |
A logical, whether to choose between sigmoid and impulse function to fit rates and concentrations. In case not, always impulse function is used. (default is TRUE) |
nInit |
number of optimization to find the best functional representation of each rate (by default 10) |
nIter |
number of max iteration during optimization (default is 300) |
Dmin |
lower bondary for degradation rates in the NoNascent mode (default 1e-06) |
Dmax |
upper bondary for degradation rates in the NoNascent mode (default 10) |
seed |
A numeric, indicatindg the seed to be set for reproducible results. If NULL it is randomly selected (default NULL) |
BPPARAM |
Parallelization parameters for bplapply. By default SerialParam() |
An object of class INSPEcT with modeled rates
viewModelRates
, calculateRatePvals
, geneClass
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ## models removal nascentInspObjThreeGenes <- removeModel(nascentInspObj10[1:3]) nascentInspObjThreeGenes <- modelRates(nascentInspObjThreeGenes, seed=1, BPPARAM=SerialParam()) ## view modeled synthesis rates viewModelRates(nascentInspObjThreeGenes, 'synthesis') ## view gene classes geneClass(nascentInspObjThreeGenes) }
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ## models removal nascentInspObjThreeGenes <- removeModel(nascentInspObj10[1:3]) nascentInspObjThreeGenes <- modelRates(nascentInspObjThreeGenes, seed=1, BPPARAM=SerialParam()) ## view modeled synthesis rates viewModelRates(nascentInspObjThreeGenes, 'synthesis') ## view gene classes geneClass(nascentInspObjThreeGenes) }
This method compute confidence intervals for the rates of synthesis, degradation and processing estimated by
newINSPEcT
that will be used to estimate the variability of each rate in ratePvals
method.
modelRatesNF(object, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' modelRatesNF(object, BPPARAM = SerialParam())
modelRatesNF(object, BPPARAM = SerialParam()) ## S4 method for signature 'INSPEcT' modelRatesNF(object, BPPARAM = SerialParam())
object |
An object of class INSPEcT |
BPPARAM |
Parallelization parameters for bplapply. By default SerialParam() |
An object of class INSPEcT with modeled rates
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ## models removal nascentInspObjThreeGenes <- removeModel(nascentInspObj10[1:3]) nascentInspObjThreeGenes <- modelRatesNF(nascentInspObjThreeGenes, BPPARAM=SerialParam()) ## view modeled synthesis rates viewModelRates(nascentInspObjThreeGenes, 'synthesis') ## view gene classes geneClass(nascentInspObjThreeGenes) }
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ## models removal nascentInspObjThreeGenes <- removeModel(nascentInspObj10[1:3]) nascentInspObjThreeGenes <- modelRatesNF(nascentInspObjThreeGenes, BPPARAM=SerialParam()) ## view modeled synthesis rates viewModelRates(nascentInspObjThreeGenes, 'synthesis') ## view gene classes geneClass(nascentInspObjThreeGenes) }
Method to visualize the criteria used to assess variability of rates.
modelSelection(object) ## S4 method for signature 'INSPEcT' modelSelection(object) ## S4 method for signature 'INSPEcT_model' modelSelection(object)
modelSelection(object) ## S4 method for signature 'INSPEcT' modelSelection(object) ## S4 method for signature 'INSPEcT_model' modelSelection(object)
object |
An object of class INSPEcT or INSPEcT_model |
modelSelection 'aic' compares nested models closest to the one with lowest AIC, 'llr' compares all nested models, 'hib' is a mix between the previous two. (default 'aic')
preferPValue a logical, if TRUE (default) limit the search for best models among the ones with succeded the goodness of fit test.
padj a logical, if TRUE (default) correct the p-values for multiple testing
goodness_of_fit a numeric, the threshold for the goodness-of-fit test (default = .1)
variability a numeric, a vector with the thresholds for the variability test (one threshold for each rate, default = c('s'=05, 'p'=.05, 'd'=.05))
limitModelComplexity a logical that limits the complexity of the function used to describe dynamics to the length of the time-course (default = FALSE)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) modelSelection(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) modelSelection(nascentInspObj10)
The function newINSPEcT creates a new instance of the class INSPEcT provided the experimental time points, expression data (like RPKMs) of mature and eventually nascent RNA. For the nascent analysis, it is also requires a collecting time and the scaling factor to normalize the nascent RNA-seq libraries. This latter parameter can also be calculated by the function itself if both exonic and intronic expression data are provided; otherwise it must be given as an input and it is essential to guarantee the robustness of the analysis.
newINSPEcT( tpts, labeling_time = NULL, nascentExpressions = NULL, matureExpressions, preexisting = FALSE, BPPARAM = SerialParam(), labeledSF = NULL, simulatedData = FALSE, degDuringPulse = FALSE, Dmin = 1e-06, Dmax = 10, genesFilter = TRUE, genesFilterThreshold = 2/3, imputeNAs = TRUE )
newINSPEcT( tpts, labeling_time = NULL, nascentExpressions = NULL, matureExpressions, preexisting = FALSE, BPPARAM = SerialParam(), labeledSF = NULL, simulatedData = FALSE, degDuringPulse = FALSE, Dmin = 1e-06, Dmax = 10, genesFilter = TRUE, genesFilterThreshold = 2/3, imputeNAs = TRUE )
tpts |
A vector of time points, one for each sample |
labeling_time |
A number, lenght of the Nascent pulse |
nascentExpressions |
A list which contains exons and introns expression matrices and variances for the nascent RNA |
matureExpressions |
A list which contains exons and introns expression matrices and variances for the mature RNA |
preexisting |
A logical, indicating if the mature expression refers to the pre-exising (unlabeled) population. Not implemented yet for the "degDuringPulse" mode. |
BPPARAM |
Configuration for BiocParallel parallelization. By default is set to SerialParam() |
labeledSF |
A vector storing user defined normalization scale over Nascent RNA exons and introns quantifications |
simulatedData |
A logical, set to TRUE in case the analysis is on simulated data |
degDuringPulse |
A logical, set to TRUE in case of a long labelling time. Also degradation of newly synthesized transcripts will be taken into account |
Dmin |
A numerical, it is the lower bound of the degradation rate domain for the prior optimization |
Dmax |
A numerical, it is the upper bound of the degradation rate domain for the prior optimization |
genesFilter |
A logical, if TRUE, filters out genes which have no signal in at least a given fraction (2/3 by default) of the observations |
genesFilterThreshold |
A number, threshold to use for genes filtering (2/3 by default) |
imputeNAs |
A logical, if TRUE the rates first guess which are not finite are imputed from the neighbours. |
An object of class INSPEcT with a first estimation of the rates which can be accessed by the method ratesFirstGuess
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') matureCounts<-allcounts$mature tpts <- c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16) expDes<-rep(tpts,3) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj<-newINSPEcT(tpts=tpts ,labeling_time=NULL ,nascentExpressions=NULL ,matureExpressions=matExp_DESeq2)
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') matureCounts<-allcounts$mature tpts <- c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16) expDes<-rep(tpts,3) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj<-newINSPEcT(tpts=tpts ,labeling_time=NULL ,nascentExpressions=NULL ,matureExpressions=matExp_DESeq2)
A method to obtain the number of the genes associated with the object of class INSPEcT
nGenes(object) ## S4 method for signature 'INSPEcT' nGenes(object)
nGenes(object) ## S4 method for signature 'INSPEcT' nGenes(object)
object |
An object of class INSPEcT |
A numeric that indicates the number of genes within the object
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nGenes(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nGenes(nascentInspObj10)
A method to obtain the number of the tpts associated with the object of class INSPEcT
nTpts(object) ## S4 method for signature 'INSPEcT' nTpts(object)
nTpts(object) ## S4 method for signature 'INSPEcT' nTpts(object)
object |
An object of class INSPEcT |
A numeric that indicates the number of time points contained the object
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nTpts(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nTpts(nascentInspObj10)
A method to see the shapes of the estimated synthesis, degradation and processing rates, pre-RNA and total RNA concentrations (solid thin lines) their variances (dashed lines) and the modeled rates and concentrations (ticker solid line) of a single gene.
plotGene( object, ix, relative_expression = FALSE, fix.yaxis = FALSE, priors = TRUE, constantModel = FALSE ) ## S4 method for signature 'INSPEcT' plotGene( object, ix, relative_expression = FALSE, fix.yaxis = FALSE, priors = TRUE, constantModel = FALSE )
plotGene( object, ix, relative_expression = FALSE, fix.yaxis = FALSE, priors = TRUE, constantModel = FALSE ) ## S4 method for signature 'INSPEcT' plotGene( object, ix, relative_expression = FALSE, fix.yaxis = FALSE, priors = TRUE, constantModel = FALSE )
object |
An object of class INSPEcT |
ix |
Eiher a rowname or a row number to select one single gene |
relative_expression |
A logical, indicating whether expressions are rates should be plotted relative to their initial value (Default=FALSE). |
fix.yaxis |
A logical, indicating whether the limits for y-axis of degradation and processing rates should be fixed relative to their distributions |
priors |
A logical, if true the priors of the rates are plotted |
constantModel |
A logical, if true the constant model for the + nascent modeling are shown |
A list containing total RNA levels and their confidence interval (levels plus and minus one standard deviation), pre-RNA lelevs and their confidence intervals, synthsis rates and their confidence intervals, degradation rates and processing rates of the selected gene.
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) plotGene(nascentInspObj10, 1)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) plotGene(nascentInspObj10, 1)
Visualize the comparison between the rates calculated from two different INSPEcT objects profiled in steady-state conditions.
## S4 method for signature 'INSPEcT_diffsteady' plotMA(object, ...)
## S4 method for signature 'INSPEcT_diffsteady' plotMA(object, ...)
object |
An object of calss INSPEcT_diffsteady |
... |
Additional parameters, see Details section |
Possible arguments to "plotMA":
"rate" - A character, which represent the rate to be visualized, either "synthesis", "processing" or "degradation". By default, "synthesis" is chosen.
"padj" - A numeric, The p-adjusted threshold for significance. Genes with p-adjusted lower than the threshold will be depicted as orange triangles. By default set to -Inf, meaning that no genes will be highlighted.
"xlim" - A numeric vector of length 2, limits of x-axis, by default the range of the data.
"xlab" - A character, the label of x-axis, by default "log2 geometric mean"
"ylim" - A numeric vector of length 2, limits of y-axis, by default the range of the data.
"ylab" - A character, the label of y-axis, by default "log2 fold change"
"main" - A character, the title of the plot, by default the name of the visualized rate.
http://en.wikipedia.org/wiki/MA_plot
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT(tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) plotMA(diffrates, padj=.01) }
if( Sys.info()["sysname"] != "Windows" ) { data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-letters[1:11] expDes<-rep(conditions,3) tL<-1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=nascentCounts ,libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) nasFullObj <- newINSPEcT(tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) diffrates = compareSteady(nasFullObj[,c(1,11)]) plotMA(diffrates, padj=.01) }
Plot the premature and mature expressions of a specific gene in the different samples of the dataset along with the null model and the log2 fold change threshold. Individal observations that fall outside of the dashed lines are considered post-transcriptional events.
plotPMgene(object, gene_id, samples_colors = 1) ## S4 method for signature 'INSPEcT_steadyNoNascent' plotPMgene(object, gene_id, samples_colors = 1)
plotPMgene(object, gene_id, samples_colors = 1) ## S4 method for signature 'INSPEcT_steadyNoNascent' plotPMgene(object, gene_id, samples_colors = 1)
object |
An object of class INSPEcT_steadyNoNascent |
gene_id |
A numeric that indicated the index of the gene to be plotted |
samples_colors |
The color code relative to the samples |
Plot the null model estimated for the specific dataset. The null model is the trend between premature and mature expression, which is usually linear in the log-log scale and generally points to an increase in the ratio between premature and mature RNA at increased levels of expression
plotPMtrend(inspectIds) ## S4 method for signature 'INSPEcT_steadyNoNascent' plotPMtrend(inspectIds)
plotPMtrend(inspectIds) ## S4 method for signature 'INSPEcT_steadyNoNascent' plotPMtrend(inspectIds)
inspectIds |
An object of class INSPEcT_steadyNoNascent |
Extract premature RNA expressions
premature(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' premature(object)
premature(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' premature(object)
object |
An object of class INSPEcT_steadyNoNascent |
A matrix containing premature RNA expressions
Extract premature RNA expressions variances
prematureVar(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' prematureVar(object)
prematureVar(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' prematureVar(object)
object |
An object of class INSPEcT_steadyNoNascent |
A matrix containing premature RNA expressions variances
These functions calculates the tau and delta metrics for all genes with introns and exons in an oblect of class INSPEcT. If the INSPEcT dataset was obtained with nascent RNA the metrics are caluclated using RNA dynamics and solving numerically the system of equations. If the INSPEcT dataset was obtained without nascent RNA the metrics are approximated using premature and mature levels.
processingDelay( inspectIds, tauThreshold = 1.2, deltaThreshold = 1, silent = TRUE ) calculateDelta(inspectIds, silent = FALSE) calculateTau(inspectIds, silent = FALSE) ## S4 method for signature 'INSPEcT' processingDelay( inspectIds, tauThreshold = 1.2, deltaThreshold = 1, silent = TRUE ) ## S4 method for signature 'INSPEcT' calculateTau(inspectIds, silent = FALSE) ## S4 method for signature 'INSPEcT' calculateDelta(inspectIds, silent = FALSE)
processingDelay( inspectIds, tauThreshold = 1.2, deltaThreshold = 1, silent = TRUE ) calculateDelta(inspectIds, silent = FALSE) calculateTau(inspectIds, silent = FALSE) ## S4 method for signature 'INSPEcT' processingDelay( inspectIds, tauThreshold = 1.2, deltaThreshold = 1, silent = TRUE ) ## S4 method for signature 'INSPEcT' calculateTau(inspectIds, silent = FALSE) ## S4 method for signature 'INSPEcT' calculateDelta(inspectIds, silent = FALSE)
inspectIds |
An object of class INSPEcT. |
tauThreshold |
A numeric representing the tau threshold to define a gene affected by processing. Default: 1.2 |
deltaThreshold |
A numeric representing the delta threshold to define a gene affected by processing. Default: 1.0 |
silent |
A logical indicating whether informaiton about the procedure should be printed or not. |
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16) expDes<-rep(conditions,3) tL <- 1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj <- newINSPEcT( tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) procDelay<- processingDelay(inspectIds=matureInspObj ,tauThreshold=1.2 ,deltaThreshold=1.0) head(procDelay) table(procDelay) head(calculateTau(matureInspObj)) head(calculateDelta(matureInspObj))
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature conditions<-c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16) expDes<-rep(conditions,3) tL <- 1/6 nasExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts( allcounts=matureCounts ,libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,experimentalDesign=expDes) matureInspObj <- newINSPEcT( tpts=conditions ,labeling_time=tL ,nascentExpressions=nasExp_DESeq2 ,matureExpressions=matExp_DESeq2) procDelay<- processingDelay(inspectIds=matureInspObj ,tauThreshold=1.2 ,deltaThreshold=1.0) head(procDelay) table(procDelay) head(calculateTau(matureInspObj)) head(calculateDelta(matureInspObj))
Extract the ratio between mature and premature RNAs
PTratio(object, infToNA = TRUE) ## S4 method for signature 'INSPEcT_steadyNoNascent' PTratio(object, infToNA = TRUE)
PTratio(object, infToNA = TRUE) ## S4 method for signature 'INSPEcT_steadyNoNascent' PTratio(object, infToNA = TRUE)
object |
An object of class INSPEcT_steadyNoNascent |
infToNA |
A logical indicating whether infinite values (originating from zero valued premature expressions) should be set artificially to NA or not |
A matrix containing the PTratios
Extract the post-transcriptional regulation matrix
PTreg(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' PTreg(object)
PTreg(object) ## S4 method for signature 'INSPEcT_steadyNoNascent' PTreg(object)
object |
An object of class INSPEcT_steadyNoNascent |
A matrix containing the post-transcriptional regulated genes. This matrix is generated by the method compareSteadyStateNoNascent. It generally report 1 for regulated genes in specific samples, 0 for non regulated genes and NA for genes that do not passed the expression threshold. In case the argument returnNormScores was set to TRUE, instead of discretes values, the deviations from the expected model normalized by the experimental standard deviation is reported.
Given a TranscriptDb object and a list of BAM or SAM files "quantifyExpressionsFormBAMs" evaluates exons and introns expressions and the associated variances per each gene.
quantifyExpressionsFromBAMs( txdb, BAMfiles, experimentalDesign, by = c("gene", "tx"), countMultiMappingReads = FALSE, allowMultiOverlap = FALSE, prioritizeExons = TRUE, libsize = c("assigned", "all"), strandSpecific = 0, isPairedEnd = FALSE, DESeq2 = TRUE, varSamplingCondition = NULL, BPPARAM = SerialParam() )
quantifyExpressionsFromBAMs( txdb, BAMfiles, experimentalDesign, by = c("gene", "tx"), countMultiMappingReads = FALSE, allowMultiOverlap = FALSE, prioritizeExons = TRUE, libsize = c("assigned", "all"), strandSpecific = 0, isPairedEnd = FALSE, DESeq2 = TRUE, varSamplingCondition = NULL, BPPARAM = SerialParam() )
txdb |
A TranscriptDB object |
BAMfiles |
A vector of paths |
experimentalDesign |
A numerical which reports the desing of the experiment in terms of time points and replicates. Time points must be ordered according to the sequence of files submitted for the analysis, these labels characterize different files as replicates of a given condition. |
by |
A character, either "gene" or "tx", indicating if expressions and counts should be summarized at the levels of genes or transcripts. "gene" by default. In case "tx" is selected, we suggest to set argument "allowMultiOverlap" to TRUE, otherwise the reads mapping to overlapping transcripts of the same gene will remain unassigned. |
countMultiMappingReads |
A logical, if multimapping reads should be counted, FALSE by default. Multimap reads are identified using the tag "NH" in the bam/sam file. |
allowMultiOverlap |
A logical, indicating if a read is allowed to be assigned to more than one feature, FALSE by default |
prioritizeExons |
A logical, indicating whether reads assigned to exon shold not be accounted for intron counts. If set to FALSE, reads with shared overlap between an exon and the following intron will be assigned also to introns. This could improve intronic quantification in experimental settings (including polyA library preparation) or compact genomes were intronic reads are sampled at a very low rate compared to exonic reads. By default, TRUE. |
libsize |
A character, either "assigned" or "all", indicating whether the libsize for expression normalization should include all mapped reads or only the reads assigned to any of the features. By default, "assigned" is selected. |
strandSpecific |
Numeric, 0 if no strand-specific read counting should be performed, 1 stranded, 2 reversely-stranded. 0 by default |
isPairedEnd |
A logical, if paired-end reads are used, FALSE by default |
DESeq2 |
A logical, if TRUE exons and introns variances are evaluated through the package DESeq2, if FALSE through plgem |
varSamplingCondition |
A character reporting which experimental condition should be used to sample the variance if DESeq2 = FALSE. |
BPPARAM |
Parallelization parameters for bplapply. By default SerialParam() By default, the first element of "experimentalDesign" with replicates. |
A list containing expressions and associated variances for exons and introns.
if( Sys.info()["sysname"] != "Windows" ) { require(TxDb.Mmusculus.UCSC.mm9.knownGene) txdb<-TxDb.Mmusculus.UCSC.mm9.knownGene expDes<-c(0,0,1,1) paths_total<-system.file('extdata/', c('bamRep1.bam' ,'bamRep2.bam' ,'bamRep3.bam' ,'bamRep4.bam') ,package='INSPEcT') matExp<-quantifyExpressionsFromBAMs(txdb=txdb ,BAMfiles=paths_total ,experimentalDesign=expDes) }
if( Sys.info()["sysname"] != "Windows" ) { require(TxDb.Mmusculus.UCSC.mm9.knownGene) txdb<-TxDb.Mmusculus.UCSC.mm9.knownGene expDes<-c(0,0,1,1) paths_total<-system.file('extdata/', c('bamRep1.bam' ,'bamRep2.bam' ,'bamRep3.bam' ,'bamRep4.bam') ,package='INSPEcT') matExp<-quantifyExpressionsFromBAMs(txdb=txdb ,BAMfiles=paths_total ,experimentalDesign=expDes) }
Given a TranscriptDb object and a list of bigWig (BW) files "quantifyExpressionsFormBWs" evaluates exons and introns expressions and the associated variances per each gene.
quantifyExpressionsFromBWs( txdb, BWfiles, experimentalDesign, readLength = 50, by = c("gene", "tx"), libsize = c("assigned", "all"), DESeq2 = TRUE, varSamplingCondition = NULL, BPPARAM = SerialParam() )
quantifyExpressionsFromBWs( txdb, BWfiles, experimentalDesign, readLength = 50, by = c("gene", "tx"), libsize = c("assigned", "all"), DESeq2 = TRUE, varSamplingCondition = NULL, BPPARAM = SerialParam() )
txdb |
A TranscriptDB object |
BWfiles |
A vector of paths |
experimentalDesign |
A numerical which reports the desing of the experiment in terms of time points and replicates. Time points must be ordered according to the sequence of files submitted for the analysis, these labels characterize different files as replicates of a given condition. |
readLength |
A numerical that indicates the read length of the RNA-seq experiment. Used to normalize the coverage. By default, 50. |
by |
A character, either "gene" or "tx", indicating if expressions and counts should be summarized at the levels of genes or transcripts. "gene" by default. In case "tx" is selected, we suggest to set argument "allowMultiOverlap" to TRUE, otherwise the reads mapping to overlapping transcripts of the same gene will remain unassigned. |
libsize |
A character, either "assigned" or "all", indicating whether the libsize for expression normalization should include all mapped reads or only the reads assigned to any of the features. By default, "assigned" is selected. |
DESeq2 |
A logical, if TRUE exons and introns variances are evaluated through the package DESeq2, if FALSE through plgem |
varSamplingCondition |
A character reporting which experimental condition should be used to sample the variance if DESeq2 = FALSE. |
BPPARAM |
Parallelization parameters for bplapply. By default SerialParam() By default, the first element of "experimentalDesign" with replicates. |
A list containing expressions and associated variances for exons and introns.
Given introns and exons abundances (for example RPKMs) this method returns their variances evaluated thorugh plgem.
quantifyExpressionsFromTrAbundance( trAbundaces, experimentalDesign, varSamplingCondition = NULL, simulatedData = FALSE )
quantifyExpressionsFromTrAbundance( trAbundaces, experimentalDesign, varSamplingCondition = NULL, simulatedData = FALSE )
trAbundaces |
A a list with elements "exonsAbundances" and "intronsAbundances". |
experimentalDesign |
A numerical which reports the desing of the experiment in terms of time points and replicates. The time points must be ordered according to the columns of the count matrices submitted for the analysis; these labels define conditions and replicates. |
varSamplingCondition |
A character reporting which experimental condition should be used to sample the variance if DESeq2 = FALSE. |
simulatedData |
A boolean which is TRUE if the data under analysis are simulated. |
A list containing RPKMs and associated variances for exons and introns.
Evaluates introns and exons RPKMs, per gene, from counts data.
quantifyExpressionsFromTrCounts( allcounts, experimentalDesign, exonsWidths, intronsWidths, libsize = NULL, DESeq2 = TRUE, varSamplingCondition = NULL )
quantifyExpressionsFromTrCounts( allcounts, experimentalDesign, exonsWidths, intronsWidths, libsize = NULL, DESeq2 = TRUE, varSamplingCondition = NULL )
allcounts |
A named list containing "exonsCounts" and "intronsCounts". |
experimentalDesign |
A numerical which reports the desing of the experiment in terms of time points and replicates. Time points must be ordered according to the sequence of files submitted for the analysis, these labels characterize different files as replicates of a given condition. |
exonsWidths |
A numeric containing the exons widths. |
intronsWidths |
A numeric containing the intorns widths. |
libsize |
A numeric containing the library size. |
DESeq2 |
A logical, if TRUE the RPKMs variances are evaluated through the package DESeq2, if FALSE plgem is used. |
varSamplingCondition |
A character reporting which experimental condition should be used to sample the variance if DESeq2 = FALSE. By default, the first element of "experimentalDesign" with replicates. |
A list containing RPKMs and associated variances for exons and introns.
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature expDes<-rep(c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16),3) nasExp_DESeq2<-quantifyExpressionsFromTrCounts(libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=nascentCounts ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts(libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=matureCounts ,experimentalDesign=expDes) nasExp_plgem<-quantifyExpressionsFromTrCounts(libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=nascentCounts ,DESeq2=FALSE ,experimentalDesign=expDes) matExp_plgem<-quantifyExpressionsFromTrCounts(libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=matureCounts ,DESeq2=FALSE ,experimentalDesign=expDes)
data('allcounts', package='INSPEcT') data('featureWidths', package='INSPEcT') data('libsizes', package='INSPEcT') nascentCounts<-allcounts$nascent matureCounts<-allcounts$mature expDes<-rep(c(0,1/6,1/3,1/2,1,1.5,2,4,8,12,16),3) nasExp_DESeq2<-quantifyExpressionsFromTrCounts(libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=nascentCounts ,experimentalDesign=expDes) matExp_DESeq2<-quantifyExpressionsFromTrCounts(libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=matureCounts ,experimentalDesign=expDes) nasExp_plgem<-quantifyExpressionsFromTrCounts(libsize=nascentLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=nascentCounts ,DESeq2=FALSE ,experimentalDesign=expDes) matExp_plgem<-quantifyExpressionsFromTrCounts(libsize=totalLS ,exonsWidths=exWdths ,intronsWidths=intWdths ,allcounts=matureCounts ,DESeq2=FALSE ,experimentalDesign=expDes)
This method is used to retrieve all the p-values relative to the variability of synthesis, processing and degradation rates.
ratePvals(object) ## S4 method for signature 'INSPEcT' ratePvals(object)
ratePvals(object) ## S4 method for signature 'INSPEcT' ratePvals(object)
object |
An object of class INSPEcT |
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratePvals(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratePvals(nascentInspObj10)
This method allow to access to the estimated synthesis, degradation, processing rates and pre mRNA and total mRNA
concentrations the way they were calculated by the constructor function newINSPEcT
.
ratesFirstGuess(object, feature) ## S4 method for signature 'INSPEcT' ratesFirstGuess(object, feature)
ratesFirstGuess(object, feature) ## S4 method for signature 'INSPEcT' ratesFirstGuess(object, feature)
object |
An object of class INSPEcT |
feature |
A character indicating the feature to retireve, "synthesis", "degradation", "processing" for rates, "total" for total mRNA concentrations or "preMRNA" for premature mRNA concentrations |
A numeric matrix containing the values for the selected feature
newINSPEcT
, ratesFirstGuessVar
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratesFirstGuess(nascentInspObj10, 'total') ratesFirstGuess(nascentInspObj10, 'preMRNA') ratesFirstGuess(nascentInspObj10, 'synthesis')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratesFirstGuess(nascentInspObj10, 'total') ratesFirstGuess(nascentInspObj10, 'preMRNA') ratesFirstGuess(nascentInspObj10, 'synthesis')
This method allow to access to the estimated variance of synthesis rates and pre mRNA and total mRNA
concentrations the way they were calculated by the constructor function newINSPEcT
.
ratesFirstGuessVar(object, feature) ## S4 method for signature 'INSPEcT' ratesFirstGuessVar(object, feature)
ratesFirstGuessVar(object, feature) ## S4 method for signature 'INSPEcT' ratesFirstGuessVar(object, feature)
object |
An object of class INSPEcT |
feature |
A character indicating the feature to retireve, "synthesis", "degradation", "processing" for rates, "total" for total mRNA concentrations or "preMRNA" for premature mRNA concentrations |
A numeric vector containing the values for the selected feature
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratesFirstGuessVar(nascentInspObj10, 'total') ratesFirstGuessVar(nascentInspObj10, 'preMRNA') ratesFirstGuessVar(nascentInspObj10, 'synthesis')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) ratesFirstGuessVar(nascentInspObj10, 'total') ratesFirstGuessVar(nascentInspObj10, 'preMRNA') ratesFirstGuessVar(nascentInspObj10, 'synthesis')
Remove the model from an INSPEcT object. It is required when subsetting an INSPEcT object per time points because when removing time points the modeling is not valid anymore.
removeModel(object) ## S4 method for signature 'INSPEcT' removeModel(object)
removeModel(object) ## S4 method for signature 'INSPEcT' removeModel(object)
object |
An Object of class INSPEcT |
An Object of class INSPEcT
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_5genes <- nascentInspObj10[1:5] ## This will turn out into an error: ## Not run: nascentInspObj10_5genes_5tpts <- nascentInspObj10[1:5, 1:5] ## Before subsetting time points, the model should be removed: nascentInspObj10_5genes_5tpts <- removeModel(nascentInspObj10)[1:5, 1:5] ## Also this will turn out into an error: ## Not run: nascentInspObj10 <- modelRates(nascentInspObj10) ## Before running the model again, or changing modeling parameters, ## the previous model should be removed: nascentInspObj10_old <- nascentInspObj10 nascentInspObj10_new <- removeModel(nascentInspObj10) ## Not run: nascentInspObj10_new <- modelRates(nascentInspObj10_new, useSigmoidFun = FALSE)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) nascentInspObj10_5genes <- nascentInspObj10[1:5] ## This will turn out into an error: ## Not run: nascentInspObj10_5genes_5tpts <- nascentInspObj10[1:5, 1:5] ## Before subsetting time points, the model should be removed: nascentInspObj10_5genes_5tpts <- removeModel(nascentInspObj10)[1:5, 1:5] ## Also this will turn out into an error: ## Not run: nascentInspObj10 <- modelRates(nascentInspObj10) ## Before running the model again, or changing modeling parameters, ## the previous model should be removed: nascentInspObj10_old <- nascentInspObj10 nascentInspObj10_new <- removeModel(nascentInspObj10) ## Not run: nascentInspObj10_new <- modelRates(nascentInspObj10_new, useSigmoidFun = FALSE)
A method to visualize the performance in the classification of synthesis, degradation
and processing rates based on the comparison of the original simulated rates and the one
obtained by the function modelRates
. For each rate, classification performance is measured
in terms of sensitivity and specificity using a ROC curve analysis. False negatives (FN) represent cases
where the rate is identified as constant while it was simulated as varying. False positives (FP) represent
cases where INSPEcT identified a rate as varying while it was simulated as constant. On the contrary,
true positives (TP) and negatives (TN) are cases of correct classification of varying and constant rates, respectively.
Consequently, sensitivity and specificity are computed using increasing thresholds for the brown p-values,
and the ability of correctly classifying a rate is measured through the area under the curve (AUC) for each rate.
rocCurve(object, object2, plot = TRUE, comparative = FALSE) ## S4 method for signature 'INSPEcT_model,INSPEcT' rocCurve(object, object2, plot = TRUE, comparative = FALSE)
rocCurve(object, object2, plot = TRUE, comparative = FALSE) ## S4 method for signature 'INSPEcT_model,INSPEcT' rocCurve(object, object2, plot = TRUE, comparative = FALSE)
object |
An object of class INSPEcT_model, with true rates |
object2 |
An modeled object of class INSPEcT |
plot |
A logical indicating whether ROC curves should be plotted or not |
comparative |
A logical indicating whether the cross-prediction should be visualized. When this mode is selected, the p-values assigned to the variability of one rate (e.g. synthesis) are tested against the variability the other rates (e.g. processing and degradation). Cross-prediction ROC curves are plotted with dashed lines. |
A list of objects of class pROC with summary of each roc curve
makeSimModel
, makeSimDataset
, rocThresholds
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) # newTpts<-simRates@params$tpts # nascentSim2replicates<-makeSimDataset(object=simRates # ,tpts=newTpts # ,nRep=3 # ,NoNascent=FALSE # ,seed=1) # nascentSim2replicates<-modelRates(nascentSim2replicates[1:100] # ,seed=1) # (not evaluated to save computational time) data("nascentSim2replicates",package='INSPEcT') rocCurve(simRates[1:100],nascentSim2replicates) title("3rep. 11t.p. Total and nascent RNA", line=3) }
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) # newTpts<-simRates@params$tpts # nascentSim2replicates<-makeSimDataset(object=simRates # ,tpts=newTpts # ,nRep=3 # ,NoNascent=FALSE # ,seed=1) # nascentSim2replicates<-modelRates(nascentSim2replicates[1:100] # ,seed=1) # (not evaluated to save computational time) data("nascentSim2replicates",package='INSPEcT') rocCurve(simRates[1:100],nascentSim2replicates) title("3rep. 11t.p. Total and nascent RNA", line=3) }
A method to visualize the performance in the classification of synthesis, degradation
and processing rates based on the comparison of the original simulated rates and the one
obtained by the function modelRates
. For each rate, classification performance is measured
in terms of sensitivity and specificity using a ROC curve analysis. False negatives (FN) represent cases
where the rate is identified as constant while it was simulated as varying. False positives (FP) represent
cases where INSPEcT identified a rate as varying while it was simulated as constant. On the contrary,
true positives (TP) and negatives (TN) are cases of correct classification of varying and constant rates, respectively.
Consequently, at increasing brown p-values different sensitivity and specificity can be achieved.
rocThresholds(object, object2, xlim = c(1e-05, 1), plot = TRUE) ## S4 method for signature 'INSPEcT_model,INSPEcT' rocThresholds(object, object2, xlim = c(1e-05, 1), plot = TRUE)
rocThresholds(object, object2, xlim = c(1e-05, 1), plot = TRUE) ## S4 method for signature 'INSPEcT_model,INSPEcT' rocThresholds(object, object2, xlim = c(1e-05, 1), plot = TRUE)
object |
An object of class INSPEcT_model, with true rates |
object2 |
An object of class INSPEcT or INSPEcT_model, with modeled rates |
xlim |
A numeric representing limits for the x-axis (default is c(1-e-5,1)) |
plot |
A logical that indicates whether to plot or not. (default=TRUE) |
The thresholds that maximize both sensitivity and specificity
makeSimModel
, makeSimDataset
, rocCurve
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) # newTpts<-simRates@params$tpts # nascentSim2replicates<-makeSimDataset(object=simRates # ,tpts=newTpts # ,nRep=3 # ,NoNascent=FALSE # ,seed=1) # nascentSim2replicates<-modelRates(nascentSim2replicates[1:100] # ,seed=1) # (not evaluated to save computational time) data("nascentSim2replicates",package='INSPEcT') rocThresholds(simRates[1:100],nascentSim2replicates) }
if( Sys.info()["sysname"] != "Windows" ) { nascentInspObj <- readRDS(system.file(package='INSPEcT', 'nascentInspObj.rds')) simRates<-makeSimModel(nascentInspObj, 1000, seed=1) # newTpts<-simRates@params$tpts # nascentSim2replicates<-makeSimDataset(object=simRates # ,tpts=newTpts # ,nRep=3 # ,NoNascent=FALSE # ,seed=1) # nascentSim2replicates<-modelRates(nascentSim2replicates[1:100] # ,seed=1) # (not evaluated to save computational time) data("nascentSim2replicates",package='INSPEcT') rocThresholds(simRates[1:100],nascentSim2replicates) }
This function is used to set the confidence intervals in the nascent RNA mode.
setConfidenceIntervals(object, confidenceIntervals) ## S4 method for signature 'INSPEcT' setConfidenceIntervals(object, confidenceIntervals)
setConfidenceIntervals(object, confidenceIntervals) ## S4 method for signature 'INSPEcT' setConfidenceIntervals(object, confidenceIntervals)
object |
An object of class INSPEcT_model |
confidenceIntervals |
list of confidence intervals. |
An object of class ExpressionSet containing the confidence intervals.
A dataset containing the rates and concentrations obtained from the dataset simRates; 3 replicates and time points corresponding to: 0,1/6,1/3,1/2,1,1.5,2,4,8,12,16 hours.
An INSPEcT object
A dataset containing the rates and concentrations obtained from the dataset simRates with 1 replicates and time points corresponding to: 0, 1/6, 1/3, 1/2, 1, 2, 4, 8, 16 hours. On this dataset rates and concentrations have been modeled with the method modelRates
An INSPEcT object
A dataset containing the rates and concentrations obtained from the dataset simRates; 3 replicates and time points corresponding to: 0,1/6,1/3,1/2,1,1.25,1.5,2,3,4,6,8,10,12,16 hours.
An INSPEcT object
A dataset containing the rates and concentrations obtained from the dataset simRates with 1 replicates and time points corresponding to: 0, 1/6, 1/3, 1/2, 1, 2, 4, 8, 16 hours. On this dataset rates and concentrations have been modeled with the method modelRates
An INSPEcT object
Divides the INSPEcT object into the groups defined by 'f',
## S4 method for signature 'INSPEcT,ANY' split(x, f, drop = FALSE, ...)
## S4 method for signature 'INSPEcT,ANY' split(x, f, drop = FALSE, ...)
x |
An object of class INSPEcT |
f |
A vector of length equal to the number of genes in x which defines the groups |
drop |
A logical belonging to the generic funciton, useless in this context. |
... |
Additional arguments to match the generic function |
A list containing objects of class INSPEcT
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) splitIdx <- c(1,1,1,2,2,2,3,3,3,4) nascentInspObj10Split <- split(nascentInspObj10, splitIdx)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) splitIdx <- c(1,1,1,2,2,2,3,3,3,4) nascentInspObj10Split <- split(nascentInspObj10, splitIdx)
Accessor to obtain the tpts associated with the object of class INSPEcT
tpts(object) ## S4 method for signature 'INSPEcT' tpts(object)
tpts(object) ## S4 method for signature 'INSPEcT' tpts(object)
object |
An object of class INSPEcT |
A numeric that indicates time points contained the object
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) tpts(nascentInspObj10)
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) tpts(nascentInspObj10)
A method to access the modeld confidence intervals computed via the method computeConfidenceIntervals
viewConfidenceIntervals(object, feature) ## S4 method for signature 'INSPEcT' viewConfidenceIntervals(object, feature)
viewConfidenceIntervals(object, feature) ## S4 method for signature 'INSPEcT' viewConfidenceIntervals(object, feature)
object |
An object of class INSPEcT |
feature |
A character indicating the feature to retireve: "synthesis", "degradation", "processing". |
A numeric matrix containing the values for the selected feature
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewConfidenceIntervals(nascentInspObj10, 'synthesis')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewConfidenceIntervals(nascentInspObj10, 'synthesis')
A method to access the modeld rates via the method modelRates
viewModelRates(object, feature) ## S4 method for signature 'INSPEcT' viewModelRates(object, feature)
viewModelRates(object, feature) ## S4 method for signature 'INSPEcT' viewModelRates(object, feature)
object |
An object of class INSPEcT |
feature |
A character indicating the feature to retireve, "synthesis", "degradation", "processing" for rates, "total" for total mRNA concentrations or "preMRNA" for premature mRNA concentrations |
A numeric matrix containing the values for the selected feature
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewModelRates(nascentInspObj10, 'synthesis')
nascentInspObj10 <- readRDS(system.file(package='INSPEcT', 'nascentInspObj10.rds')) viewModelRates(nascentInspObj10, 'synthesis')