msqrobLmer: Function to fit msqrob models with ridge regression and/or...

View source: R/msqrob.R

msqrobLmerR Documentation

Function to fit msqrob models with ridge regression and/or random effects using lme4

Description

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.

Usage

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))
)

Arguments

y

A matrix with the quantified feature intensities. The features are along the rows and samples along the columns.

formula

Model formula. The model is built based on the covariates in the data object.

data

A DataFrame with information on the design. It has the same number of rows as the number of columns (samples) of y.

rowdata

A DataFrame with the rowData information of the SummarizedExperiment. It has the same number of rows as the number of rows (features) of y.

tol

numeric(1) indicating the tolerance for declaring convergence of the M-estimation loop.

robust

boolean(1) to indicate if robust regression is performed to account for outliers. Default is TRUE. If FALSE an OLS fit is performed.

ridge

boolean(1) to indicate if ridge regression is performed. Default is FALSE. If TRUE the fixed effects are estimated via penalized regression and shrunken to zero.

maxitRob

numeric(1) indicating the maximum iterations in the IRWLS algorithm used in the M-estimation step of the robust regression.

doQR

boolean(1) to indicate if QR decomposition is used when adopting ridge regression. Default is TRUE. If FALSE the predictors of the fixed effects are not transformed, and the degree of shrinkage can depend on the encoding.

featureGroups

vector of type character or vector of type factor indicating how to aggregate the features. Is only used when multiple features are used to build the model, e.g. when starting from peptide data and modelling the fold change at the protein level. The default is NULL

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 list(control = lmerControl(calc.derivs = FALSE))

Value

A list of objects of the StatModel class.

Author(s)

Lieven Clement, Oliver M. Crook

Examples


# 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]])

statOmics/msqrob2 documentation built on May 8, 2024, 4:38 p.m.