We report a novel method to identify specific transcriptomic phenotypes based on an elitist non-dominated sorting genetic algorithm that combines the advantages of clustering methods and the exploratory properties of genetic algorithms to discover biologically and clinically relevant molecular subtypes in different cancers.
In the new era of omics data, precision medicine has become the new paradigm of cancer treatment. Among all available omics techniques, gene expression profiling, in particular, has been increasingly used to classify tumor subtypes with different biological behavior. Cancer subtype discovery is usually approached from two possible perspectives:
-Using the molecular data alone with unsupervised techniques such as clustering analysis. -Using supervised techniques focusing entirely on survival data.
The problem of finding patients subgroups with survival differences while maintaining cluster consistency could be viewed as a bi-objective problem, where there is a trade-off between the separability of the different groups and the ability of a given signature to consistently distinguish patients with different clinical outcomes. This gives rise to a set of optimal solutions, also known as Pareto-optimal solutions. To overcome these issues, we combined the advantages of clustering methods for grouping heterogeneous omics data and the search properties of genetic algorithms in GSgalgoR: A flexible yet robust multi-objective meta-heuristic for disease subtype discovery based on an elitist non-dominated sorting genetic algorithm (NSGA-II), driven by the underlying premise of maximizing survival differences between groups while getting high consistency and robustness of the clusters obtained.
In the GSgalgoR package, the NSGA-II framework was used for finding multiple Pareto-optimal solutions to classify patients according to their gene expression patterns. Basically, NSGA-II starts with a population of competing individuals which are evaluated under a set of fitness functions that estimate the survival differences and cohesiveness of the different transcriptomic groups. Then, solutions are ranked and sorted according to their non-domination level which will affect the way they are chosen to be submitted to the so-called “evolutionary operators” such as crossover and mutation. Once a set of well-suited solutions are selected and reproduced, a new offspring of individuals composed of a mixture of the “genetic information” of the parents is obtained. Parents and offspring are pooled and the best-ranked solutions are selected and passed to the next generation which will start over the same process again.
To install GSgalgoR package, start R and enter:
if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("GSgalgoR")
library(GSgalgoR)Alternatively you can install GSgalgoR from github using the devtool package
devtools::install_github("https://github.com/harpomaxx/GSgalgoR")
library(GSgalgoR)To standardize the structure of genomic data, we use the
ExpressionSet
structure for the examples given in this guide. The ExpressionSet objects are
formed mainly by:
To start testing GSgalgoR, we will use two Breast Cancer datasets. Namely, the UPP and the TRANSBIG datasets. Additionally, we will use PAM50 centroids to perform breast cancer sample classification. The datasets can be accessed from the following Bioconductor packages:
BiocManager::install("breastCancerUPP",version = "devel")
BiocManager::install("breastCancerTRANSBIG",version = "devel")
library(breastCancerTRANSBIG)
library(breastCancerUPP)Also, some basic packages are needed to run the example in this vignette
library(GSgalgoR)
library(Biobase)
library(genefu)
library(survival)
library(survminer)
library(ggplot2)
data(pam50)To access the ExpressionSets we use:
data(upp)
Train<- upp
rm(upp)
data(transbig)
Test<- transbig
rm(transbig)
#To access gene expression data
train_expr<- exprs(Train)
test_expr<- exprs(Test)
#To access feature data
train_features<- fData(Train)
test_features<- fData(Test)
#To access clinical data
train_clinic <- pData(Train) 
test_clinic <- pData(Test) Galgo can accept any numeric data, like probe intensity from microarray experiments or RNAseq normalized counts, nevertheless, features are expected to be scaled across the dataset before being plugged in into the Galgo Framework. For PAM50 classification, Gene Symbols are expected, so probesets are mapped into their respective gene symbols. Probesets mapping for multiple genes are expanded while Genes mapped to multiple probes are collapsed selecting the probes with the highest variance for each duplicated gene.
#Custom function to drop duplicated genes (keep genes with highest variance)
DropDuplicates<- function(eset, map= "Gene.symbol"){
    #Drop NA's
    drop <- which(is.na(fData(eset)[,map]))
    eset <- eset[-drop,]
    #Drop duplicates
    drop <- NULL
    Dup <- as.character(unique(fData(eset)[which(duplicated
            (fData(eset)[,map])),map]))
    Var <- apply(exprs(eset),1,var)
    for(j in Dup){
        pos <- which(fData(eset)[,map]==j)
        drop <- c(drop,pos[-which.max(Var[pos])])
    }
    eset <- eset[-drop,]
    featureNames(eset) <- fData(eset)[,map]
    return(eset)
}
# Custom function to expand probesets mapping to multiple genes
expandProbesets <- function (eset, sep = "///", map="Gene.symbol"){
    x <- lapply(featureNames(eset), function(x) strsplit(x, sep)[[1]])
    y<- lapply(as.character(fData(eset)[,map]), function(x) strsplit(x,sep))
    eset <- eset[order(sapply(x, length)), ]
    x <- lapply(featureNames(eset), function(x) strsplit(x, sep)[[1]])
    y<- lapply(as.character(fData(eset)[,map]), function(x) strsplit(x,sep))
    idx <- unlist(sapply(1:length(x), function(i) rep(i,length(x[[i]]))))
    idy <- unlist(sapply(1:length(y), function(i) rep(i,length(y[[i]]))))
    xx <- !duplicated(unlist(x))
    idx <- idx[xx]
    idy <- idy[xx]
    x <- unlist(x)[xx]
    y <- unlist(y)[xx]
    eset <- eset[idx, ]
    featureNames(eset) <- x
    fData(eset)[,map] <- x
    fData(eset)$gene <- y
    return(eset)
}Train=DropDuplicates(Train)
Train=expandProbesets(Train)
#Drop NAs in survival
Train <- Train[,!is.na(
    survival::Surv(time=pData(Train)$t.rfs,event=pData(Train)$e.rfs))] 
