eval_Coxmos_models {Coxmos} | R Documentation |
eval_Coxmos_models
Description
The eval_Coxmos_models
function facilitates the comprehensive evaluation of multiple Coxmos
models in a concurrent manner. It is designed to provide a detailed assessment of the models'
performance by calculating the Area Under the Curve (AUC) for each model at specified time points.
The results generated by this function are primed for visualization using the plot_evaluation()
function.
Usage
eval_Coxmos_models(
lst_models,
X_test,
Y_test,
pred.method = "cenROC",
pred.attr = "mean",
times = NULL,
PARALLEL = FALSE,
max_time_points = 15,
verbose = FALSE,
progress_bar = TRUE
)
Arguments
lst_models |
List of Coxmos models. Each object of the list must be named. |
X_test |
Numeric matrix or data.frame. Explanatory variables for test data (raw format). Qualitative variables must be transform into binary variables. |
Y_test |
Numeric matrix or data.frame. Response variables for test data. Object must have two columns named as "time" and "event". For event column, accepted values are: 0/1 or FALSE/TRUE for censored and event observations. |
pred.method |
Character. AUC evaluation algorithm method for evaluate the model performance. Must be one of the following: "risksetROC", "survivalROC", "cenROC", "nsROC", "smoothROCtime_C", "smoothROCtime_I" (default: "cenROC"). |
pred.attr |
Character. Way to evaluate the metric selected. Must be one of the following: "mean" or "median" (default: "mean"). |
times |
Numeric vector. Time points where the AUC will be evaluated. If NULL, a maximum of 'max_time_points' points will be selected equally distributed (default: NULL). |
PARALLEL |
Logical. Run the cross validation with multicore option. As many cores as your total cores - 1 will be used. It could lead to higher RAM consumption (default: FALSE). |
max_time_points |
Numeric. Maximum number of time points to use for evaluating the model (default: 15). |
verbose |
Logical. If verbose = TRUE, extra messages could be displayed (default: FALSE). |
progress_bar |
Logical. If progress_bar = TRUE, progress bar is shown (default = TRUE). |
Details
The function begins by validating the names of the models provided in the lst_models
list and
ensures that there are at least two events present in the dataset. It then checks for the
availability of the specified evaluation method and ensures that the test times are consistent
with the training times of the models.
The core of the function revolves around the evaluation of each model. Depending on the user's
preference, the evaluations can be executed in parallel, which can significantly expedite the
process, especially when dealing with a large number of models. The function employs various
evaluation methods, as specified by the pred.method
parameter, to compute the AUC values. These
methods include but are not limited to "risksetROC", "survivalROC", and "cenROC".
Post-evaluation, the function collates the results, including training times, AIC values, c-index, Brier scores, and AUC values for each time point. The results are then transformed into a structured data frame, making it conducive for further analysis and visualization. It's worth noting that potential issues in AUC computation, often arising from sparse samples, are flagged to the user for further inspection.
Value
A list of four objects.
df
: A data.frame which the global predictions for all models. This data.frame is used to
plot the information by the function plot_evaluation()
.
lst_AUC
: A list of models where the user can check the linear predictors computed, the
global AUC, the AUC per time point and the predicted time points selected.
lst_BRIER
: A list of models where the user can check the predicted time points selected,
the Brier Score per time point and the Integrative Brier score (computed by survcomp::sbrier.score2proba
).
time
: Time used for evaluation process.
Author(s)
Pedro Salguero Garcia. Maintainer: pedsalga@upv.edu.es
References
Harrell FE, Califf RM, Pryor DB, Lee KL, Rosati RA (1982). “Evaluating the Yield of Medical Tests.” JAMA, 247. doi:10.1001/jama.1982.03320430047030, https://jamanetwork.com/journals/jama. MS S, AC C, J Q, B H (2011). “survcomp: an R/Bioconductor package for performance assessment and comparison of survival models.” Bioinformatics, 27(22), 3206-3208. Heagerty PJ, Lumley T, Pepe MS (2000). “Time-Dependent ROC Curves for Censored Survival Data and a Diagnostic Marker.” Biometrics. Heagerty PJ, Zheng Y (2005). “Survival Model Predictive Accuracy and ROC Curves.” Biometrics, 61, 92-105. doi:10.1111/j.0006-341x.2005.030814.x. Beyene KM, Ghouch AE (2020). “Smoothed time-dependent receiver operating characteristic curve for right censored survival data.” Statistics in Medicine, 39(24), 3373-3396. ISSN 10970258, https://pubmed.ncbi.nlm.nih.gov/32687225/. Pérez-Fernández S, Martínez-Camblor P, Filzmoser P, Corral N (2018). “nsROC: An R package for Non-Standard ROC Curve Analysis.” The R Journal. doi:10.1007/s00180-020-00955-7. Díaz-Coto S, Martínez-Camblor P, Pérez-Fernández S (2020). “smoothROCtime: an R package for time-dependent ROC curve estimation.” Computational Statistics, 35(3), 1231-1251. ISSN 16139658, doi:10.1007/s00180-020-00955-7.
Examples
data("X_proteomic")
data("Y_proteomic")
set.seed(123)
index_train <- caret::createDataPartition(Y_proteomic$event, p = .5, list = FALSE, times = 1)
X_train <- X_proteomic[index_train,1:50]
Y_train <- Y_proteomic[index_train,]
X_test <- X_proteomic[-index_train,1:50]
Y_test <- Y_proteomic[-index_train,]
model_icox <- splsicox(X_train, Y_train, n.comp = 2)
model_drcox <- splsdrcox(X_train, Y_train, n.comp = 2)
lst_models <- list("splsicox" = model_icox, "splsdrcox" = model_drcox)
eval_Coxmos_models(lst_models, X_test, Y_test, pred.method = "cenROC")