Limma-voom contrast analysis with LOD imputation facade
Source:R/ContrastsFacades.R
ContrastsLimmaVoomImputeFacade.RdLimma-voom contrast analysis with LOD imputation facade
Limma-voom contrast analysis with LOD imputation facade
Details
Encapsulates the pipeline: strategy_limma ->
build_model_limma_voom_impute -> ContrastsLimma.
Combines vooma precision weights with LOD imputation for proteins whose fit produces NA coefficients (typically from entire missing groups).
See also
Other modelling:
AnovaExtractor,
Contrasts,
ContrastsDEqMSFacade,
ContrastsDEqMSVoomFacade,
ContrastsFirth,
ContrastsFirthFacade,
ContrastsLMFacade,
ContrastsLMImputeFacade,
ContrastsLMMissingFacade,
ContrastsLimma,
ContrastsLimmaFacade,
ContrastsLimmaImputeFacade,
ContrastsLimmaVoomFacade,
ContrastsLimpaFacade,
ContrastsLmerFacade,
ContrastsMissing,
ContrastsModerated,
ContrastsModeratedDEqMS,
ContrastsPlotter,
ContrastsRLMFacade,
ContrastsROPECA,
ContrastsROPECAFacade,
ContrastsTable,
INTERNAL_FUNCTIONS_BY_FAMILY,
LR_test(),
Model,
ModelFirth,
ModelLimma,
StrategyLM,
StrategyLimma,
StrategyLimpa,
StrategyLmer,
StrategyLogistf,
StrategyRLM,
build_contrast_analysis(),
build_model(),
build_model_glm_peptide(),
build_model_glm_protein(),
build_model_impute(),
build_model_limma(),
build_model_limma_impute(),
build_model_limma_voom(),
build_model_limma_voom_impute(),
build_model_limpa(),
build_model_logistf(),
compute_borrowed_variance(),
compute_borrowed_variance_limma(),
compute_contrast(),
compute_lmer_contrast(),
contrasts_fisher_exact(),
get_anova_df(),
get_complete_model_fit(),
get_p_values_pbeta(),
group_label(),
impute_refit_singular(),
isSingular_lm(),
linfct_all_possible_contrasts(),
linfct_factors_contrasts(),
linfct_from_model(),
linfct_matrix_contrasts(),
merge_contrasts_results(),
model_analyse(),
model_summary(),
moderated_p_deqms(),
moderated_p_deqms_long(),
moderated_p_limma(),
moderated_p_limma_long(),
new_lm_imputed(),
pivot_model_contrasts_2_Wide(),
plot_lmer_peptide_predictions(),
sim_build_models_lm(),
sim_build_models_lmer(),
sim_build_models_logistf(),
sim_make_model_lm(),
sim_make_model_lmer(),
strategy_limma(),
strategy_limpa(),
strategy_logistf(),
summary_ROPECA_median_p.scaled()
Public fields
modelModelLimma object (with imputed proteins)
contrastContrastsLimma object
.lfqdatastored reference to input LFQData
.contrast_namesnames of the requested contrasts
Methods
Method new()
initialize
Usage
ContrastsLimmaVoomImputeFacade$new(
lfqdata,
modelstr,
contrasts,
lod = NULL,
df_method = c("observed", "borrowed"),
weights = lfqdata$config$nr_children,
span = 0.5,
plot = FALSE,
...
)Arguments
lfqdataLFQData object (aggregated to protein level)
modelstrmodel formula string (e.g. "~ group_")
contrastsnamed character vector of contrasts
lodnumeric limit of detection; if NULL, auto-computed from data
df_method"observed" uses max(n_observed - p, 1); "borrowed" uses median df from successful fits
weightscolumn name for per-observation weights (default:
lfqdata$config$nr_children). PassNULLfor unweighted.spanlowess smoother span for vooma trend (default 0.5)
plotlogical; if TRUE, plot the mean-variance trend
...passed to
strategy_limma(e.g. trend, robust)
Examples
istar <- sim_lfq_data_protein_config(Nprot = 30, weight_missing = 0.5)
#> creating sampleName from file_name column
#> completing cases
#> completing cases done
#> setup done
lfqdata <- LFQData$new(istar$data, istar$config)
lfqdata$rename_response("transformedIntensity")
contrasts <- c("A_vs_Ctrl" = "group_A - group_Ctrl")
fa <- ContrastsLimmaVoomImputeFacade$new(lfqdata, "~ group_", contrasts)
#> Warning: Partial NA coefficients for 4 probe(s)
#> completing cases
head(fa$get_contrasts())
#> # A tibble: 6 × 14
#> facade modelName protein_Id contrast diff FDR std.error statistic p.value
#> <chr> <chr> <chr> <chr> <dbl> <dbl> <dbl> <dbl> <dbl>
#> 1 limma… limma 0EfVhX~29… A_vs_Ct… 1.24 0.325 0.617 2.00 0.0650
#> 2 limma… limma 0m5WN4~67… A_vs_Ct… -0.0361 0.973 1.03 -0.0352 0.973
#> 3 limma… limma 7QuTub~61… A_vs_Ct… 0.909 0.694 0.813 1.12 0.301
#> 4 limma… limma 7cbcrd~26… A_vs_Ct… 0.612 0.855 0.968 0.633 0.541
#> 5 limma… limma 9VUkAq~34… A_vs_Ct… 0.768 0.855 1.15 0.671 0.514
#> 6 limma… limma At886V~77… A_vs_Ct… -1.86 0.240 0.822 -2.26 0.0400
#> # ℹ 5 more variables: sigma <dbl>, df <dbl>, conf.low <dbl>, conf.high <dbl>,
#> # avgAbd <dbl>
fa$to_wide()
#> # A tibble: 30 × 5
#> protein_Id diff.A_vs_Ctrl p.value.A_vs_Ctrl FDR.A_vs_Ctrl statistic.A_vs_Ctrl
#> <chr> <dbl> <dbl> <dbl> <dbl>
#> 1 0EfVhX~29… 1.24 0.0650 0.325 2.00
#> 2 0m5WN4~67… -0.0361 0.973 0.973 -0.0352
#> 3 7QuTub~61… 0.909 0.301 0.694 1.12
#> 4 7cbcrd~26… 0.612 0.541 0.855 0.633
#> 5 9VUkAq~34… 0.768 0.514 0.855 0.671
#> 6 At886V~77… -1.86 0.0400 0.240 -2.26
#> 7 BEJI92~27… -1.30 0.154 0.462 -1.51
#> 8 CGzoYe~08… 0.196 0.849 0.951 0.194
#> 9 CtOJ9t~91… -0.0717 0.931 0.973 -0.0882
#> 10 DoWup2~28… -1.82 0.00796 0.0865 -3.09
#> # ℹ 20 more rows