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.
brca <- curatedTCGAData(diseaseCode = "BRCA", assays = "RNASeq2GeneNorm",
version = "1.1.38", dry.run = FALSE)
# keep only solid tumour (code: 01)
brca.primary.solid.tumor <- TCGAutils::TCGAsplitAssays(brca, '01')
xdata.raw <- t(assay(brca.primary.solid.tumor[[1]]))
# Get survival information
ydata.raw <- colData(brca.primary.solid.tumor) %>% 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(ydata.raw) <- ydata.raw$patientID
# Get matches between survival and assay data
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)), ]
# Using only a subset of genes previously selected to keep this short example.
set.seed(params$seed)
small.subset <- c('CD5', 'CSF2RB', 'IRGC', 'NEUROG2', 'NLRC4', 'PDE11A',
'PTEN', 'TP53', 'BRAF',
'PIK3CB', 'QARS', 'RFC3', 'RPGRIP1L', 'SDC1', 'TMEM31',
'YME1L1', 'ZBTB11', sample(colnames(xdata.raw), 100)) %>%
unique
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.
set.seed(params$seed)
fitted <- cv.glmHub(xdata, Surv(ydata$time, ydata$status),
family = 'cox',
lambda = buildLambda(1),
network = 'correlation',
network.options = networkOptions(cutoff = .6,
min.degree = .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
coefs.v <- coef(fitted, s = 'lambda.min')[,1] %>% { .[. != 0]}
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
coefs.v %>% {
data.frame(gene.name = names(.),
coefficient = .,
stringsAsFactors = FALSE)
} %>%
arrange(gene.name) %>%
knitr::kable()
| gene.name | coefficient | |
|---|---|---|
| CD5 | CD5 | -0.16632 |
names(coefs.v) %>% { hallmarks(.)$heatmap }
## Error in curl::curl_fetch_memory(url, handle = handle): Timeout was reached: [chat.lionproject.net] Connection timeout after 10005 ms
## Request failed [ERROR]. Retrying in 1.6 seconds...
## Error in curl::curl_fetch_memory(url, handle = handle): Failed to connect to chat.lionproject.net port 443 after 2544 ms: Connection refused
## Request failed [ERROR]. Retrying in 3 seconds...
## Cannot call Hallmark API, please try again later.
## NULL
separate2GroupsCox(as.vector(coefs.v),
xdata[, names(coefs.v)],
ydata,
plot.title = 'Full dataset', legend.outside = FALSE)
## $pvalue
## [1] 0.001237802
##
## $plot
##
## $km
## Call: survfit(formula = survival::Surv(time, status) ~ group, data = prognostic.index.df)
##
## n events median 0.95LCL 0.95UCL
## Low risk 540 58 3959 3492 NA
## High risk 540 94 3738 3262 4456
sessionInfo()
## R version 4.3.0 RC (2023-04-18 r84287 ucrt)
## Platform: x86_64-w64-mingw32/x64 (64-bit)
## Running under: Windows Server 2022 x64 (build 20348)
##
## Matrix products: default
##
##
## locale:
## [1] LC_COLLATE=C
## [2] LC_CTYPE=English_United States.utf8
## [3] LC_MONETARY=English_United States.utf8
## [4] LC_NUMERIC=C
## [5] LC_TIME=English_United States.utf8
##
## 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.19.0
## [5] glmnet_4.1-7 Matrix_1.5-4
## [7] TCGAutils_1.21.1 curatedTCGAData_1.23.2
## [9] MultiAssayExperiment_1.27.0 SummarizedExperiment_1.31.1
## [11] Biobase_2.61.0 GenomicRanges_1.53.1
## [13] GenomeInfoDb_1.37.1 IRanges_2.35.1
## [15] S4Vectors_0.39.1 BiocGenerics_0.47.0
## [17] MatrixGenerics_1.13.0 matrixStats_0.63.0
## [19] futile.logger_1.4.3 survival_3.5-5
## [21] ggplot2_3.4.2 dplyr_1.1.2
## [23] BiocStyle_2.29.0
##
## loaded via a namespace (and not attached):
## [1] jsonlite_1.8.4 shape_1.4.6
## [3] magrittr_2.0.3 magick_2.7.4
## [5] GenomicFeatures_1.53.0 farver_2.1.1
## [7] rmarkdown_2.21 BiocIO_1.11.0
## [9] zlibbioc_1.47.0 vctrs_0.6.2
## [11] memoise_2.0.1 Rsamtools_2.17.0
## [13] RCurl_1.98-1.12 rstatix_0.7.2
## [15] htmltools_0.5.5 S4Arrays_1.1.2
## [17] BiocBaseUtils_1.3.0 progress_1.2.2
## [19] AnnotationHub_3.9.1 lambda.r_1.2.4
## [21] curl_5.0.0 broom_1.0.4
## [23] pROC_1.18.0 SparseArray_1.1.2
## [25] sass_0.4.6 bslib_0.4.2
## [27] plyr_1.8.8 zoo_1.8-12
## [29] futile.options_1.0.1 cachem_1.0.8
## [31] GenomicAlignments_1.37.0 mime_0.12
## [33] lifecycle_1.0.3 iterators_1.0.14
## [35] pkgconfig_2.0.3 R6_2.5.1
## [37] fastmap_1.1.1 GenomeInfoDbData_1.2.10
## [39] shiny_1.7.4 digest_0.6.31
## [41] colorspace_2.1-0 AnnotationDbi_1.63.1
## [43] ExperimentHub_2.9.0 RSQLite_2.3.1
## [45] ggpubr_0.6.0 filelock_1.0.2
## [47] labeling_0.4.2 km.ci_0.5-6
## [49] fansi_1.0.4 abind_1.4-5
## [51] httr_1.4.5 compiler_4.3.0
## [53] bit64_4.0.5 withr_2.5.0
## [55] backports_1.4.1 BiocParallel_1.35.0
## [57] carData_3.0-5 DBI_1.1.3
## [59] highr_0.10 ggsignif_0.6.4
## [61] biomaRt_2.57.0 rappdirs_0.3.3
## [63] DelayedArray_0.27.2 rjson_0.2.21
## [65] tools_4.3.0 interactiveDisplayBase_1.39.0
## [67] httpuv_1.6.9 glue_1.6.2
## [69] restfulr_0.0.15 promises_1.2.0.1
## [71] generics_0.1.3 gtable_0.3.3
## [73] KMsurv_0.1-5 tzdb_0.3.0
## [75] tidyr_1.3.0 survminer_0.4.9
## [77] data.table_1.14.8 hms_1.1.3
## [79] car_3.1-2 xml2_1.3.4
## [81] utf8_1.2.3 XVector_0.41.1
## [83] BiocVersion_3.18.0 foreach_1.5.2
## [85] pillar_1.9.0 stringr_1.5.0
## [87] later_1.3.1 splines_4.3.0
## [89] BiocFileCache_2.9.0 lattice_0.21-8
## [91] rtracklayer_1.61.0 bit_4.0.5
## [93] tidyselect_1.2.0 Biostrings_2.69.0
## [95] knitr_1.42 gridExtra_2.3
## [97] bookdown_0.33 xfun_0.39
## [99] stringi_1.7.12 yaml_2.3.7
## [101] evaluate_0.20 codetools_0.2-19
## [103] tibble_3.2.1 BiocManager_1.30.20
## [105] cli_3.6.1 xtable_1.8-4
## [107] munsell_0.5.0 jquerylib_0.1.4
## [109] survMisc_0.5.6 Rcpp_1.0.10
## [111] GenomicDataCommons_1.25.0 dbplyr_2.3.2
## [113] png_0.1-8 XML_3.99-0.14
## [115] ellipsis_0.3.2 readr_2.1.4
## [117] blob_1.2.4 prettyunits_1.1.1
## [119] bitops_1.0-7 scales_1.2.1
## [121] purrr_1.0.1 crayon_1.5.2
## [123] rlang_1.1.1 KEGGREST_1.41.0
## [125] rvest_1.0.3 formatR_1.14