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): OpenSSL SSL_connect: SSL_ERROR_SYSCALL in connection to chat.lionproject.net:443
## Request failed [ERROR]. Retrying in 1.6 seconds...
## Error in curl::curl_fetch_memory(url, handle = handle): Timeout was reached: [chat.lionproject.net] Operation timed out after 10001 milliseconds with 0 out of 0 bytes received
## 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.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
##
## 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
##
## 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_0.5.2 glmSparseNet_1.16.0
## [5] glmnet_4.1-4 Matrix_1.5-1
## [7] TCGAutils_1.18.0 curatedTCGAData_1.19.2
## [9] MultiAssayExperiment_1.24.0 SummarizedExperiment_1.28.0
## [11] Biobase_2.58.0 GenomicRanges_1.50.0
## [13] GenomeInfoDb_1.34.0 IRanges_2.32.0
## [15] S4Vectors_0.36.0 BiocGenerics_0.44.0
## [17] MatrixGenerics_1.10.0 matrixStats_0.62.0
## [19] futile.logger_1.4.3 survival_3.4-0
## [21] ggplot2_3.3.6 dplyr_1.0.10
## [23] BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] backports_1.4.1 AnnotationHub_3.6.0
## [3] BiocFileCache_2.6.0 plyr_1.8.7
## [5] splines_4.2.1 BiocParallel_1.32.0
## [7] digest_0.6.30 foreach_1.5.2
## [9] htmltools_0.5.3 magick_2.7.3
## [11] fansi_1.0.3 magrittr_2.0.3
## [13] memoise_2.0.1 tzdb_0.3.0
## [15] Biostrings_2.66.0 readr_2.1.3
## [17] prettyunits_1.1.1 colorspace_2.0-3
## [19] blob_1.2.3 rvest_1.0.3
## [21] rappdirs_0.3.3 xfun_0.34
## [23] crayon_1.5.2 RCurl_1.98-1.9
## [25] jsonlite_1.8.3 zoo_1.8-11
## [27] iterators_1.0.14 glue_1.6.2
## [29] survminer_0.4.9 GenomicDataCommons_1.22.0
## [31] gtable_0.3.1 zlibbioc_1.44.0
## [33] XVector_0.38.0 DelayedArray_0.24.0
## [35] car_3.1-1 shape_1.4.6
## [37] abind_1.4-5 scales_1.2.1
## [39] futile.options_1.0.1 DBI_1.1.3
## [41] rstatix_0.7.0 Rcpp_1.0.9
## [43] xtable_1.8-4 progress_1.2.2
## [45] bit_4.0.4 km.ci_0.5-6
## [47] httr_1.4.4 ellipsis_0.3.2
## [49] pkgconfig_2.0.3 XML_3.99-0.12
## [51] farver_2.1.1 sass_0.4.2
## [53] dbplyr_2.2.1 utf8_1.2.2
## [55] tidyselect_1.2.0 labeling_0.4.2
## [57] rlang_1.0.6 later_1.3.0
## [59] AnnotationDbi_1.60.0 munsell_0.5.0
## [61] BiocVersion_3.16.0 tools_4.2.1
## [63] cachem_1.0.6 cli_3.4.1
## [65] generics_0.1.3 RSQLite_2.2.18
## [67] ExperimentHub_2.6.0 broom_1.0.1
## [69] evaluate_0.17 stringr_1.4.1
## [71] fastmap_1.1.0 yaml_2.3.6
## [73] knitr_1.40 bit64_4.0.5
## [75] survMisc_0.5.6 purrr_0.3.5
## [77] KEGGREST_1.38.0 mime_0.12
## [79] formatR_1.12 xml2_1.3.3
## [81] biomaRt_2.54.0 compiler_4.2.1
## [83] filelock_1.0.2 curl_4.3.3
## [85] png_0.1-7 interactiveDisplayBase_1.36.0
## [87] ggsignif_0.6.4 tibble_3.1.8
## [89] bslib_0.4.0 stringi_1.7.8
## [91] highr_0.9 GenomicFeatures_1.50.0
## [93] lattice_0.20-45 KMsurv_0.1-5
## [95] vctrs_0.5.0 pillar_1.8.1
## [97] lifecycle_1.0.3 BiocManager_1.30.19
## [99] jquerylib_0.1.4 data.table_1.14.4
## [101] bitops_1.0-7 httpuv_1.6.6
## [103] rtracklayer_1.58.0 R6_2.5.1
## [105] BiocIO_1.8.0 bookdown_0.29
## [107] promises_1.2.0.1 gridExtra_2.3
## [109] codetools_0.2-18 lambda.r_1.2.4
## [111] assertthat_0.2.1 rjson_0.2.21
## [113] withr_2.5.0 GenomicAlignments_1.34.0
## [115] Rsamtools_2.14.0 GenomeInfoDbData_1.2.9
## [117] hms_1.1.2 tidyr_1.2.1
## [119] rmarkdown_2.17 carData_3.0-5
## [121] ggpubr_0.4.0 pROC_1.18.0
## [123] shiny_1.7.3 restfulr_0.0.15