Title: | R Normalization and Inference of Time Series data |
---|---|
Description: | R/Bioconductor package for normalization, curve registration and inference in time course gene expression data. |
Authors: | Dipen P. Sangurdekar <[email protected]> |
Maintainer: | Dipen P. Sangurdekar <[email protected]> |
License: | GPL-3 |
Version: | 1.41.0 |
Built: | 2024-12-18 03:55:44 UTC |
Source: | https://github.com/bioc/Rnits |
This function takes high-dimensional expression data as a RGList, creates a
Rnits
object for subsequent filtering and normalization
build.Rnits(obj, probedata = NULL, phenodata = NULL, filter = NULL, normalize = NULL, normmethod = NULL, plot = FALSE, center = FALSE, background = NULL, threshold = 0.8, logscale = FALSE)
build.Rnits(obj, probedata = NULL, phenodata = NULL, filter = NULL, normalize = NULL, normmethod = NULL, plot = FALSE, center = FALSE, background = NULL, threshold = 0.8, logscale = FALSE)
obj |
Raw expression data in |
probedata |
A data frame containing the probe names that should match the probe names in raw data (optional) |
phenodata |
A data frame with information about sample names. The rownames of the data frame must match column names of the expression values. If input data is data frame of log ratios, this is required. |
filter |
An argument to perform background filtering of probes. If |
normalize |
Character string specifying the normalization method for raw data.
If |
normmethod |
Normalization method for input data. Default |
background |
Only for AffyBatch data. If |
center |
If |
plot |
If |
threshold |
Default |
logscale |
Default |
See the Limma User's Guide for more details on read.maimages
, normalizeBetweenArrays
,
normalizeWithinArrays
and RGList
. For importing microarray raw data,
use the 'Targets file' to specify experimental design. The target file has columns
SlideNumber, FileName, Cy3 (description of Cy3 channel ref/control/treatment), Cy5
(description of Cy3 channel ref/control/treatment) and Time. Time values should be
identical for control and treatment.
An object of S4 class Rnits
(which is
derived from class exprSet
), containing the probe data,
design data, expression data, phenotypical data (i.e. Time).
ExpressionSet
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between')
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between')
Calculate the optimal B-spline model using generalized cross-validation
calculateGCV(object, topcomp = 5) ## S4 method for signature 'Rnits' calculateGCV(object, topcomp = 5)
calculateGCV(object, topcomp = 5) ## S4 method for signature 'Rnits' calculateGCV(object, topcomp = 5)
object |
|
topcomp |
The number of top eigenvectors to be used for computation |
The optimal B-spline model is chosen as the largest model that minimizes the cross validation error of the top N eigenvectors of each time series data.
A list object with fields 'degree', 'df' for each time series data set.
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') opt_model <- calculateGCV(rnitsobj) ## Not run: rnitsobj <- fit(rnitsobj, gene.level = TRUE, model = opt.model) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') opt_model <- calculateGCV(rnitsobj) ## Not run: rnitsobj <- fit(rnitsobj, gene.level = TRUE, model = opt.model) ## End(Not run)
Fit a model comparing time series data set Rnits
objects
fit(object, cluster = TRUE, B = 100, verbatim = FALSE, nclus = NULL, modelhistplot = FALSE, seed = 123, gene.level = TRUE, clusterallsamples = FALSE, model = NULL) ## S4 method for signature 'Rnits' fit(object, cluster = TRUE, B = 100, verbatim = FALSE, nclus = NULL, modelhistplot = FALSE, seed = 123, gene.level = TRUE, clusterallsamples = FALSE, model = NULL)
fit(object, cluster = TRUE, B = 100, verbatim = FALSE, nclus = NULL, modelhistplot = FALSE, seed = 123, gene.level = TRUE, clusterallsamples = FALSE, model = NULL) ## S4 method for signature 'Rnits' fit(object, cluster = TRUE, B = 100, verbatim = FALSE, nclus = NULL, modelhistplot = FALSE, seed = 123, gene.level = TRUE, clusterallsamples = FALSE, model = NULL)
object |
|
cluster |
if |
B |
Default |
verbatim |
If |
nclus |
Default |
modelhistplot |
If |
seed |
Random seed for bootstrap iterations |
gene.level |
If |
clusterallsamples |
If |
model |
A data frame with fields 'degree' and 'df' indicating a specific B-spline model to be used. If provided, model selection is not run. |
The function compares multiple time-series expression data sets by i) (optional)
summarizing probes into gene-level information ii) (optional) identifying a set
of co-expressed genes by clustering iii) For each cluster (or for all genes
/probes), fit a series of B-splines with varying curvature and degrees of
freedom. Under the null hypothesis H_0
, a single model is fit for
all data sets, while under H_1
, each data set is fit separately.
P-values from the hypothesis test are then plotted and the least complex
spline parameters that result in uniformly distributed null p-values are
automatically chosen.
An object of S4 class Rnits
with fitted results
data containing cluster information, ratio statistics and p-values.
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) ## End(Not run)
Rnits
Retrieve cluster IDs of probes/genes from fitted Rnits
object after fit
has been run.
getCID(object) ## S4 method for signature 'Rnits' getCID(object)
getCID(object) ## S4 method for signature 'Rnits' getCID(object)
object |
If cluster = False
during fitting, a vector of 1s
will be returned.
A vector of cluster IDs corresponding to gene/probe names
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get cluster IDs from fitted model cid <- getCID(rnitsobj) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get cluster IDs from fitted model cid <- getCID(rnitsobj) ## End(Not run)
Rnits
objectRetrieve model fit data from Rnits
object after fit
has been run.
getFitModel(object) ## S4 method for signature 'Rnits' getFitModel(object)
getFitModel(object) ## S4 method for signature 'Rnits' getFitModel(object)
object |
Contains Ratio statistic, p-value and cluster ID data
A data frame containing the model fit results for all genes
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # P-values, ratio statistics and cluster ID's can be retrieved for all genes together fitdata <- getFitModel(rnitsobj) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # P-values, ratio statistics and cluster ID's can be retrieved for all genes together fitdata <- getFitModel(rnitsobj) ## End(Not run)
Extract normalized log-ratios from Rnits
object
getLR(object, impute = FALSE) ## S4 method for signature 'Rnits' getLR(object, impute = FALSE)
getLR(object, impute = FALSE) ## S4 method for signature 'Rnits' getLR(object, impute = FALSE)
object |
|
impute |
If |
A matrix of normalized log-ratios.
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get logratios lr <- getLR(rnitsobj) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get logratios lr <- getLR(rnitsobj) ## End(Not run)
For two color data, extract normalized channel data from Rnits
object
getNormTwoChannel(object) ## S4 method for signature 'Rnits' getNormTwoChannel(object)
getNormTwoChannel(object) ## S4 method for signature 'Rnits' getNormTwoChannel(object)
object |
|
A list containing R and G fields for normalized Red and Green channel data respectively.
Extract p-values from fitted Rnits
object
getPval(object) ## S4 method for signature 'Rnits' getPval(object)
getPval(object) ## S4 method for signature 'Rnits' getPval(object)
object |
|
An vector of p-values
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) #Get pvalues from fitted model pval <- getPval(rnitsobj) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) #Get pvalues from fitted model pval <- getPval(rnitsobj) ## End(Not run)
Extract ratio statistics from fitted Rnits
object
getStat(object) ## S4 method for signature 'Rnits' getStat(object)
getStat(object) ## S4 method for signature 'Rnits' getStat(object)
object |
|
An vector of ratio statistics
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get ratio statistics from fitted model stat <- getStat(rnitsobj) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get ratio statistics from fitted model stat <- getStat(rnitsobj) ## End(Not run)
After fit
has been applied on Rnits
object, plot the profiles
of N top ranking genes/probes.
plotResults(object, id = NULL, fdr = NULL, top = 48, pdf = FALSE, sort.by = "p-value", filename = "TopPlots.pdf", scale_y = NULL) ## S4 method for signature 'Rnits' plotResults(object, id = NULL, fdr = NULL, top = 48, pdf = FALSE, sort.by = "p-value", filename = "TopPlots.pdf", scale_y = NULL)
plotResults(object, id = NULL, fdr = NULL, top = 48, pdf = FALSE, sort.by = "p-value", filename = "TopPlots.pdf", scale_y = NULL) ## S4 method for signature 'Rnits' plotResults(object, id = NULL, fdr = NULL, top = 48, pdf = FALSE, sort.by = "p-value", filename = "TopPlots.pdf", scale_y = NULL)
object |
|
id |
Names of specific genes or probes to be plotted. Overrides |
fdr |
FDR cut-off plotting top probes or genes. Overrides |
top |
Number of top genes or probes whose profile is to be plotted. Default |
pdf |
Save plot as pdf? Default |
sort.by |
Criteria for sorting top genes or probes. Default |
filename |
Name of pdf file. Default |
scale_y |
If 'free', use free scales for plots. Default NULL. |
... |
Optional arguments to plot |
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Plot top results plotResults(rnitsobj, top = 16) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Plot top results plotResults(rnitsobj, top = 16) ## End(Not run)
The code utilizes the probe-gene mapping from features file to summarize probe-level log ratios to gene level ratios.
summarizeProbes(object) ## S4 method for signature 'Rnits' summarizeProbes(object)
summarizeProbes(object) ## S4 method for signature 'Rnits' summarizeProbes(object)
object |
|
Tukey's biweight is used to compute gene level summary
An object of class Rnits
with gene level log ratios, which can be retrieved by getLR(object)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') # Summarize gene-level data rnitsobj <- summarizeProbes(rnitsobj)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') # Summarize gene-level data rnitsobj <- summarizeProbes(rnitsobj)
Summarize top genes or probes from Rnits fit
method
## S4 method for signature 'Rnits' summary(object, top = 48, fdr = NULL, plot = FALSE, sort.by = "p-value")
## S4 method for signature 'Rnits' summary(object, top = 48, fdr = NULL, plot = FALSE, sort.by = "p-value")
object |
|
top |
Display results for top N genes/probes. Default |
fdr |
Display results for genes/probes less than FDR (%) cutoff (if provided). Overrides |
plot |
If |
sort.by |
Sort top results by either |
A table of top genes/profiles listing the ratio statistics, p-values, q-values and cluster information.
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get summary of top genes summary(rnitsobj, FDR = 5) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) # Get summary of top genes summary(rnitsobj, FDR = 5) ## End(Not run)
Align multiple time series to the average seris
timeAlign(object, iterMax = 5, seed = 123, null.frac = 0.75, anchor = NULL, rerun = FALSE, plot = FALSE) ## S4 method for signature 'Rnits' timeAlign(object, iterMax = 5, seed = 123, null.frac = 0.75, anchor = NULL, rerun = FALSE, plot = FALSE)
timeAlign(object, iterMax = 5, seed = 123, null.frac = 0.75, anchor = NULL, rerun = FALSE, plot = FALSE) ## S4 method for signature 'Rnits' timeAlign(object, iterMax = 5, seed = 123, null.frac = 0.75, anchor = NULL, rerun = FALSE, plot = FALSE)
object |
|
iterMax |
Maximum iterations to be performed |
seed |
Random seed |
null.frac |
Fraction of genes that are considered as null |
anchor |
Sample to be considerded as base for aligning time series. If not provided, the average is used |
rerun |
If |
plot |
If |
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') # Do curve-registration on data rnitsobj <- timeAlign(rnitsobj)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') # Do curve-registration on data rnitsobj <- timeAlign(rnitsobj)
Extract expression data for top genes/probes
topData(object, id = NULL, fdr = NULL, top = 16, sort.by = "p-value") ## S4 method for signature 'Rnits' topData(object, id = NULL, fdr = NULL, top = 16, sort.by = "p-value")
topData(object, id = NULL, fdr = NULL, top = 16, sort.by = "p-value") ## S4 method for signature 'Rnits' topData(object, id = NULL, fdr = NULL, top = 16, sort.by = "p-value")
object |
|
id |
Names of probes or genes |
top |
Display results for top N genes/probes. Default |
fdr |
Display results for genes/probes less than FDR cutoff (if provided). Overrides |
sort.by |
Sort top results by either |
A table of expression values of top genes/profiles
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) #Get data for top genes td <- topData(rnitsobj, FDR = 5) ## End(Not run)
# load pre-compiled expressionSet object for Ronen and Botstein yeast chemostat data data(yeastchemostat) rnitsobj = build.Rnits(yeastchemostat, logscale = TRUE, normmethod = 'Between') ## Not run: # Fit model using gene-level summarization rnitsobj <- fit(rnitsobj, gene.level = TRUE, clusterallsamples = FALSE) #Get data for top genes td <- topData(rnitsobj, FDR = 5) ## End(Not run)
(From author's GEO submission) Transcriptional response of steady-state yeast cultures to transient perturbations in carbon source
yeastchemostat
yeastchemostat
An ExpressionSet
object with containing 'Sample' and 'Time' columns and replicates removed.
http://www.ncbi.nlm.nih.gov/geo/query/acc.cgi?acc=GSE4158