View source: R/testDA_censoredGLMM.R
testDA_censoredGLMM | R Documentation |
Calculate tests for differential abundance of cell populations using method 'censcyt-DA-censored-GLMM'
testDA_censoredGLMM( d_counts, formula, contrast, mi_reps = 10, imputation_method = c("km", "km_exp", "km_wei", "km_os", "rs", "mrl", "cc", "pmm"), min_cells = 3, min_samples = NULL, normalize = FALSE, norm_factors = "TMM", BPPARAM = BiocParallel::SerialParam(), verbose = FALSE )
d_counts |
|
formula |
Model formula object, see |
contrast |
Contrast matrix, created with |
mi_reps |
number of imputations in multiple imputation. default = 10. |
imputation_method |
which method should be used in the imputation step. One of 'km','km_exp','km_wei','km_os', 'rs', 'mrl', 'cc', 'pmm'. See details. default = 'km'. |
min_cells |
Filtering parameter. Default = 3. Clusters are kept for differential
testing if they have at least |
min_samples |
Filtering parameter. Default = |
normalize |
Whether to include optional normalization factors to adjust for composition effects (see details). Default = FALSE. |
norm_factors |
Normalization factors to use, if |
BPPARAM |
specify parallelization option as one of
|
verbose |
Logical. |
Calculates tests for differential abundance of clusters, using generalized linear mixed models (GLMMs) where a covariate is subject to right censoring.
The same underlying testing as described in testDA_GLMM
is
applied here. The main difference is that multiple imputation is used to
handle a censored covariate. In short, multiple imputation consists of three
steps: imputation, analysis and pooling. In the imputation step multiple complete
data sets are generated by imputation. The imputed data is then analysed in
the second step and the results are combined in the third step. See also pool
.
The imputation in the first step is specific for censored data in contrast to
the 'normal' use of multiple imputation where data is missing.
Alternatively the samples with censored data can be removed (complete case analysis)
or the censored values can be treated as missing (predictive mean matching).
Possible imputation methods in argument 'imputation_method' are:
Kaplan Meier imputation is similar to 'rs' (Risk set imputation) but the random draw is according to the survival function of the respective risk set. The largest value is treated as observed to obtain a complete survival function. (Taylor et al. 2002)
The same as 'km' but if the largest value is censored the tail of the survival function is modeled as an exponential distribution where the rate parameter is obtained by fixing the distribution to the last observed value. See (Moeschberger and Klein, 1985).
The same as 'km' but if the largest value is censored the tail of the survival function is modeled as an weibull distribution where the parameters are obtained by MLE fitting on the whole data. See (Moeschberger and Klein, 1985).
The same as 'km' but if the largest value is censored the tail of the survival function is modeled by order statistics. See (Moeschberger and Klein, 1985).
Risk Set imputation replaces the censored values with a random draw from the risk set of the respective censored value. (Taylor et al. 2002)
Mean Residual Life (Conditional multiple imputation, See Atem et al. 2017) is a multiple imputation procedure that bootstraps the data and imputes the censored values by replacing them with their respective mean residual life.
complete case (listwise deletion) analysis removes incomlete samples.
predictive mean matching treats censored values as missing and
uses predictive mean matching from mice
.
Returns a new SummarizedExperiment
object, with differential test
results stored in the rowData
slot. Results include raw p-values
(p_val
) and adjusted p-values (p_adj
), which can be used to rank
clusters by evidence for differential abundance. The results can be accessed with the
rowData
accessor function.
A Comparison of Several Methods of Estimating the Survival Function When There is Extreme Right Censoring (M. L. Moeschberger and John P. Klein, 1985)
Improved conditional imputation for linear regression with a randomly censored predictor (Atem et al. 2017)
Survival estimation and testing via multiple imputation (Taylor et al. 2002)
# create small data set with 2 differential clusters with 10 samples. d_counts <- simulate_multicluster(alphas = runif(10,1e4,1e5), sizes = runif(10,1e4,1e5), nr_diff = 2, group=2, return_summarized_experiment = TRUE)$counts # extract covariates data.frame experiment_info <- SummarizedExperiment::colData(d_counts) # add censoring experiment_info$status <- sample(c(0,1),size=10,replace = TRUE,prob = c(0.3,0.7)) experiment_info$covariate[experiment_info$status == 0] <- runif(10-sum(experiment_info$status), min=0, max=experiment_info$covariate[experiment_info$status == 0]) # create model formula object da_formula <- createFormula(experiment_info, cols_fixed = c("covariate", "group_covariate"), cols_random = "sample",event_indicator = "status") # create contrast matrix contrast <- diffcyt::createContrast(c(0, 1, 0)) # run testing with imputation method 'km' outs <- testDA_censoredGLMM(d_counts = d_counts, formula = da_formula, contrast = contrast, mi_reps = 2, imputation_method = "km") diffcyt::topTable(outs) # differential clusters: which(!is.na(SummarizedExperiment::rowData(d_counts)$paired))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.