# To install scGPS from github (Depending on the configuration of the local
# computer or HPC, possible custom C++ compilation may be required - see
# installation trouble-shootings below)
devtools::install_github("IMB-Computational-Genomics-Lab/scGPS")
# for C++ compilation trouble-shooting, manual download and installation can be
# done from github
git clone https://github.com/IMB-Computational-Genomics-Lab/scGPS
# then check in scGPS/src if any of the precompiled (e.g.  those with *.so and
# *.o) files exist and delete them before recompiling
# then with the scGPS as the R working directory, manually install and load
# using devtools functionality
# Install the package
devtools::install()
#load the package to the workspace 
library(scGPS)The purpose of this workflow is to solve the following task:
# load mixed population 1 (loaded from day_2_cardio_cell_sample dataset, 
# named it as day2)
library(scGPS)
day2 <- day_2_cardio_cell_sample
mixedpop1 <- new_scGPS_object(ExpressionMatrix = day2$dat2_counts,
    GeneMetadata = day2$dat2geneInfo, CellMetadata = day2$dat2_clusters)
# load mixed population 2 (loaded from day_5_cardio_cell_sample dataset, 
# named it as day5)
day5 <- day_5_cardio_cell_sample
mixedpop2 <- new_scGPS_object(ExpressionMatrix = day5$dat5_counts,
    GeneMetadata = day5$dat5geneInfo, CellMetadata = day5$dat5_clusters)
# select a subpopulation
c_selectID <- 1
# load gene list (this can be any lists of user selected genes)
genes <- training_gene_sample
genes <- genes$Merged_unique
# load cluster information 
cluster_mixedpop1 <- colData(mixedpop1)[,1]
cluster_mixedpop2 <- colData(mixedpop2)[,1]
#run training (running nboots = 3 here, but recommend to use nboots = 50-100)
LSOLDA_dat <- bootstrap_prediction(nboots = 3, mixedpop1 = mixedpop1, 
    mixedpop2 = mixedpop2, genes = genes, c_selectID  = c_selectID,
    listData = list(), cluster_mixedpop1 = cluster_mixedpop1, 
    cluster_mixedpop2 = cluster_mixedpop2, trainset_ratio = 0.7)
names(LSOLDA_dat)
#> [1] "Accuracy"          "ElasticNetGenes"   "Deviance"         
#> [4] "ElasticNetFit"     "LDAFit"            "predictor_S1"     
#> [7] "ElasticNetPredict" "LDAPredict"        "cell_results"# summary results LDA
sum_pred_lda <- summary_prediction_lda(LSOLDA_dat = LSOLDA_dat, nPredSubpop = 4)
# summary results Lasso to show the percent of cells
# classified as cells belonging 
sum_pred_lasso <- summary_prediction_lasso(LSOLDA_dat = LSOLDA_dat,
    nPredSubpop = 4)
# plot summary results 
plot_sum <-function(sum_dat){
    sum_dat_tf <- t(sum_dat)
    sum_dat_tf <- na.omit(sum_dat_tf)
    sum_dat_tf <- apply(sum_dat[, -ncol(sum_dat)],1,
        function(x){as.numeric(as.vector(x))})
    sum_dat$names <- gsub("ElasticNet for subpop","sp",  sum_dat$names )
    sum_dat$names <- gsub("in target mixedpop","in p",  sum_dat$names) 
    sum_dat$names <- gsub("LDA for subpop","sp",  sum_dat$names )
    sum_dat$names <- gsub("in target mixedpop","in p",  sum_dat$names)
    colnames(sum_dat_tf) <- sum_dat$names
    boxplot(sum_dat_tf, las=2)
}
plot_sum(sum_pred_lasso)# summary accuracy to check the model accuracy in the leave-out test set 
summary_accuracy(object = LSOLDA_dat)
#> [1] 68.39623 64.59330 64.95327
# summary maximum deviance explained by the model 
summary_deviance(object = LSOLDA_dat)
#> $allDeviance
#> [1] "6.69" "6.56" "8.57"
#> 
#> $DeviMax
#>          dat_DE$Dfd          Deviance           DEgenes
#> 1                 0              8.57    genes_cluster1
#> 2                 1              8.57    genes_cluster1
#> 3                 3              8.57    genes_cluster1
#> 4 remaining DEgenes remaining DEgenes remaining DEgenes
#> 
#> $LassoGenesMax
#> NULLThe purpose of this workflow is to solve the following task:
(skip this step if clusters are known)
# find clustering information in an expresion data using CORE
day5 <- day_5_cardio_cell_sample
cellnames <- colnames(day5$dat5_counts)
cluster <-day5$dat5_clusters
cellnames <-data.frame("Cluster"=cluster, "cellBarcodes" = cellnames)
mixedpop2 <-new_scGPS_object(ExpressionMatrix = day5$dat5_counts,
                    GeneMetadata = day5$dat5geneInfo, CellMetadata = cellnames)
CORE_cluster <- CORE_clustering(mixedpop2, remove_outlier = c(0), PCA=FALSE)
# to update the clustering information, users can ...
key_height <- CORE_cluster$optimalClust$KeyStats$Height
optimal_res <- CORE_cluster$optimalClust$OptimalRes
optimal_index = which(key_height == optimal_res)
clustering_after_outlier_removal <- unname(unlist(
 CORE_cluster$Cluster[[optimal_index]]))
