Contents

0.1 Instalation

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

1 Required Packages

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

2 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 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.

skcm <- curatedTCGAData(
    diseaseCode = "SKCM", assays = "RNASeq2GeneNorm",
    version = "1.1.38", dry.run = FALSE
)

Build the survival data from the clinical columns.

skcmMetastatic <- TCGAutils::TCGAsplitAssays(skcm, "06")
xdataRaw <- t(assay(skcmMetastatic[[1]]))

# Get survival information
ydataRaw <- colData(skcmMetastatic) |>
    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()

# Get survival information
ydataRaw <- colData(skcm) |>
    as.data.frame() |>
    # Find max time between all days (ignoring missings)
    dplyr::filter(
        !is.na(days_to_last_followup) | !is.na(days_to_death)
    ) |>
    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(ydataRaw) <- ydataRaw$patientID

# keep only features that have standard deviation > 0
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)), ]

set.seed(params$seed)
smallSubset <- c(
    "FOXL2", "KLHL5", "PCYT2", "SLC6A10P", "STRAP", "TMEM33",
    "WT1-AS", sample(colnames(xdataRaw), 100)
)

xdata <- xdataRaw[, smallSubset[smallSubset %in% colnames(xdataRaw)]]
ydata <- ydataRaw |> dplyr::select(time, status)

3 Fit models

Fit model model penalizing by the hubs using the cross-validation function by cv.glmHub.

fitted <- cv.glmHub(
    xdata,
    Surv(ydata$time, ydata$status),
    family = "cox",
    foldid = glmSparseNet:::balancedCvFolds(ydata$status)$output,
    network = "correlation",
    options = networkOptions(
        cutoff = .6,
        minDegree = .2
    )
)

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

4.1 Coefficients of selected model from Cross-Validation

Taking the best model described by lambda.min

coefsCV <- Filter(function(.x) .x != 0, coef(fitted, s = "lambda.min")[, 1])
data.frame(
    ensembl.id = names(coefsCV),
    gene.name = geneNames(names(coefsCV))$external_gene_name,
    coefficient = coefsCV,
    stringsAsFactors = FALSE
) |>
    arrange(gene.name) |>
    knitr::kable()
ensembl.id gene.name coefficient
AMICA1 AMICA1 AMICA1 -0.2758400
C4orf49 C4orf49 C4orf49 -0.0059089
PCYT2 PCYT2 PCYT2 0.0646641

4.2 Survival curves and Log rank test

separate2GroupsCox(as.vector(coefsCV),
    xdata[, names(coefsCV)],
    ydata,
    plotTitle = "Full dataset", legendOutside = FALSE
)
## $pvalue
## [1] 0.0001269853
## 
## $plot

## 
## $km
## Call: survfit(formula = survival::Surv(time, status) ~ group, data = prognosticIndexDf)
## 
##                 n events median 0.95LCL 0.95UCL
## Low risk - 1  180     79   4000    2927    6164
## High risk - 1 179    114   2005    1524    2829

5 Session Info