Test=DropDuplicates(Test)
Test=expandProbesets(Test)
#Drop NAs in survival
Test <- 
    Test[,!is.na(survival::Surv(
        time=pData(Test)$t.rfs,event=pData(Test)$e.rfs))] 
#Determine common probes (Genes)
Int= intersect(rownames(Train),rownames(Test))
Train= Train[Int,]
Test= Test[Int,]
identical(rownames(Train),rownames(Test))
#> [1] TRUEFor simplicity and speed, we will create a reduced expression matrix for the examples.
#First we will get PAM50 centroids from genefu package
PAM50Centroids <- pam50$centroids
PAM50Genes <- pam50$centroids.map$probe
PAM50Genes<- featureNames(Train)[ featureNames(Train) %in% PAM50Genes]
#Now we sample 200 random genes from expression matrix
Non_PAM50Genes<- featureNames(Train)[ !featureNames(Train) %in% PAM50Genes]
Non_PAM50Genes <- sample(Non_PAM50Genes,200, replace=FALSE)
reduced_set <- c(PAM50Genes, Non_PAM50Genes)
#Now we get the reduced training and test sets
Train<- Train[reduced_set,]
Test<- Test[reduced_set,]Apply robust linear scaling as proposed in paper reference
exprs(Train) <- t(apply(exprs(Train),1,genefu::rescale,na.rm=TRUE,q=0.05))
exprs(Test) <- t(apply(exprs(Test),1,genefu::rescale,na.rm=TRUE,q=0.05))
train_expr <- exprs(Train)
test_expr <- exprs(Test)The ‘Surv’ object is created by the Surv() function of the survival package.
This uses phenotypic data that are contained in the corresponding datasets,
accessed by the pData command.
train_clinic <- pData(Train) 
test_clinic <- pData(Test)
train_surv <- survival::Surv(time=train_clinic$t.rfs,event=train_clinic$e.rfs)
test_surv <- survival::Surv(time=test_clinic$t.rfs,event=test_clinic$e.rfs)The main function in this package is galgo(). It accepts an expression matrix
and survival object to find robust gene expression signatures related to a given
outcome. This function contains some parameters that can be modified, according
to the characteristics of the analysis to be performed.
The principal parameters are:
# For testing reasons it is set to a low number but ideally should be above 100
population <- 30 
# For testing reasons it is set to a low number but ideally should be above 150
generations <-15
nCV <- 5                      
distancetype <- "pearson"     
TournamentSize <- 2
period <- 3650set.seed(264)
output <- GSgalgoR::galgo(generations = generations, 
                        population = population, 
                        prob_matrix = train_expr, 
                        OS = train_surv,
                        nCV = nCV, 
                        distancetype = distancetype,
                        TournamentSize = TournamentSize, 
                        period = period)print(class(output))
