msqrobHurdle,SummarizedExperiment-method | R Documentation |
Fitting a hurdle msqrob model with an intensity component for assessing differential abundance and an count component that is modeling peptide counts using quasibinomial glm for differential detection of the number of features that were not missing and used for aggregation
## S4 method for signature 'SummarizedExperiment'
msqrobHurdle(
object,
formula,
modelColumnName = "msqrobHurdle",
overwrite = FALSE,
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
tol = 1e-06,
doQR = TRUE,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE)),
priorCount = 0.1,
binomialBound = TRUE
)
## S4 method for signature 'QFeatures'
msqrobHurdle(
object,
i,
formula,
modelColumnName = "msqrobHurdle",
overwrite = FALSE,
robust = TRUE,
ridge = FALSE,
maxitRob = 1,
tol = 1e-06,
doQR = TRUE,
lmerArgs = list(control = lmerControl(calc.derivs = FALSE)),
priorCount = 0.1,
binomialBound = TRUE
)
object |
|
formula |
Model formula. Both model components are 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 |
priorCount |
A 'numeric(1)', which is a prior count to be added to the observations to shrink the estimated odds ratios of the count component towards zero. Default is 0.1. |
binomialBound |
logical, if ‘TRUE’ then the quasibinomial variance estimator will be never smaller than 1 (no underdispersion). Default is TRUE. |
i |
|
SummarizedExperiment or QFeatures instance
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 to protein expression values
pe <- aggregateFeatures(pe, i = "peptide", fcol = "Proteins", name = "protein")
# Fit Hurdle MSqrob model
# For summarized SummarizedExperiment
se <- pe[["protein"]]
se
colData(se) <- colData(pe)
se <- msqrobHurdle(se, formula = ~condition)
getCoef(rowData(se)$msqrobHurdleIntensity[[1]])
getCoef(rowData(se)$msqrobHurdleCount[[1]])
# For features object
pe <- msqrobHurdle(pe, i = "protein", formula = ~condition)
getCoef(rowData(pe[["protein"]])$msqrobHurdleIntensity[[1]])
getCoef(rowData(pe[["protein"]])$msqrobHurdleCount[[1]])
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.