sessionInfo()
## R Under development (unstable) (2024-10-26 r87273 ucrt)
## Platform: x86_64-w64-mingw32/x64
## 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] glmnet_4.1-8                VennDiagram_1.7.3          
##  [3] reshape2_1.4.4              forcats_1.0.0              
##  [5] Matrix_1.7-1                glmSparseNet_1.25.0        
##  [7] TCGAutils_1.27.0            curatedTCGAData_1.27.1     
##  [9] MultiAssayExperiment_1.33.0 SummarizedExperiment_1.37.0
## [11] Biobase_2.67.0              GenomicRanges_1.59.0       
## [13] GenomeInfoDb_1.43.0         IRanges_2.41.0             
## [15] S4Vectors_0.45.0            BiocGenerics_0.53.1        
## [17] generics_0.1.3              MatrixGenerics_1.19.0      
## [19] matrixStats_1.4.1           futile.logger_1.4.3        
## [21] survival_3.7-0              ggplot2_3.5.1              
## [23] dplyr_1.1.4                 BiocStyle_2.35.0           
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.9            shape_1.4.6.1            
##   [3] magrittr_2.0.3            magick_2.8.5             
##   [5] GenomicFeatures_1.59.0    farver_2.1.2             
##   [7] rmarkdown_2.28            BiocIO_1.17.0            
##   [9] zlibbioc_1.53.0           vctrs_0.6.5              
##  [11] memoise_2.0.1             Rsamtools_2.23.0         
##  [13] RCurl_1.98-1.16           rstatix_0.7.2            
##  [15] tinytex_0.54              progress_1.2.3           
##  [17] htmltools_0.5.8.1         S4Arrays_1.7.1           
##  [19] BiocBaseUtils_1.9.0       AnnotationHub_3.15.0     
##  [21] lambda.r_1.2.4            curl_5.2.3               
##  [23] broom_1.0.7               Formula_1.2-5            
##  [25] pROC_1.18.5               SparseArray_1.7.0        
##  [27] sass_0.4.9                bslib_0.8.0              
##  [29] plyr_1.8.9                httr2_1.0.5              
##  [31] zoo_1.8-12                futile.options_1.0.1     
##  [33] cachem_1.1.0              GenomicAlignments_1.43.0 
##  [35] mime_0.12                 lifecycle_1.0.4          
##  [37] iterators_1.0.14          pkgconfig_2.0.3          
##  [39] R6_2.5.1                  fastmap_1.2.0            
##  [41] GenomeInfoDbData_1.2.13   digest_0.6.37            
##  [43] colorspace_2.1-1          AnnotationDbi_1.69.0     
##  [45] ps_1.8.1                  ExperimentHub_2.15.0     
##  [47] RSQLite_2.3.7             ggpubr_0.6.0             
##  [49] labeling_0.4.3            filelock_1.0.3           
##  [51] km.ci_0.5-6               fansi_1.0.6              
##  [53] httr_1.4.7                abind_1.4-8              
##  [55] compiler_4.5.0            bit64_4.5.2              
##  [57] withr_3.0.2               backports_1.5.0          
##  [59] BiocParallel_1.41.0       carData_3.0-5            
##  [61] DBI_1.2.3                 highr_0.11               
##  [63] ggsignif_0.6.4            biomaRt_2.63.0           
##  [65] rappdirs_0.3.3            DelayedArray_0.33.1      
##  [67] rjson_0.2.23              tools_4.5.0              
##  [69] chromote_0.3.1            glue_1.8.0               
##  [71] restfulr_0.0.15           promises_1.3.0           
##  [73] checkmate_2.3.2           gtable_0.3.6             
##  [75] KMsurv_0.1-5              tzdb_0.4.0               
##  [77] tidyr_1.3.1               survminer_0.5.0          
##  [79] websocket_1.4.2           data.table_1.16.2        
##  [81] hms_1.1.3                 car_3.1-3                
##  [83] xml2_1.3.6                utf8_1.2.4               
##  [85] XVector_0.47.0            BiocVersion_3.21.1       
##  [87] foreach_1.5.2             pillar_1.9.0             
##  [89] stringr_1.5.1             later_1.3.2              
##  [91] splines_4.5.0             BiocFileCache_2.15.0     
##  [93] lattice_0.22-6            rtracklayer_1.67.0       
##  [95] bit_4.5.0                 tidyselect_1.2.1         
##  [97] Biostrings_2.75.0         knitr_1.48               
##  [99] gridExtra_2.3             bookdown_0.41            
## [101] xfun_0.49                 stringi_1.8.4            
## [103] UCSC.utils_1.3.0          yaml_2.3.10              
## [105] evaluate_1.0.1            codetools_0.2-20         
## [107] tibble_3.2.1              BiocManager_1.30.25      
## [109] cli_3.6.3                 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.13-1             GenomicDataCommons_1.31.0
## [117] dbplyr_2.5.0              png_0.1-8                
## [119] XML_3.99-0.17             readr_2.1.5              
## [121] blob_1.2.4                prettyunits_1.2.0        
## [123] bitops_1.0-9              scales_1.3.0             
## [125] purrr_1.0.2               crayon_1.5.3             
## [127] rlang_1.1.4               KEGGREST_1.47.0          
## [129] rvest_1.0.4               formatR_1.14