#> [1] "galgo.Obj"
#> attr(,"package")
#> [1] "GSgalgoR"The output of the galgo() function is an object of type galgo.Obj that has two
slots with the elements:
Is a l x (n + 5) matrix where n is the number of features evaluated and l is the number of solutions obtained.
Is a list of length equal to the number of generations run in the algorithm. Each element is a l x 2 matrix where l is the number of solutions obtained and the columns are the SC Fitness and the Survival Fitness values respectively.
For easier interpretation of the galgo.Obj, the output can be transformed to a
list or to a data.frame objects.
This function restructurates a galgo.Obj to a more easy to understand an use
list. This output is particularly useful if one wants to select a given solution
and use its outputs in a new classifier. The output of type list has a length
equals to the number of solutions obtained by the galgo algorithm.
Basically this output is a list of lists, where each element of the output is named after the solution’s name (solution.n, where n is the number assigned to that solution), and inside of it, it has all the constituents for that given solution with the following structure:
outputList <- to_list(output)
head(names(outputList))
#> [1] "Solution.1" "Solution.2" "Solution.3" "Solution.4" "Solution.5"
#> [6] "Solution.6"To evaluate the structure of the first solution we can run:
outputList[["Solution.1"]]
#> $Genes
#>  [1] "RRM2"      "BIRC5"     "MYC"       "CDC20"     "MAPT"      "ACTR3B"   
#>  [7] "EFR3A"     "HTR7P1"    "RAB6A"     "CHAT"      "SNX11"     "ZGPAT"    
#> [13] "AKAP7"     "SIRT5"     "SAMD4A"    "CCNA2"     "SMPX"      "SLCO1B3"  
#> [19] "MAGEL2"    "INTS1"     "LRIT1"     "CUL5"      "OSTF1"     "SFTPA2"   
#> [25] "RNFT2"     "POM121L4P" "FUT9"      "ADAM23"    "ZBTB40"    "LDHB"     
#> [31] "MEA1"      "HABP2"     "RPL21"     "HMGN5"     "UBA1"      "TNMD"     
#> [37] "EIF3H"    
#> 
#> $k
#> [1] 7
#> 
#> $SC.Fit
#> [1] 0.0366157
#> 
#> $Surv.Fit
#> [1] 924.6521
#> 
#> $rank
#> [1] 1
#> 
#> $CrowD
#> [1] InfThe current function restructures a galgo.Obj to a more easy to understand
an use data.frame. The output data frame has m x n dimensions, were the
rownames (m) are the solutions obtained by the galgo algorithm. The columns has
the following structure:
outputDF <- to_dataframe(output)
head(outputDF)
#>                    Genes k     SC.Fit Surv.Fit Rank     CrowD
#> Solutions.1 RRM2, BI.... 7 0.03661570 924.6521    1       Inf
#> Solutions.2 MYBL2, B.... 3 0.17120138 192.2896    1       Inf
#> Solutions.3 RRM2, EG.... 2 0.10729955 688.8392    1 0.6748982
#> Solutions.4 BIRC5, C.... 5 0.05379444 807.6729    1 0.6745774
#> Solutions.5 BIRC5, M.... 2 0.16893797 518.2101    1 0.4648913
#> Solutions.6 BIRC5, C.... 2 0.13804813 641.3836    1 0.3006677Once we obtain the galgo.obj from the output of galgo() we can plot the
obtained Pareto front and see how it evolved trough the tested number of
generations
plot_pareto(output)Breast cancer (BRCA) is the most common neoplasm in women to date and one of the best studied cnacer types. Currently, numerous molecular alteration for this type of cancer are well known and many transcriptomic signatures have been developed for this type of cancer. In this regards, Perou et al. proposed one of the first molecular subtype classification according to transcriptomic profiles of the tumor, which recapitulates naturally-occurring gene expression patterns that encompass different functional pathways and patient outcomes. These subtypes, (LumA, LumB, Basal-like, HER2 and Normal-Like) have a strong overlap with the classical histopathological classification of BRCA tumors and might affect decision making when used to decided chemotherapy in certain cases.
To evaluate Galgo’s performance along with PAM50 classification, we will use the two already scaled and reduced BRCA gene expression datasets and will compare Galgo performance with the widely used intrinsic molecular subtype PAM50 classification. Galgo performs feature selection by design, so this step is not strictly necessary to use galgoR (although feature selection might fasten GSgalgoRruns), nevertheless, appropriate gene expression scaling is critical when running GSgalgoR.
The scaled expression values of each patient are compared with the prototypical centroids using Pearson’s correlation coefficient and the closest centroid to each patient is used to assign the corresponding labels.
#The reduced UPP dataset will be used as training set 
train_expression <- exprs(Train) 
train_clinic<- pData(Train)
train_features<- fData(Train)
train_surv<- survival::Surv(time=train_clinic$t.rfs,event=train_clinic$e.rfs)
#The reduced TRANSBIG dataset will be used as test set 
test_expression <- exprs(Test) 
test_clinic<- pData(Test)
test_features<- fData(Test)
test_surv<- survival::Surv(time=test_clinic$t.rfs,event=test_clinic$e.rfs)
#PAM50 centroids
centroids<- pam50$centroids
#Extract features from both data.frames
inBoth<- Reduce(intersect, list(rownames(train_expression),rownames(centroids)))
#Classify samples 
PAM50_train<- cluster_classify(train_expression[inBoth,],centroids[inBoth,],
                            method = "spearman")
