| getModelFit {limorhyde2} | R Documentation |
Fit linear models for rhythmicity in one or more conditions
Description
This is the first step in an analysis using limorhyde2, the second is to
moderate the fits using getPosteriorFit().
Usage
getModelFit(
y,
metadata,
period = 24,
nKnots = 3L,
degree = if (nKnots > 2) 3L else 2L,
sinusoid = FALSE,
timeColname = "time",
condColname = NULL,
covarColnames = NULL,
sampleColname = "sample",
nShifts = 3L,
method = c("trend", "voom", "deseq2"),
lmFitArgs = list(),
eBayesArgs = if (method == "trend") list(trend = TRUE) else list(),
DESeqArgs = list(),
keepLmFits = FALSE
)
Arguments
y |
Matrix-like object of measurements, with rows corresponding to features and columns to samples. |
metadata |
data.frame containing experimental design information for
each sample. Rows of |
period |
Number specifying the period for the time variable, in the same
units as the values in the |
nKnots |
Number of internal knots for the periodic spline for the time variable. |
degree |
Integer indicating degree of the piecewise polynomial for the spline. |
sinusoid |
Logical indicating whether to fit a cosinor-based model instead of a spline-based model. |
timeColname |
String indicating the column in |
condColname |
String indicating the column in |
covarColnames |
Character vector indicating the columns in |
sampleColname |
String indicating the column in |
nShifts |
Number of shifted models to fit. Only used for periodic splines, not for cosinor. Do not change from the default unless you know what you're doing. |
method |
String indicating method to estimate model coefficients. For microarray data, use 'trend'. For RNA-seq count data, use 'voom' or 'deseq2'. |
lmFitArgs |
List of arguments passed to |
eBayesArgs |
List of arguments passed to |
DESeqArgs |
List of arguments passed to |
keepLmFits |
Logical indicating whether to keep the complete fit objects
from |
Value
A limorhyde2 object with elements:
-
metadata: As supplied above, converted to adata.table. -
timeColname: As supplied above. -
condColname: As supplied above. -
covarColnames: As supplied above. -
coefficients: Matrix with rows corresponding to features and columns to model terms, including all shifted models. -
shifts: Numeric vector indicating amount by which timepoints were shifted for each shifted model. -
period: As supplied above. -
conds: IfcondColnameis notNULL, a vector of unique values of the condition variable. -
nKnots: Number of knots. -
degree: As supplied above. -
sinusoid: As supplied above. -
nConds: Number of conditions. -
nCovs: Number of covariates. -
lmFits: IfkeepLmFitsisTRUE, a list of objects fromlimmaorDESeq2, with length equal to length of theshiftselement.
See Also
Examples
library('data.table')
# rhythmicity in one condition
y = GSE54650$y
metadata = GSE54650$metadata
fit = getModelFit(y, metadata)
fit = getPosteriorFit(fit)
rhyStats = getRhythmStats(fit, features = c('13170', '13869'))
# rhythmicity and differential rhythmicity in multiple conditions
y = GSE34018$y
metadata = GSE34018$metadata
fit = getModelFit(y, metadata, nKnots = 3L, condColname = 'cond')
fit = getPosteriorFit(fit)
rhyStats = getRhythmStats(fit, features = c('13170', '12686'))
diffRhyStats = getDiffRhythmStats(fit, rhyStats)