corresponding_cells_after_outlier_removal <- CORE_cluster$cellsForClustering
original_cells_before_removal <- colData(mixedpop2)[,2]
corresponding_index <- match(corresponding_cells_after_outlier_removal,
                            original_cells_before_removal )
# check the matching
identical(as.character(original_cells_before_removal[corresponding_index]),
         corresponding_cells_after_outlier_removal)
#> [1] TRUE
# create new object with the new clustering after removing outliers
mixedpop2_post_clustering <- mixedpop2[,corresponding_index]
colData(mixedpop2_post_clustering)[,1] <- clustering_after_outlier_removal(skip this step if clusters are known)
(SCORE aims to get stable subpopulation results by introducing bagging aggregation and bootstrapping to the CORE algorithm)
# find clustering information in an expresion data using SCORE
day5 <- day_5_cardio_cell_sample
cellnames <- colnames(day5$dat5_counts)
cluster <-day5$dat5_clusters
cellnames <-data.frame("Cluster"=cluster, "cellBarcodes" = cellnames)
mixedpop2 <-new_scGPS_object(ExpressionMatrix = day5$dat5_counts,
                    GeneMetadata = day5$dat5geneInfo, CellMetadata = cellnames )
SCORE_test <- CORE_bagging(mixedpop2, remove_outlier = c(0), PCA=FALSE,
                                bagging_run = 20, subsample_proportion = .8)dev.off()
#> null device 
#>           1
##3.2.1 plot CORE clustering
p1 <- plot_CORE(CORE_cluster$tree, CORE_cluster$Cluster, 
    color_branch = c("#208eb7", "#6ce9d3", "#1c5e39", "#8fca40", "#154975",
        "#b1c8eb"))
p1
#> $mar
#> [1] 1 5 0 1
#extract optimal index identified by CORE
key_height <- CORE_cluster$optimalClust$KeyStats$Height
optimal_res <- CORE_cluster$optimalClust$OptimalRes
optimal_index = which(key_height == optimal_res)
#plot one optimal clustering bar
plot_optimal_CORE(original_tree= CORE_cluster$tree,
                 optimal_cluster = unlist(CORE_cluster$Cluster[optimal_index]),
                 shift = -2000)
#> Ordering and assigning labels...
#> 2
#> 162335NA
#> 3
#> 162335423
#> Plotting the colored dendrogram now....
#> Plotting the bar underneath now....
##3.2.2 plot SCORE clustering
#plot all clustering bars
plot_CORE(SCORE_test$tree, list_clusters = SCORE_test$Cluster)
#plot one stable optimal clustering bar
plot_optimal_CORE(original_tree= SCORE_test$tree,
                 optimal_cluster = unlist(SCORE_test$Cluster[
                    SCORE_test$optimal_index]),
                 shift = -100)
#> Ordering and assigning labels...
#> 2
#> 24112NANANANANA
#> 3
#> 24112250NANANANA
#> 4
#> 24112250335NANANA
#> 5
#> 24112250335367NANA
#> 6
#> 24112250335367414NA
#> 7
#> 24112250335367414470
#> Plotting the colored dendrogram now....
#> Plotting the bar underneath now....t <- tSNE(expression.mat=assay(mixedpop2))
#> Preparing PCA inputs using the top 1500 genes ...
#> Computing PCA values...
#> Running tSNE ...
p2 <-plot_reduced(t, color_fac = factor(colData(mixedpop2)[,1]),
                      palletes =1:length(unique(colData(mixedpop2)[,1])))
#> Warning: The dot-dot notation (`..count..`) was deprecated in ggplot2 3.4.0.
#> ℹ Please use `after_stat(count)` instead.
#> ℹ The deprecated feature was likely used in the cowplot package.
#>   Please report the issue at <https://github.com/wilkelab/cowplot/issues>.
#> This warning is displayed once every 8 hours.
#> Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
#> generated.
#> Warning: Use of `reduced_dat_toPlot$Dim1` is discouraged.
#> ℹ Use `Dim1` instead.
#> Warning: Use of `reduced_dat_toPlot$Dim2` is discouraged.
#> ℹ Use `Dim2` instead.
p2#load gene list (this can be any lists of user-selected genes)
genes <-training_gene_sample
genes <-genes$Merged_unique
#the gene list can also be objectively identified by differential expression
#analysis cluster information is requied for find_markers. Here, we use
#CORE results.
#colData(mixedpop2)[,1] <- unlist(SCORE_test$Cluster[SCORE_test$optimal_index])
suppressMessages(library(locfit))
DEgenes <- find_markers(expression_matrix=assay(mixedpop2),
                            cluster = colData(mixedpop2)[,1],
                            selected_cluster=unique(colData(mixedpop2)[,1]))