table(PAM50_train)
#> PAM50_train
#>  1  2  3  4  5 
#> 22 30 94 73 15
PAM50_test<- cluster_classify(test_expression[inBoth,],centroids[inBoth,],
                            method = "spearman")
table(PAM50_test)
#> PAM50_test
#>  1  2  3  4  5 
#> 45 26 80 44  3
# Classify samples using genefu
#annot<- fData(Train)
#colnames(annot)[3]="Gene.Symbol"
#PAM50_train<- molecular.subtyping(sbt.model = "pam50",
#         data = t(train_expression), annot = annot,do.mapping = TRUE)Once the patients are classified according to their closest centroids, we can now evaluate the survival curves for the different types in each of the datasets
surv_formula <- 
    as.formula("Surv(train_clinic$t.rfs,train_clinic$e.rfs)~ PAM50_train")
tumortotal1 <- surv_fit(surv_formula,data=train_clinic)
tumortotal1diff <- survdiff(surv_formula)
tumortotal1pval<- pchisq(tumortotal1diff$chisq, length(tumortotal1diff$n) - 1,
                         lower.tail = FALSE) 
p<-ggsurvplot(tumortotal1,
            data=train_clinic,
            risk.table=TRUE,
            pval=TRUE,
            palette="dark2",
            title="UPP breast cancer \n PAM50 subtypes survival",
            surv.scale="percent",
            conf.int=FALSE, 
            xlab="time (days)", 
            ylab="survival(%)", 
            xlim=c(0,3650),
            break.time.by = 365, 
            ggtheme = theme_minimal(), 
            risk.table.y.text.col = TRUE, 
            risk.table.y.text = FALSE,censor=FALSE)
print(p)surv_formula <- 
    as.formula("Surv(test_clinic$t.rfs,test_clinic$e.rfs)~ PAM50_test")
tumortotal2 <- surv_fit(surv_formula,data=test_clinic)
tumortotal2diff <- survdiff(surv_formula)
tumortotal2pval<- pchisq(tumortotal2diff$chisq, length(tumortotal2diff$n) - 1,
                        lower.tail = FALSE) 
p<-ggsurvplot(tumortotal2,
            data=test_clinic,
            risk.table=TRUE,
            pval=TRUE,
            palette="dark2",
            title="TRANSBIG breast cancer \n PAM50 subtypes survival",
            surv.scale="percent",
            conf.int=FALSE,
            xlab="time (days)",
            ylab="survival(%)",
            xlim=c(0,3650),
            break.time.by = 365,
            ggtheme = theme_minimal(),
            risk.table.y.text.col = TRUE,
            risk.table.y.text = FALSE,
            censor=FALSE)
print(p)Now we run Galgo to find cohesive and clinically meaningful signatures for BRCA using UPP data as training set and TRANSBIG data as test set
population <- 15             
generations <-5             
nCV <- 5                      
distancetype <- "pearson"     
TournamentSize <- 2
period <- 3650Run Galgo on the training set
output= GSgalgoR::galgo(generations = generations,
                    population = population,
                    prob_matrix = train_expression,
                    OS=train_surv,
                    nCV= nCV, 
                    distancetype=distancetype,
                    TournamentSize=TournamentSize,
                    period=period)
