if (!require('BiocManager'))
install.packages('BiocManager')
BiocManager::install('glmSparseNet')
library(dplyr)
library(ggplot2)
library(survival)
library(futile.logger)
library(curatedTCGAData)
library(TCGAutils)
#
library(glmSparseNet)
#
# Some general options for futile.logger the debugging package
.Last.value <- flog.layout(layout.format('[~l] ~m'))
.Last.value <- glmSparseNet:::show.message(FALSE)
# Setting ggplot2 default theme as minimal
theme_set(ggplot2::theme_minimal())
The data is loaded from an online curated dataset downloaded from TCGA using
curatedTCGAData
bioconductor package and processed.
To accelerate the process we use a very reduced dataset down to 107 variables only (genes), which is stored as a data object in this package. However, the procedure to obtain the data manually is described in the following chunk.
skcm <- curatedTCGAData(diseaseCode = 'SKCM', assays = 'RNASeq2GeneNorm',
version = '1.1.38', dry.run = FALSE)
Build the survival data from the clinical columns.
xdata
and ydata
skcm.metastatic <- TCGAutils::TCGAsplitAssays(skcm, '06')
xdata.raw <- t(assay(skcm.metastatic[[1]]))
# Get survival information
ydata.raw <- colData(skcm.metastatic) %>% as.data.frame %>%
# Find max time between all days (ignoring missings)
dplyr::rowwise() %>%
dplyr::mutate(
time = max(days_to_last_followup,
days_to_death,
na.rm = TRUE)
) %>%
# Keep only survival variables and codes
dplyr::select(patientID, status = vital_status, time) %>%
# Discard individuals with survival time less or equal to 0
dplyr::filter(!is.na(time) & time > 0) %>%
as.data.frame()
# Get survival information
ydata.raw <- colData(skcm) %>% as.data.frame %>%
# Find max time between all days (ignoring missings)
dplyr::rowwise() %>%
dplyr::mutate(
time = max(days_to_last_followup, days_to_death, na.rm = TRUE)
) %>%
# Keep only survival variables and codes
dplyr::select(patientID, status = vital_status, time) %>%
# Discard individuals with survival time less or equal to 0
dplyr::filter(!is.na(time) & time > 0) %>% as.data.frame
# Set index as the patientID
rownames(ydata.raw) <- ydata.raw$patientID
# keep only features that have standard deviation > 0
xdata.raw <- xdata.raw[TCGAbarcode(rownames(xdata.raw)) %in%
rownames(ydata.raw),]
xdata.raw <- xdata.raw %>%
{ (apply(., 2, sd) != 0) } %>%
{ xdata.raw[, .] } %>%
scale
# Order ydata the same as assay
ydata.raw <- ydata.raw[TCGAbarcode(rownames(xdata.raw)), ]
set.seed(params$seed)
small.subset <- c('FOXL2', 'KLHL5', 'PCYT2', 'SLC6A10P', 'STRAP', 'TMEM33',
'WT1-AS', sample(colnames(xdata.raw), 100))
xdata <- xdata.raw[, small.subset[small.subset %in% colnames(xdata.raw)]]
ydata <- ydata.raw %>% dplyr::select(time, status)
Fit model model penalizing by the hubs using the cross-validation function by
cv.glmHub
.
fitted <- cv.glmHub(
xdata,
Surv(ydata$time, ydata$status),
family = 'cox',
foldid = glmSparseNet:::balanced.cv.folds(!!ydata$status)$output,
network = 'correlation',
network.options = networkOptions(min.degree = .2,
cutoff = .6)
)
Shows the results of 100
different parameters used to find the optimal value
in 10-fold cross-validation. The two vertical dotted lines represent the best
model and a model with less variables selected (genes), but within a standard
error distance from the best.
plot(fitted)
Taking the best model described by lambda.min
coefs.v <- coef(fitted, s = 'lambda.min')[,1] %>% { .[. != 0]}
coefs.v %>% {
data.frame(ensembl.id = names(.),
gene.name = geneNames(names(.))$external_gene_name,
coefficient = .,
stringsAsFactors = FALSE)
} %>%
arrange(gene.name) %>%
knitr::kable()
ensembl.id | gene.name | coefficient | |
---|---|---|---|
PCYT2 | PCYT2 | AMICA1 | 0.0646641 |
AMICA1 | AMICA1 | C4orf49 | -0.2758400 |
C4orf49 | C4orf49 | PCYT2 | -0.0059089 |
separate2GroupsCox(as.vector(coefs.v),
xdata[, names(coefs.v)],
ydata,
plot.title = 'Full dataset', legend.outside = FALSE)
## $pvalue
## [1] 0.0001269853
##
## $plot
##
## $km
## Call: survfit(formula = survival::Surv(time, status) ~ group, data = prognostic.index.df)
##
## n events median 0.95LCL 0.95UCL
## Low risk 180 79 4000 2927 6164
## High risk 179 114 2005 1524 2829
sessionInfo()
## R version 4.3.2 Patched (2023-11-01 r85457)
## Platform: x86_64-apple-darwin20 (64-bit)
## Running under: macOS Monterey 12.7.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRlapack.dylib; LAPACK version 3.11.0
##
## locale:
## [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
##
## time zone: America/New_York
## tzcode source: internal
##
## attached base packages:
## [1] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] VennDiagram_1.7.3 reshape2_1.4.4
## [3] forcats_1.0.0 glmSparseNet_1.20.1
## [5] glmnet_4.1-8 Matrix_1.6-5
## [7] TCGAutils_1.22.2 curatedTCGAData_1.24.0
## [9] MultiAssayExperiment_1.28.0 SummarizedExperiment_1.32.0
## [11] Biobase_2.62.0 GenomicRanges_1.54.1
## [13] GenomeInfoDb_1.38.5 IRanges_2.36.0
## [15] S4Vectors_0.40.2 BiocGenerics_0.48.1
## [17] MatrixGenerics_1.14.0 matrixStats_1.2.0
## [19] futile.logger_1.4.3 survival_3.5-7
## [21] ggplot2_3.4.4 dplyr_1.1.4
## [23] BiocStyle_2.30.0
##
## loaded via a namespace (and not attached):
## [1] jsonlite_1.8.8 shape_1.4.6
## [3] magrittr_2.0.3 magick_2.8.2
## [5] GenomicFeatures_1.54.3 farver_2.1.1
## [7] rmarkdown_2.25 BiocIO_1.12.0
## [9] zlibbioc_1.48.0 vctrs_0.6.5
## [11] memoise_2.0.1 Rsamtools_2.18.0
## [13] RCurl_1.98-1.14 rstatix_0.7.2
## [15] BiocBaseUtils_1.4.0 htmltools_0.5.7
## [17] S4Arrays_1.2.0 progress_1.2.3
## [19] AnnotationHub_3.10.0 lambda.r_1.2.4
## [21] curl_5.2.0 broom_1.0.5
## [23] SparseArray_1.2.3 sass_0.4.8
## [25] bslib_0.6.1 plyr_1.8.9
## [27] zoo_1.8-12 futile.options_1.0.1
## [29] cachem_1.0.8 GenomicAlignments_1.38.2
## [31] mime_0.12 lifecycle_1.0.4
## [33] iterators_1.0.14 pkgconfig_2.0.3
## [35] R6_2.5.1 fastmap_1.1.1
## [37] GenomeInfoDbData_1.2.11 shiny_1.8.0
## [39] digest_0.6.34 colorspace_2.1-0
## [41] AnnotationDbi_1.64.1 ExperimentHub_2.10.0
## [43] RSQLite_2.3.5 ggpubr_0.6.0
## [45] filelock_1.0.3 labeling_0.4.3
## [47] km.ci_0.5-6 fansi_1.0.6
## [49] httr_1.4.7 abind_1.4-5
## [51] compiler_4.3.2 bit64_4.0.5
## [53] withr_3.0.0 backports_1.4.1
## [55] BiocParallel_1.36.0 carData_3.0-5
## [57] DBI_1.2.1 highr_0.10
## [59] ggsignif_0.6.4 biomaRt_2.58.2
## [61] rappdirs_0.3.3 DelayedArray_0.28.0
## [63] rjson_0.2.21 tools_4.3.2
## [65] interactiveDisplayBase_1.40.0 httpuv_1.6.14
## [67] glue_1.7.0 restfulr_0.0.15
## [69] promises_1.2.1 generics_0.1.3
## [71] gtable_0.3.4 KMsurv_0.1-5
## [73] tzdb_0.4.0 tidyr_1.3.1
## [75] survminer_0.4.9 data.table_1.15.0
## [77] hms_1.1.3 car_3.1-2
## [79] xml2_1.3.6 utf8_1.2.4
## [81] XVector_0.42.0 BiocVersion_3.18.1
## [83] foreach_1.5.2 pillar_1.9.0
## [85] stringr_1.5.1 later_1.3.2
## [87] splines_4.3.2 BiocFileCache_2.10.1
## [89] lattice_0.22-5 rtracklayer_1.62.0
## [91] bit_4.0.5 tidyselect_1.2.0
## [93] Biostrings_2.70.2 knitr_1.45
## [95] gridExtra_2.3 bookdown_0.37
## [97] xfun_0.41 stringi_1.8.3
## [99] yaml_2.3.8 evaluate_0.23
## [101] codetools_0.2-19 tibble_3.2.1
## [103] BiocManager_1.30.22 cli_3.6.2
## [105] xtable_1.8-4 munsell_0.5.0
## [107] jquerylib_0.1.4 survMisc_0.5.6
## [109] Rcpp_1.0.12 GenomicDataCommons_1.26.0
## [111] dbplyr_2.4.0 png_0.1-8
## [113] XML_3.99-0.16.1 ellipsis_0.3.2
## [115] readr_2.1.5 blob_1.2.4
## [117] prettyunits_1.2.0 bitops_1.0-7
## [119] scales_1.3.0 purrr_1.0.2
## [121] crayon_1.5.2 rlang_1.1.3
## [123] KEGGREST_1.42.0 rvest_1.0.3
## [125] formatR_1.14