#the output contains dataframes for each cluster.
#the data frame contains all genes, sorted by p-values
names(DEgenes)
#> [1] "baseMean"       "log2FoldChange" "lfcSE"          "stat"          
#> [5] "pvalue"         "padj"           "id"
#you can annotate the identified clusters
DEgeneList_1vsOthers <- DEgenes$DE_Subpop1vsRemaining$id
#users need to check the format of the gene input to make sure they are
#consistent to the gene names in the expression matrix
#the following command saves the file "PathwayEnrichment.xlsx" to the
#working dir
#use 500 top DE genes
suppressMessages(library(DOSE))
suppressMessages(library(ReactomePA))
suppressMessages(library(clusterProfiler))
genes500 <- as.factor(DEgeneList_1vsOthers[seq_len(500)])
enrichment_test <- annotate_clusters(genes, pvalueCutoff=0.05, gene_symbol=TRUE)
#the enrichment outputs can be displayed by running
clusterProfiler::dotplot(enrichment_test, showCategory=10, font.size = 6)The purpose of this workflow is to solve the following task:
#select a subpopulation, and input gene list
c_selectID <- 1
#note make sure the format for genes input here is the same to the format
#for genes in the mixedpop1 and mixedpop2
genes = DEgenes$id[1:500]
#run the test bootstrap with nboots = 2 runs
cluster_mixedpop1 <- colData(mixedpop1)[,1]
cluster_mixedpop2 <- colData(mixedpop2)[,1]
LSOLDA_dat <- bootstrap_prediction(nboots = 2, mixedpop1 = mixedpop1,
                             mixedpop2 = mixedpop2, genes = genes, 
                             c_selectID  = c_selectID,
                             listData = list(),
                             cluster_mixedpop1 = cluster_mixedpop1,
                             cluster_mixedpop2 = cluster_mixedpop2)#get the number of rows for the summary matrix
row_cluster <-length(unique(colData(mixedpop2)[,1]))
#summary results LDA to to show the percent of cells classified as cells
#belonging by LDA classifier
summary_prediction_lda(LSOLDA_dat=LSOLDA_dat, nPredSubpop = row_cluster )
#>                 V1               V2                                names
#> 1 6.95187165775401 40.1069518716578 LDA for subpop 1 in target mixedpop2
#> 2 53.5714285714286 98.5714285714286 LDA for subpop 2 in target mixedpop2
#> 3 1.50375939849624 42.1052631578947 LDA for subpop 3 in target mixedpop2
#> 4               15               60 LDA for subpop 4 in target mixedpop2
#summary results Lasso to show the percent of cells classified as cells
#belonging by Lasso classifier
summary_prediction_lasso(LSOLDA_dat=LSOLDA_dat, nPredSubpop = row_cluster)
#>                 V1               V2                                      names
#> 1 36.3636363636364 64.7058823529412 ElasticNet for subpop1 in target mixedpop2
#> 2              100              100 ElasticNet for subpop2 in target mixedpop2
#> 3 83.4586466165414 86.4661654135338 ElasticNet for subpop3 in target mixedpop2
#> 4               90               90 ElasticNet for subpop4 in target mixedpop2
# summary maximum deviance explained by the model during the model training
summary_deviance(object = LSOLDA_dat)
#> $allDeviance
#> [1] "50.91" "52.15"
#> 
#> $DeviMax
#>           dat_DE$Dfd          Deviance           DEgenes
#> 1                  0             52.15    genes_cluster1
#> 2                  1             52.15    genes_cluster1
#> 3                  3             52.15    genes_cluster1
#> 4                  5             52.15    genes_cluster1
#> 5                  6             52.15    genes_cluster1
#> 6                  9             52.15    genes_cluster1
#> 7                 10             52.15    genes_cluster1
#> 8                 13             52.15    genes_cluster1
#> 9                 14             52.15    genes_cluster1
#> 10                20             52.15    genes_cluster1
#> 11                22             52.15    genes_cluster1
#> 12                26             52.15    genes_cluster1
#> 13                28             52.15    genes_cluster1
#> 14                32             52.15    genes_cluster1
#> 15                36             52.15    genes_cluster1
#> 16                41             52.15    genes_cluster1
#> 17                45             52.15    genes_cluster1
#> 18 remaining DEgenes remaining DEgenes remaining DEgenes
#> 
#> $LassoGenesMax
#> NULL
# summary accuracy to check the model accuracy in the leave-out test set
summary_accuracy(object = LSOLDA_dat)
#> [1] 68.75000 67.41071Here we look at one example use case to find relationship between clusters within one sample or between two sample
#run prediction for 3 clusters
cluster_mixedpop1 <- colData(mixedpop1)[,1]
cluster_mixedpop2 <- colData(mixedpop2)[,1]
#cluster_mixedpop2 <- as.numeric(as.vector(colData(mixedpop2)[,1]))
c_selectID <- 1
#top 200 gene markers distinguishing cluster 1
genes = DEgenes$id[1:200]
LSOLDA_dat1 <- bootstrap_prediction(nboots = 2, mixedpop1 = mixedpop2,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop2,
                        cluster_mixedpop2 = cluster_mixedpop2)
c_selectID <- 2
genes = DEgenes$id[1:200]
LSOLDA_dat2 <- bootstrap_prediction(nboots = 2,mixedpop1 = mixedpop2,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop2,
                        cluster_mixedpop2 = cluster_mixedpop2)
c_selectID <- 3
genes = DEgenes$id[1:200]
LSOLDA_dat3 <- bootstrap_prediction(nboots = 2,mixedpop1 = mixedpop2,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop2,
                        cluster_mixedpop2 = cluster_mixedpop2)