print(class(output))plot_pareto(output)
output_df<- to_dataframe(output)
NonDom_solutions<- output_df[output_df$Rank==1,]
# N of non-dominated solutions 
nrow(NonDom_solutions)
#> [1] 5
# N of partitions found
table(NonDom_solutions$k)
#> 
#>  3  6  8 10 
#>  2  1  1  1
#Average N of genes per signature
mean(unlist(lapply(NonDom_solutions$Genes,length)))
#> [1] 77.8
#SC range
range(NonDom_solutions$SC.Fit)
#> [1] 0.01534911 0.12736885
# Survival fitnesss range
range(NonDom_solutions$Surv.Fit)
#> [1] 447.8611 919.9274Now we select the best performing solutions for each number of partitions (k) according to C.Index
RESULT<- non_dominated_summary(output=output,
                            OS=train_surv, 
                            prob_matrix= train_expression,
                            distancetype =distancetype 
                            )
best_sol=NULL
for(i in unique(RESULT$k)){
    best_sol=c(
    best_sol,
    RESULT[RESULT$k==i,"solution"][which.max(RESULT[RESULT$k==i,"C.Index"])])
}
print(best_sol)
#> [1] "Solutions.1" "Solutions.5" "Solutions.3" "Solutions.4"Now we create the prototypic centroids of the selected solutions
CentroidsList <- create_centroids(output, 
                                solution_names = best_sol,
                                trainset = train_expression)We will test the Galgo signatures found with the UPP training set in an independent test set (TRANSBIG)
train_classes<- classify_multiple(prob_matrix=train_expression,
                                centroid_list= CentroidsList, 
                                distancetype = distancetype)
test_classes<- classify_multiple(prob_matrix=test_expression,
                                centroid_list= CentroidsList, 
                                distancetype = distancetype)To calculate the train and test C.Index, the risk coefficients are calculated for each subclass in the training set and then are used to predict the risk of the different groups in the test set. This is particularly important for signatures with high number of partitions, were the survival differences of different groups might overlap and change their relative order, which is of great importance in the C.Index calculation.
Prediction.models<- list()
for(i in best_sol){
    OS<- train_surv
    predicted_class<- as.factor(train_classes[,i])
    predicted_classdf <- as.data.frame(predicted_class)
    colnames(predicted_classdf)<- i
    surv_formula <- as.formula(paste0("OS~ ",i))
    coxsimple=coxph(surv_formula,data=predicted_classdf)
    Prediction.models[[i]]<- coxsimple
}
C.indexes<- data.frame(train_CI=rep(NA,length(best_sol)),
                    test_CI=rep(NA,length(best_sol)))
rownames(C.indexes)<- best_sol
for(i in best_sol){
    predicted_class_train<- as.factor(train_classes[,i])
    predicted_class_train_df <- as.data.frame(predicted_class_train)
    colnames(predicted_class_train_df)<- i
    CI_train<- 
        concordance.index(predict(Prediction.models[[i]],
                                predicted_class_train_df),
                                surv.time=train_surv[,1],
                                surv.event=train_surv[,2],
                                outx=FALSE)$c.index
    C.indexes[i,"train_CI"]<- CI_train
    predicted_class_test<- as.factor(test_classes[,i])
    predicted_class_test_df <- as.data.frame(predicted_class_test)
    colnames(predicted_class_test_df)<- i
    CI_test<- 
        concordance.index(predict(Prediction.models[[i]],
                                predicted_class_test_df),
                                surv.time=test_surv[,1],
                                surv.event=test_surv[,2],
                                outx=FALSE)$c.index
    C.indexes[i,"test_CI"]<- CI_test
    }
print(C.indexes)
#>              train_CI   test_CI
#> Solutions.1 0.6494606 0.5879290
#> Solutions.5 0.6059079 0.5595661
#> Solutions.3 0.6125618 0.5335700
#> Solutions.4 0.6423027 0.5377909
best_signature<- best_sol[which.max(C.indexes$test_CI)]
print(best_signature)
#> [1] "Solutions.1"We test best galgo signature with training and test sets
train_class <- train_classes[,best_signature]
surv_formula <- 
    as.formula("Surv(train_clinic$t.rfs,train_clinic$e.rfs)~ train_class")
