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.

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)

3 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",
    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

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])
## 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

4.2 Survival curves and Log rank test

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

5 Session Info

sessionInfo()
## R Under development (unstable) (2024-01-16 r85808)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-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.6-5                glmSparseNet_1.21.2        
##  [7] TCGAutils_1.23.3            curatedTCGAData_1.25.4     
##  [9] MultiAssayExperiment_1.29.1 SummarizedExperiment_1.33.3
## [11] Biobase_2.63.0              GenomicRanges_1.55.3       
## [13] GenomeInfoDb_1.39.6         IRanges_2.37.1             
## [15] S4Vectors_0.41.3            BiocGenerics_0.49.1        
## [17] MatrixGenerics_1.15.0       matrixStats_1.2.0          
## [19] futile.logger_1.4.3         survival_3.5-8             
## [21] ggplot2_3.4.4               dplyr_1.1.4                
## [23] BiocStyle_2.31.0           
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.8            shape_1.4.6              
##   [3] magrittr_2.0.3            magick_2.8.2             
##   [5] GenomicFeatures_1.55.3    farver_2.1.1             
##   [7] rmarkdown_2.25            BiocIO_1.13.0            
##   [9] zlibbioc_1.49.0           vctrs_0.6.5              
##  [11] memoise_2.0.1             Rsamtools_2.19.3         
##  [13] RCurl_1.98-1.14           rstatix_0.7.2            
##  [15] htmltools_0.5.7           S4Arrays_1.3.3           
##  [17] BiocBaseUtils_1.5.0       progress_1.2.3           
##  [19] AnnotationHub_3.11.1      lambda.r_1.2.4           
##  [21] curl_5.2.0                broom_1.0.5              
##  [23] pROC_1.18.5               SparseArray_1.3.4        
##  [25] sass_0.4.8                bslib_0.6.1              
##  [27] plyr_1.8.9                httr2_1.0.0              
##  [29] zoo_1.8-12                futile.options_1.0.1     
##  [31] cachem_1.0.8              GenomicAlignments_1.39.4 
##  [33] mime_0.12                 lifecycle_1.0.4          
##  [35] iterators_1.0.14          pkgconfig_2.0.3          
##  [37] R6_2.5.1                  fastmap_1.1.1            
##  [39] GenomeInfoDbData_1.2.11   digest_0.6.34            
##  [41] colorspace_2.1-0          AnnotationDbi_1.65.2     
##  [43] ps_1.7.6                  ExperimentHub_2.11.1     
##  [45] RSQLite_2.3.5             ggpubr_0.6.0             
##  [47] labeling_0.4.3            filelock_1.0.3           
##  [49] km.ci_0.5-6               fansi_1.0.6              
##  [51] httr_1.4.7                abind_1.4-5              
##  [53] compiler_4.4.0            bit64_4.0.5              
##  [55] withr_3.0.0               backports_1.4.1          
##  [57] BiocParallel_1.37.0       carData_3.0-5            
##  [59] DBI_1.2.1                 highr_0.10               
##  [61] ggsignif_0.6.4            biomaRt_2.59.1           
##  [63] rappdirs_0.3.3            DelayedArray_0.29.2      
##  [65] rjson_0.2.21              tools_4.4.0              
##  [67] chromote_0.2.0            glue_1.7.0               
##  [69] restfulr_0.0.15           promises_1.2.1           
##  [71] checkmate_2.3.1           generics_0.1.3           
##  [73] gtable_0.3.4              KMsurv_0.1-5             
##  [75] tzdb_0.4.0                tidyr_1.3.1              
##  [77] survminer_0.4.9           websocket_1.4.1          
##  [79] data.table_1.15.0         hms_1.1.3                
##  [81] car_3.1-2                 xml2_1.3.6               
##  [83] utf8_1.2.4                XVector_0.43.1           
##  [85] BiocVersion_3.19.1        foreach_1.5.2            
##  [87] pillar_1.9.0              stringr_1.5.1            
##  [89] later_1.3.2               splines_4.4.0            
##  [91] BiocFileCache_2.11.1      lattice_0.22-5           
##  [93] rtracklayer_1.63.0        bit_4.0.5                
##  [95] tidyselect_1.2.0          Biostrings_2.71.2        
##  [97] knitr_1.45                gridExtra_2.3            
##  [99] bookdown_0.37             xfun_0.42                
## [101] stringi_1.8.3             yaml_2.3.8               
## [103] evaluate_0.23             codetools_0.2-19         
## [105] tibble_3.2.1              BiocManager_1.30.22      
## [107] cli_3.6.2                 xtable_1.8-4             
## [109] munsell_0.5.0             processx_3.8.3           
## [111] jquerylib_0.1.4           survMisc_0.5.6           
## [113] Rcpp_1.0.12               GenomicDataCommons_1.27.1
## [115] dbplyr_2.4.0              png_0.1-8                
## [117] XML_3.99-0.16.1           readr_2.1.5              
## [119] blob_1.2.4                prettyunits_1.2.0        
## [121] bitops_1.0-7              scales_1.3.0             
## [123] purrr_1.0.2               crayon_1.5.2             
## [125] rlang_1.1.3               KEGGREST_1.43.0          
## [127] rvest_1.0.4               formatR_1.14