c_selectID <- 4
genes = DEgenes$id[1:200]
LSOLDA_dat4 <- bootstrap_prediction(nboots = 2,mixedpop1 = mixedpop2,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop2,
                        cluster_mixedpop2 = cluster_mixedpop2)
#prepare table input for sankey plot
LASSO_C1S2  <- reformat_LASSO(c_selectID=1, mp_selectID = 2,
                             LSOLDA_dat=LSOLDA_dat1,
                             nPredSubpop = length(unique(colData(mixedpop2)
                                [,1])),
                             Nodes_group ="#7570b3")
LASSO_C2S2  <- reformat_LASSO(c_selectID=2, mp_selectID =2,
                             LSOLDA_dat=LSOLDA_dat2,
                             nPredSubpop = length(unique(colData(mixedpop2)
                                [,1])),
                             Nodes_group ="#1b9e77")
LASSO_C3S2  <- reformat_LASSO(c_selectID=3, mp_selectID =2,
                             LSOLDA_dat=LSOLDA_dat3,
                             nPredSubpop = length(unique(colData(mixedpop2)
                                [,1])),
                             Nodes_group ="#e7298a")
LASSO_C4S2  <- reformat_LASSO(c_selectID=4, mp_selectID =2,
                             LSOLDA_dat=LSOLDA_dat4,
                             nPredSubpop = length(unique(colData(mixedpop2)
                                [,1])),
                             Nodes_group ="#00FFFF")
combined <- rbind(LASSO_C1S2,LASSO_C2S2,LASSO_C3S2, LASSO_C4S2 )
combined <- combined[is.na(combined$Value) != TRUE,]
nboots = 2
#links: source, target, value
#source: node, nodegroup
combined_D3obj <-list(Nodes=combined[,(nboots+3):(nboots+4)],
                     Links=combined[,c((nboots+2):(nboots+1),ncol(combined))])
library(networkD3)
Node_source <- as.vector(sort(unique(combined_D3obj$Links$Source)))
Node_target <- as.vector(sort(unique(combined_D3obj$Links$Target)))
Node_all <-unique(c(Node_source, Node_target))
#assign IDs for Source (start from 0)
Source <-combined_D3obj$Links$Source
Target <- combined_D3obj$Links$Target
for(i in 1:length(Node_all)){
   Source[Source==Node_all[i]] <-i-1
   Target[Target==Node_all[i]] <-i-1
}
# 
combined_D3obj$Links$Source <- as.numeric(Source)
combined_D3obj$Links$Target <- as.numeric(Target)
combined_D3obj$Links$LinkColor <- combined$NodeGroup
#prepare node info
node_df <-data.frame(Node=Node_all)
node_df$id <-as.numeric(c(0, 1:(length(Node_all)-1)))
suppressMessages(library(dplyr))
Color <- combined %>% count(Node, color=NodeGroup) %>% select(2)
node_df$color <- Color$color
suppressMessages(library(networkD3))
p1<-sankeyNetwork(Links =combined_D3obj$Links, Nodes = node_df,
                 Value = "Value", NodeGroup ="color", LinkGroup = "LinkColor", 
                 NodeID="Node", Source="Source", Target="Target", fontSize = 22)
p1Here we look at one example use case to find relationship between clusters within one sample or between two sample
#run prediction for 3 clusters
cluster_mixedpop1 <- colData(mixedpop1)[,1]
cluster_mixedpop2 <- colData(mixedpop2)[,1]
row_cluster <-length(unique(colData(mixedpop2)[,1]))
c_selectID <- 1
#top 200 gene markers distinguishing cluster 1
genes = DEgenes$id[1:200]
LSOLDA_dat1 <- bootstrap_prediction(nboots = 2, mixedpop1 = mixedpop1,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop1,
                        cluster_mixedpop2 = cluster_mixedpop2)
c_selectID <- 2
genes = DEgenes$id[1:200]
LSOLDA_dat2 <- bootstrap_prediction(nboots = 2,mixedpop1 = mixedpop1,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop1,
                        cluster_mixedpop2 = cluster_mixedpop2)
c_selectID <- 3
genes = DEgenes$id[1:200]
LSOLDA_dat3 <- bootstrap_prediction(nboots = 2,mixedpop1 = mixedpop1,
                        mixedpop2 = mixedpop2, genes=genes, c_selectID, 
                        listData =list(),
                        cluster_mixedpop1 = cluster_mixedpop1,
                        cluster_mixedpop2 = cluster_mixedpop2)
#prepare table input for sankey plot
LASSO_C1S1  <- reformat_LASSO(c_selectID=1, mp_selectID = 1,
                             LSOLDA_dat=LSOLDA_dat1, nPredSubpop = row_cluster, 
                             Nodes_group = "#7570b3")
LASSO_C2S1  <- reformat_LASSO(c_selectID=2, mp_selectID = 1,
                             LSOLDA_dat=LSOLDA_dat2, nPredSubpop = row_cluster, 
                             Nodes_group = "#1b9e77")
