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.4.0 RC (2024-04-16 r86468)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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-8 VennDiagram_1.7.3
## [3] reshape2_1.4.4 forcats_1.0.0
## [5] Matrix_1.7-0 glmSparseNet_1.23.0
## [7] TCGAutils_1.25.0 curatedTCGAData_1.25.4
## [9] MultiAssayExperiment_1.31.0 SummarizedExperiment_1.35.0
## [11] Biobase_2.65.0 GenomicRanges_1.57.0
## [13] GenomeInfoDb_1.41.0 IRanges_2.39.0
## [15] S4Vectors_0.43.0 BiocGenerics_0.51.0
## [17] MatrixGenerics_1.17.0 matrixStats_1.3.0
## [19] futile.logger_1.4.3 survival_3.6-4
## [21] ggplot2_3.5.1 dplyr_1.1.4
## [23] BiocStyle_2.33.0
##
## loaded via a namespace (and not attached):
## [1] jsonlite_1.8.8 shape_1.4.6.1
## [3] magrittr_2.0.3 magick_2.8.3
## [5] GenomicFeatures_1.57.0 farver_2.1.1
## [7] rmarkdown_2.26 BiocIO_1.15.0
## [9] zlibbioc_1.51.0 vctrs_0.6.5
## [11] memoise_2.0.1 Rsamtools_2.21.0
## [13] RCurl_1.98-1.14 rstatix_0.7.2
## [15] tinytex_0.50 progress_1.2.3
## [17] htmltools_0.5.8.1 S4Arrays_1.5.0
## [19] BiocBaseUtils_1.7.0 AnnotationHub_3.13.0
## [21] lambda.r_1.2.4 curl_5.2.1
## [23] broom_1.0.5 pROC_1.18.5
## [25] SparseArray_1.5.0 sass_0.4.9
## [27] bslib_0.7.0 plyr_1.8.9
## [29] httr2_1.0.1 zoo_1.8-12
## [31] futile.options_1.0.1 cachem_1.0.8
## [33] GenomicAlignments_1.41.0 mime_0.12
## [35] lifecycle_1.0.4 iterators_1.0.14
## [37] pkgconfig_2.0.3 R6_2.5.1
## [39] fastmap_1.1.1 GenomeInfoDbData_1.2.12
## [41] digest_0.6.35 colorspace_2.1-0
## [43] AnnotationDbi_1.67.0 ps_1.7.6
## [45] ExperimentHub_2.13.0 RSQLite_2.3.6
## [47] ggpubr_0.6.0 labeling_0.4.3
## [49] filelock_1.0.3 km.ci_0.5-6
## [51] fansi_1.0.6 httr_1.4.7
## [53] abind_1.4-5 compiler_4.4.0
## [55] bit64_4.0.5 withr_3.0.0
## [57] backports_1.4.1 BiocParallel_1.39.0
## [59] carData_3.0-5 DBI_1.2.2
## [61] highr_0.10 ggsignif_0.6.4
## [63] biomaRt_2.61.0 rappdirs_0.3.3
## [65] DelayedArray_0.31.0 rjson_0.2.21
## [67] tools_4.4.0 chromote_0.2.0
## [69] glue_1.7.0 restfulr_0.0.15
## [71] promises_1.3.0 checkmate_2.3.1
## [73] generics_0.1.3 gtable_0.3.5
## [75] KMsurv_0.1-5 tzdb_0.4.0
## [77] tidyr_1.3.1 survminer_0.4.9
## [79] websocket_1.4.1 data.table_1.15.4
## [81] hms_1.1.3 car_3.1-2
## [83] xml2_1.3.6 utf8_1.2.4
## [85] XVector_0.45.0 BiocVersion_3.20.0
## [87] foreach_1.5.2 pillar_1.9.0
## [89] stringr_1.5.1 later_1.3.2
## [91] splines_4.4.0 BiocFileCache_2.13.0
## [93] lattice_0.22-6 rtracklayer_1.65.0
## [95] bit_4.0.5 tidyselect_1.2.1
## [97] Biostrings_2.73.0 knitr_1.46
## [99] gridExtra_2.3 bookdown_0.39
## [101] xfun_0.43 stringi_1.8.3
## [103] UCSC.utils_1.1.0 yaml_2.3.8
## [105] evaluate_0.23 codetools_0.2-20
## [107] tibble_3.2.1 BiocManager_1.30.22
## [109] cli_3.6.2 xtable_1.8-4
## [111] munsell_0.5.1 processx_3.8.4
## [113] jquerylib_0.1.4 survMisc_0.5.6
## [115] Rcpp_1.0.12 GenomicDataCommons_1.29.0
## [117] dbplyr_2.5.0 png_0.1-8
## [119] XML_3.99-0.16.1 readr_2.1.5
## [121] blob_1.2.4 prettyunits_1.2.0
## [123] bitops_1.0-7 scales_1.3.0
## [125] purrr_1.0.2 crayon_1.5.2
## [127] rlang_1.1.3 KEGGREST_1.45.0
## [129] rvest_1.0.4 formatR_1.14