msqrobAggregate,SummarizedExperiment-method | R Documentation |
Parameter estimation of msqrob models for QFeatures
instance.
The method aggregates features within the model e.g. from peptides to proteins.
It provides fold change estimates and their associated uncertainty at the aggregated
level (e.g. protein level) while correcting for the peptide species that are observed
in each sample. It also addresses the correlation in the data, e.g. the peptide data
for the same protein in a sample are correlate because they originate from the same
protein pool. The method however does not return aggregated expression values for each sample.
For visualisation purposes aggregated expression values are provide by the aggregateFeatures
function from the QFeatures
Package
## S4 method for signature 'SummarizedExperiment'
msqrobAggregate(
object,
formula,
fcol,
aggregateFun = MsCoreUtils::robustSummary,
modelColumnName = "msqrobModels",
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
tol = 1e-06,
doQR = TRUE,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE))
)
## S4 method for signature 'QFeatures'
msqrobAggregate(
object,
formula,
i,
fcol,
name = "msqrobAggregate",
aggregateFun = MsCoreUtils::robustSummary,
modelColumnName = "msqrobModels",
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. |
fcol |
The feature variable of assay ‘i’ defining how to summarise the features. |
aggregateFun |
A function used for quantitative feature aggregation.
Details can be found in the documentation of the |
modelColumnName |
|
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 |
|
name |
A ‘character(1)’ naming the new assay. Default is ‘newAssay’. Note that the function will fail if there's already an assay with ‘name’. |
A ‘QFeatures’ object with an additional assay.
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)
# 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.
colData(pe)$samples <- rownames(colData(pe))
pe <- msqrobAggregate(pe, i = "peptide", fcol = "Proteins",
formula = ~condition + (1|samples) + (1|Sequence),
ridge = TRUE)
getCoef(rowData(pe[["msqrobAggregate"]])$msqrobModels[["P00956"]])
## Same but on a SummarizedExperiment object
se <- getWithColData(pe, "peptide")
se <- msqrobAggregate(se, fcol = "Proteins",
formula = ~condition + (1|samples) + (1|Sequence),
ridge = TRUE)
getCoef(rowData(se)$msqrobModels[["P00956"]])
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.