if (!require("BiocManager")) {
install.packages("BiocManager")
}
BiocManager::install("glmSparseNet")
library(dplyr)
library(ggplot2)
library(survival)
library(futile.logger)
library(curatedTCGAData)
library(TCGAutils)
library(MultiAssayExperiment)
#
library(glmSparseNet)
#
# Some general options for futile.logger the debugging package
flog.layout(layout.format("[~l] ~m"))
options(
"glmSparseNet.show_message" = FALSE,
"glmSparseNet.base_dir" = withr::local_tempdir()
)
# 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.
brca <- curatedTCGAData(
diseaseCode = "BRCA", assays = "RNASeq2GeneNorm",
version = "1.1.38", dry.run = FALSE
)
# keep only solid tumour (code: 01)
brcaPrimarySolidTumor <- TCGAutils::TCGAsplitAssays(brca, "01")
xdataRaw <- t(assay(brcaPrimarySolidTumor[[1]]))
# Get survival information
ydataRaw <- colData(brcaPrimarySolidTumor) |>
as.data.frame() |>
# Keep only data relative to survival or samples
dplyr::select(
patientID, vital_status,
Days.to.date.of.Death, Days.to.Date.of.Last.Contact,
days_to_death, days_to_last_followup,
Vital.Status
) |>
# Convert days to integer
dplyr::mutate(Days.to.date.of.Death = as.integer(Days.to.date.of.Death)) |>
dplyr::mutate(
Days.to.Last.Contact = as.integer(Days.to.Date.of.Last.Contact)
) |>
# Find max time between all days (ignoring missings)
dplyr::rowwise() |>
dplyr::mutate(
time = max(days_to_last_followup, Days.to.date.of.Death,
Days.to.Last.Contact, 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(ydataRaw) <- ydataRaw$patientID
# Get matches between survival and assay data
xdataRaw <- xdataRaw[
TCGAbarcode(rownames(xdataRaw)) %in% rownames(ydataRaw),
]
xdataRaw <- xdataRaw[, apply(xdataRaw, 2, sd) != 0] |>
scale()
# Order ydata the same as assay
ydataRaw <- ydataRaw[TCGAbarcode(rownames(xdataRaw)), ]
# Using only a subset of genes previously selected to keep this short example.
set.seed(params$seed)
smallSubset <- c(
"CD5", "CSF2RB", "IRGC", "NEUROG2", "NLRC4", "PDE11A",
"PTEN", "TP53", "BRAF",
"PIK3CB", "QARS", "RFC3", "RPGRIP1L", "SDC1", "TMEM31",
"YME1L1", "ZBTB11", sample(colnames(xdataRaw), 100)
) |>
unique()
xdata <- xdataRaw[, smallSubset[smallSubset %in% colnames(xdataRaw)]]
ydata <- ydataRaw |> dplyr::select(time, status)
Fit model model penalizing by the hubs using the cross-validation function by
cv.glmHub
.
set.seed(params$seed)
fitted <- cv.glmHub(xdata, Surv(ydata$time, ydata$status),
family = "cox",
lambda = buildLambda(1),
network = "correlation",
options = networkOptions(
cutoff = .6,
minDegree = .2
)
)
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
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
coefsCV <- Filter(function(.x) .x != 0, coef(fitted, s = "lambda.min")[, 1])
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
data.frame(
gene.name = names(coefsCV),
coefficient = coefsCV,
stringsAsFactors = FALSE
) |>
arrange(gene.name) |>
knitr::kable()
gene.name | coefficient | |
---|---|---|
CD5 | CD5 | -0.16632 |
separate2GroupsCox(as.vector(coefsCV),
xdata[, names(coefsCV)],
ydata,
plotTitle = "Full dataset", legendOutside = FALSE
)
## $pvalue
## [1] 0.001237802
##
## $plot
##
## $km
## Call: survfit(formula = survival::Surv(time, status) ~ group, data = prognosticIndexDf)
##
## n events median 0.95LCL 0.95UCL
## Low risk - 1 540 58 3959 3492 NA
## High risk - 1 540 94 3738 3262 4456
sessionInfo()
## R version 4.5.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] glmnet_4.1-10 VennDiagram_1.7.3
## [3] reshape2_1.4.4 forcats_1.0.1
## [5] Matrix_1.7-4 glmSparseNet_1.27.0
## [7] TCGAutils_1.29.5 curatedTCGAData_1.31.2
## [9] MultiAssayExperiment_1.35.9 SummarizedExperiment_1.39.2
## [11] Biobase_2.69.1 GenomicRanges_1.61.5
## [13] Seqinfo_0.99.2 IRanges_2.43.5
## [15] S4Vectors_0.47.4 BiocGenerics_0.55.1
## [17] generics_0.1.4 MatrixGenerics_1.21.0
## [19] matrixStats_1.5.0 futile.logger_1.4.3
## [21] survival_3.8-3 ggplot2_4.0.0
## [23] dplyr_1.1.4 BiocStyle_2.37.1
##
## loaded via a namespace (and not attached):
## [1] RColorBrewer_1.1-3 jsonlite_2.0.0
## [3] shape_1.4.6.1 magrittr_2.0.4
## [5] magick_2.9.0 GenomicFeatures_1.61.6
## [7] farver_2.1.2 rmarkdown_2.30
## [9] BiocIO_1.19.0 vctrs_0.6.5
## [11] memoise_2.0.1 Rsamtools_2.25.3
## [13] RCurl_1.98-1.17 rstatix_0.7.2
## [15] tinytex_0.57 htmltools_0.5.8.1
## [17] S4Arrays_1.9.1 BiocBaseUtils_1.11.2
## [19] progress_1.2.3 AnnotationHub_3.99.6
## [21] lambda.r_1.2.4 curl_7.0.0
## [23] broom_1.0.10 Formula_1.2-5
## [25] pROC_1.19.0.1 SparseArray_1.9.1
## [27] sass_0.4.10 bslib_0.9.0
## [29] plyr_1.8.9 httr2_1.2.1
## [31] zoo_1.8-14 futile.options_1.0.1
## [33] cachem_1.1.0 GenomicAlignments_1.45.4
## [35] lifecycle_1.0.4 iterators_1.0.14
## [37] pkgconfig_2.0.3 R6_2.6.1
## [39] fastmap_1.2.0 digest_0.6.37
## [41] AnnotationDbi_1.71.1 ps_1.9.1
## [43] ExperimentHub_2.99.5 RSQLite_2.4.3
## [45] ggpubr_0.6.1 filelock_1.0.3
## [47] labeling_0.4.3 km.ci_0.5-6
## [49] httr_1.4.7 abind_1.4-8
## [51] compiler_4.5.1 bit64_4.6.0-1
## [53] withr_3.0.2 S7_0.2.0
## [55] backports_1.5.0 BiocParallel_1.43.4
## [57] carData_3.0-5 DBI_1.2.3
## [59] ggsignif_0.6.4 biomaRt_2.65.16
## [61] rappdirs_0.3.3 DelayedArray_0.35.3
## [63] rjson_0.2.23 tools_4.5.1
## [65] chromote_0.5.1 glue_1.8.0
## [67] restfulr_0.0.16 promises_1.3.3
## [69] checkmate_2.3.3 gtable_0.3.6
## [71] KMsurv_0.1-6 tzdb_0.5.0
## [73] tidyr_1.3.1 survminer_0.5.1
## [75] websocket_1.4.4 data.table_1.17.8
## [77] hms_1.1.3 car_3.1-3
## [79] xml2_1.4.0 XVector_0.49.1
## [81] BiocVersion_3.22.0 foreach_1.5.2
## [83] pillar_1.11.1 stringr_1.5.2
## [85] later_1.4.4 splines_4.5.1
## [87] BiocFileCache_2.99.6 lattice_0.22-7
## [89] rtracklayer_1.69.1 bit_4.6.0
## [91] tidyselect_1.2.1 Biostrings_2.77.2
## [93] knitr_1.50 gridExtra_2.3
## [95] bookdown_0.45 xfun_0.53
## [97] stringi_1.8.7 UCSC.utils_1.5.0
## [99] yaml_2.3.10 evaluate_1.0.5
## [101] codetools_0.2-20 tibble_3.3.0
## [103] BiocManager_1.30.26 cli_3.6.5
## [105] xtable_1.8-4 processx_3.8.6
## [107] jquerylib_0.1.4 survMisc_0.5.6
## [109] dichromat_2.0-0.1 Rcpp_1.1.0
## [111] GenomeInfoDb_1.45.12 GenomicDataCommons_1.33.1
## [113] dbplyr_2.5.1 png_0.1-8
## [115] XML_3.99-0.19 readr_2.1.5
## [117] blob_1.2.4 prettyunits_1.2.0
## [119] bitops_1.0-9 scales_1.4.0
## [121] purrr_1.1.0 crayon_1.5.3
## [123] rlang_1.1.6 KEGGREST_1.49.1
## [125] rvest_1.0.5 formatR_1.14