LASSO_C3S1  <- reformat_LASSO(c_selectID=3, mp_selectID = 1,
                             LSOLDA_dat=LSOLDA_dat3, nPredSubpop = row_cluster, 
                             Nodes_group = "#e7298a")
combined <- rbind(LASSO_C1S1,LASSO_C2S1,LASSO_C3S1)
nboots = 2
#links: source, target, value
#source: node, nodegroup
combined_D3obj <-list(Nodes=combined[,(nboots+3):(nboots+4)],
                     Links=combined[,c((nboots+2):(nboots+1),ncol(combined))])
combined <- combined[is.na(combined$Value) != TRUE,]
library(networkD3)
Node_source <- as.vector(sort(unique(combined_D3obj$Links$Source)))
Node_target <- as.vector(sort(unique(combined_D3obj$Links$Target)))
Node_all <-unique(c(Node_source, Node_target))
#assign IDs for Source (start from 0)
Source <-combined_D3obj$Links$Source
Target <- combined_D3obj$Links$Target
for(i in 1:length(Node_all)){
   Source[Source==Node_all[i]] <-i-1
   Target[Target==Node_all[i]] <-i-1
}
combined_D3obj$Links$Source <- as.numeric(Source)
combined_D3obj$Links$Target <- as.numeric(Target)
combined_D3obj$Links$LinkColor <- combined$NodeGroup
#prepare node info
node_df <-data.frame(Node=Node_all)
node_df$id <-as.numeric(c(0, 1:(length(Node_all)-1)))
suppressMessages(library(dplyr))
n <- length(unique(node_df$Node))
getPalette = colorRampPalette(RColorBrewer::brewer.pal(9, "Set1"))
Color = getPalette(n)
node_df$color <- Color
suppressMessages(library(networkD3))
p1<-sankeyNetwork(Links =combined_D3obj$Links, Nodes = node_df,
                 Value = "Value", NodeGroup ="color", LinkGroup = "LinkColor",
                 NodeID="Node", Source="Source", Target="Target", fontSize = 22)