tumortotal1 <- surv_fit(surv_formula,data=train_clinic)
tumortotal1diff <- survdiff(surv_formula)
tumortotal1pval<- pchisq(tumortotal1diff$chisq,
                        length(tumortotal1diff$n) - 1,
                        lower.tail = FALSE) 
p<-ggsurvplot(tumortotal1,
            data=train_clinic,
            risk.table=TRUE,pval=TRUE,palette="dark2",
            title="UPP breast cancer \n Galgo subtypes survival",
            surv.scale="percent",
            conf.int=FALSE, xlab="time (days)", 
            ylab="survival(%)", xlim=c(0,3650),
            break.time.by = 365,
            ggtheme = theme_minimal(), 
            risk.table.y.text.col = TRUE, 
            risk.table.y.text = FALSE,censor=FALSE)
print(p)
test_class <- test_classes[,best_signature]
surv_formula <- 
    as.formula("Surv(test_clinic$t.rfs,test_clinic$e.rfs)~ test_class")
tumortotal1 <- surv_fit(surv_formula,data=test_clinic)
tumortotal1diff <- survdiff(surv_formula)
tumortotal1pval<- pchisq(tumortotal1diff$chisq,
                        length(tumortotal1diff$n) - 1,
                        lower.tail = FALSE) 
p<-ggsurvplot(tumortotal1,
            data=test_clinic,
            risk.table=TRUE,
            pval=TRUE,palette="dark2",
            title="TRANSBIG breast cancer \n Galgo subtypes survival",
            surv.scale="percent",
            conf.int=FALSE, 
            xlab="time (days)",
            ylab="survival(%)",
            xlim=c(0,3650),
            break.time.by = 365, 
            ggtheme = theme_minimal(), 
            risk.table.y.text.col = TRUE,
            risk.table.y.text = FALSE,
            censor=FALSE)
print(p)Compare PAM50 classification vs Galgo classification in the TRANSBIG (test) dataset
surv_formula1 <- 
    as.formula("Surv(test_clinic$t.rfs,test_clinic$e.rfs)~ test_class")
tumortotal1 <- surv_fit(surv_formula1,data=test_clinic)
tumortotal1diff <- survdiff(surv_formula1)
tumortotal1pval<- pchisq(tumortotal1diff$chisq,
                        length(tumortotal1diff$n) - 1,
                        lower.tail = FALSE) 
surv_formula2 <- 
    as.formula("Surv(test_clinic$t.rfs,test_clinic$e.rfs)~ PAM50_test")
tumortotal2 <- surv_fit(surv_formula2,data=test_clinic)
tumortotal2diff <- survdiff(surv_formula2)
tumortotal2pval<- pchisq(tumortotal1diff$chisq,
                        length(tumortotal2diff$n) - 1,
                        lower.tail = FALSE) 
SURV=list(GALGO=tumortotal1,PAM50=tumortotal2 )
COLS=c(1:8,10)
par(cex=1.35, mar=c(3.8, 3.8, 2.5, 2.5) + 0.1)
p=ggsurvplot(SURV,
            combine=TRUE,
            data=test_clinic,
            risk.table=TRUE,
            pval=TRUE,
            palette="dark2",
            title="Galgo vs. PAM50 subtypes \n BRCA survival comparison",
            surv.scale="percent",
            conf.int=FALSE,
            xlab="time (days)",
            ylab="survival(%)",
            xlim=c(0,period),
            break.time.by = 365, 
            ggtheme = theme_minimal(),
            risk.table.y.text.col = TRUE,
            risk.table.y.text = FALSE,
            censor=FALSE)
