msqrob,SummarizedExperiment-method | R Documentation |
Parameter estimation of msqrob models for QFeatures
and SummarizedExperiment
instance.
## S4 method for signature 'SummarizedExperiment'
msqrob(
object,
formula,
modelColumnName = "msqrobModels",
overwrite = FALSE,
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
tol = 1e-06,
doQR = TRUE,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE))
)
## S4 method for signature 'QFeatures'
msqrob(
object,
i,
formula,
modelColumnName = "msqrobModels",
overwrite = FALSE,
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
tol = 1e-06,
doQR = TRUE,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE))
)
object |
|
formula |
Model formula. The model is built based on the covariates in the data object. |
modelColumnName |
|
overwrite |
|
robust |
|
ridge |
|
maxitRob |
|
tol |
|
doQR |
|
lmerArgs |
a list (of correct class, resulting from ‘lmerControl()’
containing control parameters, including the nonlinear optimizer to be used
and parameters to be passed through to the nonlinear optimizer, see the
‘lmerControl’ documentation of the lme4 package for more details.
Default is |
i |
|
A SummarizedExperiment or a QFeatures
instance with the models.
Lieven Clement
# Load example data
# The data are a Feature object with containing
# a SummarizedExperiment named "peptide" with MaxQuant peptide intensities
# The data are a subset of spike-in the human-ecoli study
# The variable condition in the colData of the Feature object
# contains information on the spike in condition a-e (from low to high)
data(pe)
# Aggregate peptide intensities in protein expression values
pe <- aggregateFeatures(pe, i = "peptide", fcol = "Proteins", name = "protein")
# Fit MSqrob model using robust linear regression upon summarization of
# peptide intensities into protein expression values.
# For summarized SummarizedExperiment
se <- pe[["protein"]]
se
colData(se) <- colData(pe)
se <- msqrob(se, formula = ~condition, modelColumnName = "rlm")
getCoef(rowData(se)$rlm[[1]])
# For features object
pe <- msqrob(pe, i = "protein", formula = ~condition, modelColumnName = "rlm")
# with ridge regression (slower)
pe <- msqrob(pe, i = "protein", formula = ~condition, ridge = TRUE, modelColumnName = "ridge")
# compare for human protein (no DE)==> large shrinkage to zero
cbind(getCoef(rowData(pe[["protein"]])$rlm[[1]]), getCoef(rowData(pe[["protein"]])$ridge[[1]]))
# compare for ecoli protein (DE)==> almost no shrinkage to zero
cbind(
getCoef(rowData(pe[["protein"]])$rlm[["P00956"]]),
getCoef(rowData(pe[["protein"]])$ridge[["P00956"]])
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.