Contents

1 Instalation

if (!require("BiocManager"))
  install.packages("BiocManager")
BiocManager::install("glmSparseNet")

2 Required Packages

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())

3 Load data

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 around 100 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.

prad <- curatedTCGAData(diseaseCode = "PRAD", assays = "RNASeq2GeneNorm",
                        version = '1.1.38', dry.run = FALSE))

Build the survival data from the clinical columns.

# keep only solid tumour (code: 01)
prad.primary.solid.tumor <- TCGAutils::TCGAsplitAssays(prad, '01')
xdata.raw <- t(assay(prad.primary.solid.tumor[[1]]))

# Get survival information
ydata.raw <- colData(prad.primary.solid.tumor) %>% 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(geneNames(c('ENSG00000103091', 'ENSG00000064787', 
                              'ENSG00000119915', 'ENSG00000120158', 
                              'ENSG00000114491', 'ENSG00000204176', 
                              'ENSG00000138399'))$external_gene_name, 
                  sample(colnames(xdata.raw), 100)) %>% unique %>% sort

xdata <- xdata.raw[, small.subset[small.subset %in% colnames(xdata.raw)]]
ydata <- ydata.raw %>% dplyr::select(time, status)

4 Fit models

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',
                    nlambda = 1000,
                    network = 'correlation', 
                    network.options = networkOptions(cutoff = .6, 
                                                     min.degree = .2))

5 Results of Cross Validation

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)

5.1 Coefficients of selected model from Cross-Validation

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
AKAP9 AKAP9 AKAP9 0.2616307
ALPK2 ALPK2 ALPK2 -0.0714527
ATP5G2 ATP5G2 ATP5G2 -0.2575987
C22orf32 C22orf32 C22orf32 -0.2119992
CSNK2A1P CSNK2A1P CSNK2A1P -1.4875518
MYST3 MYST3 MYST3 -1.6177076
NBPF10 NBPF10 NBPF10 0.4507147
PFN1 PFN1 PFN1 0.4161846
SCGB2A2 SCGB2A2 SCGB2A2 0.0749064
SLC25A1 SLC25A1 SLC25A1 -0.8484827
STX4 STX4 STX4 -0.1690185
SYP SYP SYP 0.2425939
TMEM141 TMEM141 TMEM141 -0.8273147
UMPS UMPS UMPS 0.2214068
ZBTB26 ZBTB26 ZBTB26 0.3696515

5.2 Hallmarks of Cancer

geneNames(names(coefs.v)) %>% { hallmarks(.$external_gene_name)$heatmap }
## Error in curl::curl_fetch_memory(url, handle = handle): error:0A000126:SSL routines::unexpected eof while reading
## Request failed [ERROR]. Retrying in 2 seconds...
## Error in curl::curl_fetch_memory(url, handle = handle): error:0A000126:SSL routines::unexpected eof while reading
## Request failed [ERROR]. Retrying in 2 seconds...
## Cannot call Hallmark API, please try again later.
## NULL

5.3 Survival curves and Log rank test

separate2GroupsCox(as.vector(coefs.v), 
                   xdata[, names(coefs.v)], 
                   ydata, 
                   plot.title = 'Full dataset', legend.outside = FALSE)
## $pvalue
## [1] 0.001155155
## 
## $plot

## 
## $km
## Call: survfit(formula = survival::Surv(time, status) ~ group, data = prognostic.index.df)
## 
##             n events median 0.95LCL 0.95UCL
## Low risk  249      0     NA      NA      NA
## High risk 248     10   3502    3467      NA

6 Session Info

sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.17-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] VennDiagram_1.7.3           reshape2_1.4.4             
##  [3] forcats_1.0.0               glmSparseNet_1.18.0        
##  [5] glmnet_4.1-7                Matrix_1.5-4               
##  [7] TCGAutils_1.20.0            curatedTCGAData_1.21.3     
##  [9] MultiAssayExperiment_1.26.0 SummarizedExperiment_1.30.0
## [11] Biobase_2.60.0              GenomicRanges_1.52.0       
## [13] GenomeInfoDb_1.36.0         IRanges_2.34.0             
## [15] S4Vectors_0.38.0            BiocGenerics_0.46.0        
## [17] MatrixGenerics_1.12.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.28.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.52.0        farver_2.1.1                 
##   [7] rmarkdown_2.21                BiocIO_1.10.0                
##   [9] zlibbioc_1.46.0               vctrs_0.6.2                  
##  [11] memoise_2.0.1                 Rsamtools_2.16.0             
##  [13] RCurl_1.98-1.12               rstatix_0.7.2                
##  [15] htmltools_0.5.5               progress_1.2.2               
##  [17] AnnotationHub_3.8.0           lambda.r_1.2.4               
##  [19] curl_5.0.0                    broom_1.0.4                  
##  [21] pROC_1.18.0                   sass_0.4.5                   
##  [23] bslib_0.4.2                   plyr_1.8.8                   
##  [25] zoo_1.8-12                    futile.options_1.0.1         
##  [27] cachem_1.0.7                  GenomicAlignments_1.36.0     
##  [29] mime_0.12                     lifecycle_1.0.3              
##  [31] iterators_1.0.14              pkgconfig_2.0.3              
##  [33] R6_2.5.1                      fastmap_1.1.1                
##  [35] GenomeInfoDbData_1.2.10       shiny_1.7.4                  
##  [37] digest_0.6.31                 colorspace_2.1-0             
##  [39] AnnotationDbi_1.62.0          ExperimentHub_2.8.0          
##  [41] RSQLite_2.3.1                 ggpubr_0.6.0                 
##  [43] filelock_1.0.2                labeling_0.4.2               
##  [45] km.ci_0.5-6                   fansi_1.0.4                  
##  [47] abind_1.4-5                   httr_1.4.5                   
##  [49] compiler_4.3.0                bit64_4.0.5                  
##  [51] withr_2.5.0                   backports_1.4.1              
##  [53] BiocParallel_1.34.0           carData_3.0-5                
##  [55] DBI_1.1.3                     highr_0.10                   
##  [57] ggsignif_0.6.4                biomaRt_2.56.0               
##  [59] rappdirs_0.3.3                DelayedArray_0.26.0          
##  [61] rjson_0.2.21                  tools_4.3.0                  
##  [63] interactiveDisplayBase_1.38.0 httpuv_1.6.9                 
##  [65] glue_1.6.2                    restfulr_0.0.15              
##  [67] promises_1.2.0.1              generics_0.1.3               
##  [69] gtable_0.3.3                  KMsurv_0.1-5                 
##  [71] tzdb_0.3.0                    tidyr_1.3.0                  
##  [73] survminer_0.4.9               data.table_1.14.8            
##  [75] hms_1.1.3                     car_3.1-2                    
##  [77] xml2_1.3.3                    utf8_1.2.3                   
##  [79] XVector_0.40.0                BiocVersion_3.17.1           
##  [81] foreach_1.5.2                 pillar_1.9.0                 
##  [83] stringr_1.5.0                 later_1.3.0                  
##  [85] splines_4.3.0                 BiocFileCache_2.8.0          
##  [87] lattice_0.21-8                rtracklayer_1.60.0           
##  [89] bit_4.0.5                     tidyselect_1.2.0             
##  [91] Biostrings_2.68.0             knitr_1.42                   
##  [93] gridExtra_2.3                 bookdown_0.33                
##  [95] xfun_0.39                     stringi_1.7.12               
##  [97] yaml_2.3.7                    evaluate_0.20                
##  [99] codetools_0.2-19              tibble_3.2.1                 
## [101] BiocManager_1.30.20           cli_3.6.1                    
## [103] xtable_1.8-4                  munsell_0.5.0                
## [105] jquerylib_0.1.4               survMisc_0.5.6               
## [107] Rcpp_1.0.10                   GenomicDataCommons_1.24.0    
## [109] dbplyr_2.3.2                  png_0.1-8                    
## [111] XML_3.99-0.14                 ellipsis_0.3.2               
## [113] readr_2.1.4                   blob_1.2.4                   
## [115] prettyunits_1.1.1             bitops_1.0-7                 
## [117] scales_1.2.1                  purrr_1.0.1                  
## [119] crayon_1.5.2                  rlang_1.1.0                  
## [121] KEGGREST_1.40.0               rvest_1.0.3                  
## [123] formatR_1.14