print(p)sessionInfo()
#> R version 4.4.0 beta (2024-04-15 r86425)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.4 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_US.UTF-8        LC_COLLATE=en_US.UTF-8    
#>  [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] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#>  [1] survminer_0.4.9             ggpubr_0.6.0               
#>  [3] ggplot2_3.5.1               genefu_2.36.0              
#>  [5] AIMS_1.36.0                 e1071_1.7-14               
#>  [7] iC10_1.5                    iC10TrainingData_1.3.1     
#>  [9] impute_1.78.0               pamr_1.56.2                
#> [11] cluster_2.1.6               biomaRt_2.60.0             
#> [13] survcomp_1.54.0             prodlim_2023.08.28         
#> [15] survival_3.6-4              Biobase_2.64.0             
#> [17] BiocGenerics_0.50.0         breastCancerUPP_1.41.0     
#> [19] breastCancerTRANSBIG_1.41.0 GSgalgoR_1.14.0            
#> [21] BiocStyle_2.32.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] jsonlite_1.8.8          magrittr_2.0.3          magick_2.8.3           
#>   [4] SuppDists_1.1-9.7       farver_2.1.1            rmarkdown_2.26         
#>   [7] zlibbioc_1.50.0         vctrs_0.6.5             memoise_2.0.1          
#>  [10] tinytex_0.50            rstatix_0.7.2           htmltools_0.5.8.1      
#>  [13] progress_1.2.3          curl_5.2.1              broom_1.0.5            
#>  [16] sass_0.4.9              parallelly_1.37.1       KernSmooth_2.23-22     
#>  [19] bslib_0.7.0             httr2_1.0.1             zoo_1.8-12             
#>  [22] cachem_1.0.8            commonmark_1.9.1        lifecycle_1.0.4        
#>  [25] iterators_1.0.14        pkgconfig_2.0.3         Matrix_1.7-0           
#>  [28] R6_2.5.1                fastmap_1.1.1           GenomeInfoDbData_1.2.12
#>  [31] future_1.33.2           digest_0.6.35           colorspace_2.1-0       
#>  [34] AnnotationDbi_1.66.0    S4Vectors_0.42.0        nsga2R_1.1             
#>  [37] RSQLite_2.3.6           labeling_0.4.3          filelock_1.0.3         
#>  [40] km.ci_0.5-6             fansi_1.0.6             httr_1.4.7             
#>  [43] abind_1.4-5             compiler_4.4.0          proxy_0.4-27           
#>  [46] bit64_4.0.5             withr_3.0.0             doParallel_1.0.17      
#>  [49] backports_1.4.1         carData_3.0-5           DBI_1.2.2              
#>  [52] highr_0.10              ggsignif_0.6.4          lava_1.8.0             
#>  [55] rappdirs_0.3.3          tools_4.4.0             future.apply_1.11.2    
#>  [58] bootstrap_2019.6        glue_1.7.0              gridtext_0.1.5         
#>  [61] grid_4.4.0              generics_0.1.3          gtable_0.3.5           
#>  [64] KMsurv_0.1-5            class_7.3-22            tidyr_1.3.1            
#>  [67] data.table_1.15.4       hms_1.1.3               xml2_1.3.6             
#>  [70] car_3.1-2               utf8_1.2.4              XVector_0.44.0         
#>  [73] markdown_1.12           foreach_1.5.2           pillar_1.9.0           
#>  [76] stringr_1.5.1           limma_3.60.0            splines_4.4.0          
#>  [79] ggtext_0.1.2            dplyr_1.1.4             BiocFileCache_2.12.0   
#>  [82] lattice_0.22-6          bit_4.0.5               tidyselect_1.2.1       
#>  [85] Biostrings_2.72.0       knitr_1.46              gridExtra_2.3          
#>  [88] bookdown_0.39           IRanges_2.38.0          stats4_4.4.0           
#>  [91] xfun_0.43               statmod_1.5.0           stringi_1.8.3          
#>  [94] UCSC.utils_1.0.0        yaml_2.3.8              evaluate_0.23          
#>  [97] codetools_0.2-20        tibble_3.2.1            BiocManager_1.30.22    
#> [100] cli_3.6.2               survivalROC_1.0.3.1     xtable_1.8-4           
#> [103] munsell_0.5.1           jquerylib_0.1.4         survMisc_0.5.6         
#> [106] Rcpp_1.0.12             GenomeInfoDb_1.40.0     rmeta_3.0              
#> [109] globals_0.16.3          dbplyr_2.5.0            png_0.1-8              
#> [112] parallel_4.4.0          mco_1.16                blob_1.2.4             
#> [115] prettyunits_1.2.0       mclust_6.1.1            listenv_0.9.1          
#> [118] scales_1.3.0            purrr_1.0.2             crayon_1.5.2           
#> [121] rlang_1.1.3             KEGGREST_1.44.0