p1devtools::session_info()
#> ─ Session info ───────────────────────────────────────────────────────────────
#>  setting  value
#>  version  R version 4.3.1 (2023-06-16)
#>  os       Ubuntu 22.04.3 LTS
#>  system   x86_64, linux-gnu
#>  ui       X11
#>  language (EN)
#>  collate  C
#>  ctype    en_US.UTF-8
#>  tz       America/New_York
#>  date     2023-10-24
#>  pandoc   2.7.3 @ /usr/bin/ (via rmarkdown)
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────
#>  package                * version    date (UTC) lib source
#>  abind                    1.4-5      2016-07-21 [2] CRAN (R 4.3.1)
#>  AnnotationDbi          * 1.64.0     2023-10-24 [2] Bioconductor
#>  AnnotationHub            3.10.0     2023-10-24 [2] Bioconductor
#>  ape                      5.7-1      2023-03-13 [2] CRAN (R 4.3.1)
#>  aplot                    0.2.2      2023-10-06 [2] CRAN (R 4.3.1)
#>  Biobase                * 2.62.0     2023-10-24 [2] Bioconductor
#>  BiocFileCache            2.10.0     2023-10-24 [2] Bioconductor
#>  BiocGenerics           * 0.48.0     2023-10-24 [2] Bioconductor
#>  BiocManager              1.30.22    2023-08-08 [2] CRAN (R 4.3.1)
#>  BiocParallel             1.36.0     2023-10-24 [2] Bioconductor
#>  BiocVersion              3.18.0     2023-10-24 [2] Bioconductor
#>  Biostrings               2.70.0     2023-10-24 [2] Bioconductor
#>  bit                      4.0.5      2022-11-15 [2] CRAN (R 4.3.1)
#>  bit64                    4.0.5      2020-08-30 [2] CRAN (R 4.3.1)
#>  bitops                   1.0-7      2021-04-24 [2] CRAN (R 4.3.1)
#>  blob                     1.2.4      2023-03-17 [2] CRAN (R 4.3.1)
#>  bslib                    0.5.1      2023-08-11 [2] CRAN (R 4.3.1)
#>  cachem                   1.0.8      2023-05-01 [2] CRAN (R 4.3.1)
#>  callr                    3.7.3      2022-11-02 [2] CRAN (R 4.3.1)
#>  caret                  * 6.0-94     2023-03-21 [2] CRAN (R 4.3.1)
#>  class                    7.3-22     2023-05-03 [3] CRAN (R 4.3.1)
#>  cli                      3.6.1      2023-03-23 [2] CRAN (R 4.3.1)
#>  clusterProfiler        * 4.10.0     2023-10-24 [2] Bioconductor
#>  codetools                0.2-19     2023-02-01 [3] CRAN (R 4.3.1)
#>  colorspace               2.1-0      2023-01-23 [2] CRAN (R 4.3.1)
#>  cowplot                  1.1.1      2020-12-30 [2] CRAN (R 4.3.1)
#>  crayon                   1.5.2      2022-09-29 [2] CRAN (R 4.3.1)
#>  curl                     5.1.0      2023-10-02 [2] CRAN (R 4.3.1)
#>  data.table               1.14.8     2023-02-17 [2] CRAN (R 4.3.1)
#>  DBI                      1.1.3      2022-06-18 [2] CRAN (R 4.3.1)
#>  dbplyr                   2.3.4      2023-09-26 [2] CRAN (R 4.3.1)
#>  DelayedArray             0.28.0     2023-10-24 [2] Bioconductor
#>  dendextend               1.17.1     2023-03-25 [2] CRAN (R 4.3.1)
#>  DESeq2                 * 1.42.0     2023-10-24 [2] Bioconductor
#>  devtools                 2.4.5      2022-10-11 [2] CRAN (R 4.3.1)
#>  digest                   0.6.33     2023-07-07 [2] CRAN (R 4.3.1)
#>  DOSE                   * 3.28.0     2023-10-24 [2] Bioconductor
#>  dplyr                  * 1.1.3      2023-09-03 [2] CRAN (R 4.3.1)
#>  dynamicTreeCut         * 1.63-1     2016-03-11 [2] CRAN (R 4.3.1)
#>  e1071                    1.7-13     2023-02-01 [2] CRAN (R 4.3.1)
#>  ellipsis                 0.3.2      2021-04-29 [2] CRAN (R 4.3.1)
#>  enrichplot               1.22.0     2023-10-24 [2] Bioconductor
#>  evaluate                 0.22       2023-09-29 [2] CRAN (R 4.3.1)
#>  fansi                    1.0.5      2023-10-08 [2] CRAN (R 4.3.1)
#>  farver                   2.1.1      2022-07-06 [2] CRAN (R 4.3.1)
#>  fastcluster              1.2.3      2021-05-24 [2] CRAN (R 4.3.1)
#>  fastmap                  1.1.1      2023-02-24 [2] CRAN (R 4.3.1)
#>  fastmatch                1.1-4      2023-08-18 [2] CRAN (R 4.3.1)
#>  fgsea                    1.28.0     2023-10-24 [2] Bioconductor
#>  filelock                 1.0.2      2018-10-05 [2] CRAN (R 4.3.1)
#>  foreach                  1.5.2      2022-02-02 [2] CRAN (R 4.3.1)
#>  fs                       1.6.3      2023-07-20 [2] CRAN (R 4.3.1)
#>  future                   1.33.0     2023-07-01 [2] CRAN (R 4.3.1)
#>  future.apply             1.11.0     2023-05-21 [2] CRAN (R 4.3.1)
#>  generics                 0.1.3      2022-07-05 [2] CRAN (R 4.3.1)
#>  GenomeInfoDb           * 1.38.0     2023-10-24 [2] Bioconductor
#>  GenomeInfoDbData         1.2.11     2023-10-17 [2] Bioconductor
#>  GenomicRanges          * 1.54.0     2023-10-24 [2] Bioconductor
#>  ggforce                  0.4.1      2022-10-04 [2] CRAN (R 4.3.1)
#>  ggfun                    0.1.3      2023-09-15 [2] CRAN (R 4.3.1)
#>  ggplot2                * 3.4.4      2023-10-12 [2] CRAN (R 4.3.1)
#>  ggplotify                0.1.2      2023-08-09 [2] CRAN (R 4.3.1)
#>  ggraph                   2.1.0      2022-10-09 [2] CRAN (R 4.3.1)
#>  ggrepel                  0.9.4      2023-10-13 [2] CRAN (R 4.3.1)
#>  ggtree                   3.10.0     2023-10-24 [2] Bioconductor
#>  glmnet                   4.1-8      2023-08-22 [2] CRAN (R 4.3.1)
#>  globals                  0.16.2     2022-11-21 [2] CRAN (R 4.3.1)
#>  glue                     1.6.2      2022-02-24 [2] CRAN (R 4.3.1)
#>  GO.db                    3.18.0     2023-09-21 [2] Bioconductor
#>  GOSemSim                 2.28.0     2023-10-24 [2] Bioconductor
#>  gower                    1.0.1      2022-12-22 [2] CRAN (R 4.3.1)
#>  graph                    1.80.0     2023-10-24 [2] Bioconductor
#>  graphite                 1.48.0     2023-10-24 [2] Bioconductor
#>  graphlayouts             1.0.1      2023-09-19 [2] CRAN (R 4.3.1)
#>  gridExtra                2.3        2017-09-09 [2] CRAN (R 4.3.1)
#>  gridGraphics             0.5-1      2020-12-13 [2] CRAN (R 4.3.1)
#>  gson                     0.1.0      2023-03-07 [2] CRAN (R 4.3.1)
#>  gtable                   0.3.4      2023-08-21 [2] CRAN (R 4.3.1)
#>  hardhat                  1.3.0      2023-03-30 [2] CRAN (R 4.3.1)
#>  HDO.db                   0.99.1     2023-06-20 [2] Bioconductor
#>  HPO.db                   0.99.2     2023-10-18 [2] Bioconductor
#>  htmltools                0.5.6.1    2023-10-06 [2] CRAN (R 4.3.1)
#>  htmlwidgets              1.6.2      2023-03-17 [2] CRAN (R 4.3.1)
#>  httpuv                   1.6.12     2023-10-23 [2] CRAN (R 4.3.1)
#>  httr                     1.4.7      2023-08-15 [2] CRAN (R 4.3.1)
#>  igraph                   1.5.1      2023-08-10 [2] CRAN (R 4.3.1)
#>  interactiveDisplayBase   1.40.0     2023-10-24 [2] Bioconductor
#>  ipred                    0.9-14     2023-03-09 [2] CRAN (R 4.3.1)
#>  IRanges                * 2.36.0     2023-10-24 [2] Bioconductor
#>  iterators                1.0.14     2022-02-05 [2] CRAN (R 4.3.1)
#>  jquerylib                0.1.4      2021-04-26 [2] CRAN (R 4.3.1)
#>  jsonlite                 1.8.7      2023-06-29 [2] CRAN (R 4.3.1)
#>  KEGGREST                 1.42.0     2023-10-24 [2] Bioconductor
#>  knitr                    1.44       2023-09-11 [2] CRAN (R 4.3.1)
#>  labeling                 0.4.3      2023-08-29 [2] CRAN (R 4.3.1)
#>  later                    1.3.1      2023-05-02 [2] CRAN (R 4.3.1)
#>  lattice                * 0.22-5     2023-10-24 [3] CRAN (R 4.3.1)
#>  lava                     1.7.2.1    2023-02-27 [2] CRAN (R 4.3.1)
#>  lazyeval                 0.2.2      2019-03-15 [2] CRAN (R 4.3.1)
#>  lifecycle                1.0.3      2022-10-07 [2] CRAN (R 4.3.1)
#>  listenv                  0.9.0      2022-12-16 [2] CRAN (R 4.3.1)
#>  locfit                 * 1.5-9.8    2023-06-11 [2] CRAN (R 4.3.1)
#>  lubridate                1.9.3      2023-09-27 [2] CRAN (R 4.3.1)
#>  magrittr                 2.0.3      2022-03-30 [2] CRAN (R 4.3.1)
#>  MASS                     7.3-60     2023-05-04 [3] CRAN (R 4.3.1)
#>  Matrix                   1.6-1.1    2023-09-18 [3] CRAN (R 4.3.1)
#>  MatrixGenerics         * 1.14.0     2023-10-24 [2] Bioconductor
#>  matrixStats            * 1.0.0      2023-06-02 [2] CRAN (R 4.3.1)
#>  memoise                  2.0.1      2021-11-26 [2] CRAN (R 4.3.1)
#>  mime                     0.12       2021-09-28 [2] CRAN (R 4.3.1)
#>  miniUI                   0.1.1.1    2018-05-18 [2] CRAN (R 4.3.1)
#>  ModelMetrics             1.2.2.2    2020-03-17 [2] CRAN (R 4.3.1)
#>  MPO.db                   0.99.7     2023-10-18 [2] Bioconductor
#>  munsell                  0.5.0      2018-06-12 [2] CRAN (R 4.3.1)
#>  networkD3              * 0.4        2017-03-18 [2] CRAN (R 4.3.1)
#>  nlme                     3.1-163    2023-08-09 [3] CRAN (R 4.3.1)
#>  nnet                     7.3-19     2023-05-03 [3] CRAN (R 4.3.1)
#>  org.Hs.eg.db           * 3.18.0     2023-09-21 [2] Bioconductor
#>  parallelly               1.36.0     2023-05-26 [2] CRAN (R 4.3.1)
#>  patchwork                1.1.3      2023-08-14 [2] CRAN (R 4.3.1)
#>  pillar                   1.9.0      2023-03-22 [2] CRAN (R 4.3.1)
#>  pkgbuild                 1.4.2      2023-06-26 [2] CRAN (R 4.3.1)
#>  pkgconfig                2.0.3      2019-09-22 [2] CRAN (R 4.3.1)
#>  pkgload                  1.3.3      2023-09-22 [2] CRAN (R 4.3.1)
#>  plyr                     1.8.9      2023-10-02 [2] CRAN (R 4.3.1)
#>  png                      0.1-8      2022-11-29 [2] CRAN (R 4.3.1)
#>  polyclip                 1.10-6     2023-09-27 [2] CRAN (R 4.3.1)
#>  prettyunits              1.2.0      2023-09-24 [2] CRAN (R 4.3.1)
#>  pROC                     1.18.4     2023-07-06 [2] CRAN (R 4.3.1)
#>  processx                 3.8.2      2023-06-30 [2] CRAN (R 4.3.1)
#>  prodlim                  2023.08.28 2023-08-28 [2] CRAN (R 4.3.1)
#>  profvis                  0.3.8      2023-05-02 [2] CRAN (R 4.3.1)
#>  promises                 1.2.1      2023-08-10 [2] CRAN (R 4.3.1)
#>  proxy                    0.4-27     2022-06-09 [2] CRAN (R 4.3.1)
#>  ps                       1.7.5      2023-04-18 [2] CRAN (R 4.3.1)
#>  purrr                    1.0.2      2023-08-10 [2] CRAN (R 4.3.1)
#>  qvalue                   2.34.0     2023-10-24 [2] Bioconductor
#>  R6                       2.5.1      2021-08-19 [2] CRAN (R 4.3.1)
#>  rappdirs                 0.3.3      2021-01-31 [2] CRAN (R 4.3.1)
#>  RColorBrewer             1.1-3      2022-04-03 [2] CRAN (R 4.3.1)
#>  Rcpp                     1.0.11     2023-07-06 [2] CRAN (R 4.3.1)
#>  RcppArmadillo            0.12.6.4.0 2023-09-10 [2] CRAN (R 4.3.1)
#>  RcppParallel             5.1.7      2023-02-27 [2] CRAN (R 4.3.1)
#>  RCurl                    1.98-1.12  2023-03-27 [2] CRAN (R 4.3.1)
#>  reactome.db              1.86.0     2023-10-17 [2] Bioconductor
#>  ReactomePA             * 1.46.0     2023-10-24 [2] Bioconductor
#>  recipes                  1.0.8      2023-08-25 [2] CRAN (R 4.3.1)
#>  remotes                  2.4.2.1    2023-07-18 [2] CRAN (R 4.3.1)
#>  reshape2                 1.4.4      2020-04-09 [2] CRAN (R 4.3.1)
#>  rlang                    1.1.1      2023-04-28 [2] CRAN (R 4.3.1)
#>  rmarkdown                2.25       2023-09-18 [2] CRAN (R 4.3.1)
#>  rpart                    4.1.21     2023-10-09 [3] CRAN (R 4.3.1)
#>  RSQLite                  2.3.1      2023-04-03 [2] CRAN (R 4.3.1)
#>  Rtsne                    0.16       2022-04-17 [2] CRAN (R 4.3.1)
#>  S4Arrays                 1.2.0      2023-10-24 [2] Bioconductor
#>  S4Vectors              * 0.40.0     2023-10-24 [2] Bioconductor
#>  sass                     0.4.7      2023-07-15 [2] CRAN (R 4.3.1)
#>  scales                   1.2.1      2022-08-20 [2] CRAN (R 4.3.1)
#>  scatterpie               0.2.1      2023-06-07 [2] CRAN (R 4.3.1)
#>  scGPS                  * 1.16.0     2023-10-24 [1] Bioconductor
#>  sessioninfo              1.2.2      2021-12-06 [2] CRAN (R 4.3.1)
#>  shadowtext               0.1.2      2022-04-22 [2] CRAN (R 4.3.1)
#>  shape                    1.4.6      2021-05-19 [2] CRAN (R 4.3.1)
#>  shiny                    1.7.5.1    2023-10-14 [2] CRAN (R 4.3.1)
#>  SingleCellExperiment   * 1.24.0     2023-10-24 [2] Bioconductor
#>  SparseArray              1.2.0      2023-10-24 [2] Bioconductor
#>  stringi                  1.7.12     2023-01-11 [2] CRAN (R 4.3.1)
#>  stringr                  1.5.0      2022-12-02 [2] CRAN (R 4.3.1)
#>  SummarizedExperiment   * 1.32.0     2023-10-24 [2] Bioconductor
#>  survival                 3.5-7      2023-08-14 [3] CRAN (R 4.3.1)
#>  tibble                   3.2.1      2023-03-20 [2] CRAN (R 4.3.1)
#>  tidygraph                1.2.3      2023-02-01 [2] CRAN (R 4.3.1)
#>  tidyr                    1.3.0      2023-01-24 [2] CRAN (R 4.3.1)
#>  tidyselect               1.2.0      2022-10-10 [2] CRAN (R 4.3.1)
#>  tidytree                 0.4.5      2023-08-10 [2] CRAN (R 4.3.1)
#>  timechange               0.2.0      2023-01-11 [2] CRAN (R 4.3.1)
#>  timeDate                 4022.108   2023-01-07 [2] CRAN (R 4.3.1)
#>  treeio                   1.26.0     2023-10-24 [2] Bioconductor
#>  tweenr                   2.0.2      2022-09-06 [2] CRAN (R 4.3.1)
#>  urlchecker               1.0.1      2021-11-30 [2] CRAN (R 4.3.1)
#>  usethis                  2.2.2      2023-07-06 [2] CRAN (R 4.3.1)
#>  utf8                     1.2.4      2023-10-22 [2] CRAN (R 4.3.1)
#>  vctrs                    0.6.4      2023-10-12 [2] CRAN (R 4.3.1)
#>  viridis                  0.6.4      2023-07-22 [2] CRAN (R 4.3.1)
#>  viridisLite              0.4.2      2023-05-02 [2] CRAN (R 4.3.1)
#>  withr                    2.5.1      2023-09-26 [2] CRAN (R 4.3.1)
#>  xfun                     0.40       2023-08-09 [2] CRAN (R 4.3.1)
#>  xtable                   1.8-4      2019-04-21 [2] CRAN (R 4.3.1)
#>  XVector                  0.42.0     2023-10-24 [2] Bioconductor
#>  yaml                     2.3.7      2023-01-23 [2] CRAN (R 4.3.1)
#>  yulab.utils              0.1.0      2023-09-20 [2] CRAN (R 4.3.1)
#>  zlibbioc                 1.48.0     2023-10-24 [2] Bioconductor
#> 
#>  [1] /tmp/RtmpnJviOk/Rinst339749cf88fe9
#>  [2] /home/biocbuild/bbs-3.18-bioc/R/site-library
#>  [3] /home/biocbuild/bbs-3.18-bioc/R/library
#> 
#> ──────────────────────────────────────────────────────────────────────────────