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
)
brca <- curatedTCGAData(diseaseCode = "BRCA", assays = "RNASeq2GeneNorm",
version = "1.1.38", dry.run = FALSE)
brca <- TCGAutils::TCGAsplitAssays(brca, c('01','11'))
xdata.raw <- t(cbind(assay(brca[[1]]), assay(brca[[2]])))
# Get matches between survival and assay data
class.v <- TCGAbiospec(rownames(xdata.raw))$sample_definition %>% factor
names(class.v) <- rownames(xdata.raw)
# keep features with standard deviation > 0
xdata.raw <- xdata.raw %>%
{ (apply(., 2, sd) != 0) } %>%
{ xdata.raw[, .] } %>%
scale
set.seed(params$seed)
small.subset <- c('CD5', 'CSF2RB', 'HSF1', 'IRGC', 'LRRC37A6P', 'NEUROG2',
'NLRC4', 'PDE11A', 'PIK3CB', 'QARS', 'RPGRIP1L', 'SDC1',
'TMEM31', 'YME1L1', 'ZBTB11',
sample(colnames(xdata.raw), 100))
xdata <- xdata.raw[, small.subset[small.subset %in% colnames(xdata.raw)]]
ydata <- class.v
Fit model model penalizing by the hubs using the cross-validation function by
cv.glmHub.
fitted <- cv.glmHub(xdata, ydata,
family = 'binomial',
network = 'correlation',
nlambda = 1000,
network.options = networkOptions(cutoff = .6,
min.degree = .2))
Shows the results of 1000 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]}
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 | |
|---|---|---|---|
| (Intercept) | (Intercept) | (Intercept) | -6.8189813 |
| CD5 | CD5 | AMOTL1 | -1.1200445 |
| NLRC4 | NLRC4 | ATR | -1.4434578 |
| PIK3CB | PIK3CB | B3GALT2 | -0.3880002 |
| ZBTB11 | ZBTB11 | BAG2 | -0.3325729 |
| ATR | ATR | C16orf82 | 1.2498304 |
| IL2 | IL2 | CD5 | 0.6327083 |
| GDF11 | GDF11 | CIITA | -0.2676642 |
| DCP1A | DCP1A | DCP1A | 0.2994599 |
| AMOTL1 | AMOTL1 | FAM86B1 | 0.4430643 |
| BAG2 | BAG2 | FNIP2 | -0.1841676 |
| C16orf82 | C16orf82 | GDF11 | 0.0396368 |
| FAM86B1 | FAM86B1 | GNG11 | 0.2025463 |
| FNIP2 | FNIP2 | GREM2 | 0.6101759 |
| MS4A4A | MS4A4A | GZMB | 1.1614779 |
| B3GALT2 | B3GALT2 | HAX1 | -0.0867011 |
| GNG11 | GNG11 | IL2 | 3.0659066 |
| NDRG2 | NDRG2 | MMP28 | 1.1142519 |
| HAX1 | HAX1 | MS4A4A | -0.1516837 |
| GREM2 | GREM2 | NDRG2 | -0.2014884 |
| CIITA | CIITA | NLRC4 | 0.4256103 |
| GZMB | GZMB | PIK3CB | -2.7663574 |
| MMP28 | MMP28 | ZBTB11 | -0.8438024 |
geneNames(names(coefs.v)) %>% { hallmarks(.$external_gene_name)$heatmap }
## Error in curl::curl_fetch_memory(url, handle = handle): Failed to connect to chat.lionproject.net port 443 after 2680 ms: Connection refused
## Request failed [ERROR]. Retrying in 1.1 seconds...
## Error in curl::curl_fetch_memory(url, handle = handle): Failed to connect to chat.lionproject.net port 443 after 7340 ms: Connection refused
## Request failed [ERROR]. Retrying in 1.1 seconds...
## Cannot call Hallmark API, please try again later.
## NULL
## [INFO] Misclassified (11)
## [INFO] * False primary solid tumour: 7
## [INFO] * False normal : 4
Histogram of predicted response
ROC curve
## Setting levels: control = Primary Solid Tumor, case = Solid Tissue Normal
## Setting direction: controls < cases
## Warning: Using `size` aesthetic for lines was deprecated in ggplot2 3.4.0.
## ℹ Please use `linewidth` instead.
## This warning is displayed once every 8 hours.
## Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
## generated.
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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] glmSparseNet_1.19.0 glmnet_4.1-7
## [3] Matrix_1.5-4 TCGAutils_1.21.1
## [5] curatedTCGAData_1.23.2 MultiAssayExperiment_1.27.0
## [7] SummarizedExperiment_1.31.1 Biobase_2.61.0
## [9] GenomicRanges_1.53.1 GenomeInfoDb_1.37.1
## [11] IRanges_2.35.1 S4Vectors_0.39.1
## [13] BiocGenerics_0.47.0 MatrixGenerics_1.13.0
## [15] matrixStats_0.63.0 futile.logger_1.4.3
## [17] survival_3.5-5 ggplot2_3.4.2
## [19] dplyr_1.1.2 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 htmltools_0.5.5
## [15] S4Arrays_1.1.2 forcats_1.0.0
## [17] BiocBaseUtils_1.3.0 progress_1.2.2
## [19] AnnotationHub_3.9.1 lambda.r_1.2.4
## [21] curl_5.0.0 pROC_1.18.0
## [23] SparseArray_1.1.2 sass_0.4.6
## [25] bslib_0.4.2 plyr_1.8.8
## [27] futile.options_1.0.1 cachem_1.0.8
## [29] GenomicAlignments_1.37.0 mime_0.12
## [31] lifecycle_1.0.3 iterators_1.0.14
## [33] pkgconfig_2.0.3 R6_2.5.1
## [35] fastmap_1.1.1 GenomeInfoDbData_1.2.10
## [37] shiny_1.7.4 digest_0.6.31
## [39] colorspace_2.1-0 AnnotationDbi_1.63.1
## [41] ExperimentHub_2.9.0 RSQLite_2.3.1
## [43] filelock_1.0.2 labeling_0.4.2
## [45] fansi_1.0.4 httr_1.4.5
## [47] compiler_4.3.0 bit64_4.0.5
## [49] withr_2.5.0 BiocParallel_1.35.0
## [51] DBI_1.1.3 highr_0.10
## [53] biomaRt_2.57.0 rappdirs_0.3.3
## [55] DelayedArray_0.27.2 rjson_0.2.21
## [57] tools_4.3.0 interactiveDisplayBase_1.39.0
## [59] httpuv_1.6.9 glue_1.6.2
## [61] restfulr_0.0.15 promises_1.2.0.1
## [63] grid_4.3.0 generics_0.1.3
## [65] gtable_0.3.3 tzdb_0.3.0
## [67] hms_1.1.3 xml2_1.3.4
## [69] utf8_1.2.3 XVector_0.41.1
## [71] BiocVersion_3.18.0 foreach_1.5.2
## [73] pillar_1.9.0 stringr_1.5.0
## [75] later_1.3.1 splines_4.3.0
## [77] BiocFileCache_2.9.0 lattice_0.21-8
## [79] rtracklayer_1.61.0 bit_4.0.5
## [81] tidyselect_1.2.0 Biostrings_2.69.0
## [83] knitr_1.42 bookdown_0.33
## [85] xfun_0.39 stringi_1.7.12
## [87] yaml_2.3.7 evaluate_0.20
## [89] codetools_0.2-19 tibble_3.2.1
## [91] BiocManager_1.30.20 cli_3.6.1
## [93] xtable_1.8-4 munsell_0.5.0
## [95] jquerylib_0.1.4 Rcpp_1.0.10
## [97] GenomicDataCommons_1.25.0 dbplyr_2.3.2
## [99] png_0.1-8 XML_3.99-0.14
## [101] parallel_4.3.0 ellipsis_0.3.2
## [103] readr_2.1.4 blob_1.2.4
## [105] prettyunits_1.1.1 bitops_1.0-7
## [107] scales_1.2.1 purrr_1.0.1
## [109] crayon_1.5.2 rlang_1.1.1
## [111] KEGGREST_1.41.0 rvest_1.0.3
## [113] formatR_1.14