msqrobLmer | R Documentation |
Low-level function for parameter estimation with msqrob using the robust ridge regression. The models can be fitted for each feature (e.g. summarised protein expression values) or multiple features belonging to the same accession can be modelled simultaneously e.g. peptide-based models where all peptide intensities for the same protein are modelled simultaneously. The fold changes and uncertainty estimates are then calculated at the protein level while correcting for peptide species and within sample correlation.
msqrobLmer(
y,
formula,
data,
rowdata = NULL,
tol = 1e-06,
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
doQR = TRUE,
featureGroups = NULL,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE))
)
y |
A |
formula |
Model formula. The model is built based on the covariates in the data object. |
data |
A |
rowdata |
A |
tol |
|
robust |
|
ridge |
|
maxitRob |
|
doQR |
|
featureGroups |
vector of type |
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 |
A list of objects of the StatModel
class.
Lieven Clement, Oliver M. Crook
# 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 ridge regression upon summarization of
# peptide intensities into protein expression values
modelsRidge <- msqrobLmer(assay(pe[["protein"]]), ~condition, data = colData(pe),
ridge = TRUE)
getCoef(modelsRidge[[1]])
# Fit MSqrob model using robust ridge regression starting from peptide intensities
# The fold changes are calculated at the protein level while correcting for
# the different peptide species in each sample and the correlation between
# peptide intensities of peptides of the same protein in the same sample.
# Add the samples variable to colData
colData(pe)$samples <- rownames(colData(pe))
modelsPepBased <- msqrobLmer(assay(pe[["peptide"]]),
formula = ~condition + (1|samples) + (1|Sequence), data = colData(pe),
rowdata = rowData(pe[["peptide"]]), featureGroups = rowData(pe[["peptide"]])$Proteins,
ridge = TRUE)
getCoef(modelsPepBased[[1]])
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.