# 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.37209 62.50000 65.25822
# summary maximum deviance explained by the model 
summary_deviance(object = LSOLDA_dat)
#> $allDeviance
#> [1] "8.21" "7.9"  "6.56"
#> 
#> $DeviMax
#>          dat_DE$Dfd          Deviance           DEgenes
#> 1                 0              8.21    genes_cluster1
#> 2                 1              8.21    genes_cluster1
#> 3                 2              8.21    genes_cluster1
#> 4                 3              8.21    genes_cluster1
#> 5                 4              8.21    genes_cluster1
#> 6 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
#> 24112NANANANANANA
#> 3
#> 24112224NANANANANA
#> 4
#> 24112224299NANANANA
#> 5
#> 24112224299335NANANA
#> 6
#> 24112224299335367NANA
#> 7
#> 24112224299335367414NA
#> 8
#> 24112224299335367414470
#> 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>.
#> 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 2.67379679144385 68.4491978609626 LDA for subpop 1 in target mixedpop2
#> 2               20 58.5714285714286 LDA for subpop 2 in target mixedpop2
#> 3 3.00751879699248 64.6616541353383 LDA for subpop 3 in target mixedpop2
#> 4              7.5               75 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 55.0802139037433  46.524064171123 ElasticNet for subpop1 in target mixedpop2
#> 2 79.2857142857143 96.4285714285714 ElasticNet for subpop2 in target mixedpop2
#> 3 41.3533834586466 34.5864661654135 ElasticNet for subpop3 in target mixedpop2
#> 4               55               65 ElasticNet for subpop4 in target mixedpop2
# summary maximum deviance explained by the model during the model training
summary_deviance(object = LSOLDA_dat)
#> $allDeviance
#> [1] "98.23" "61.41"
#> 
#> $DeviMax
#>           dat_DE$Dfd          Deviance           DEgenes
#> 1                  0             98.23    genes_cluster1
#> 2                  1             98.23    genes_cluster1
#> 3                  2             98.23    genes_cluster1
#> 4                  5             98.23    genes_cluster1
#> 5                  7             98.23    genes_cluster1
#> 6                 11             98.23    genes_cluster1
#> 7                 12             98.23    genes_cluster1
#> 8                 13             98.23    genes_cluster1
#> 9                 16             98.23    genes_cluster1
#> 10                17             98.23    genes_cluster1
#> 11                19             98.23    genes_cluster1
#> 12                20             98.23    genes_cluster1
#> 13                21             98.23    genes_cluster1
#> 14                24             98.23    genes_cluster1
#> 15                27             98.23    genes_cluster1
#> 16                30             98.23    genes_cluster1
#> 17                33             98.23    genes_cluster1
#> 18                34             98.23    genes_cluster1
#> 19                36             98.23    genes_cluster1
#> 20                42             98.23    genes_cluster1
#> 21                44             98.23    genes_cluster1
#> 22                46             98.23    genes_cluster1
#> 23                48             98.23    genes_cluster1
#> 24                50             98.23    genes_cluster1
#> 25                54             98.23    genes_cluster1
#> 26                56             98.23    genes_cluster1
#> 27                57             98.23    genes_cluster1
#> 28                58             98.23    genes_cluster1
#> 29                61             98.23    genes_cluster1
#> 30                63             98.23    genes_cluster1
#> 31                64             98.23    genes_cluster1
#> 32                66             98.23    genes_cluster1
#> 33                67             98.23    genes_cluster1
#> 34                68             98.23    genes_cluster1
#> 35                69             98.23    genes_cluster1
#> 36                71             98.23    genes_cluster1
#> 37                72             98.23    genes_cluster1
#> 38                73             98.23    genes_cluster1
#> 39                76             98.23    genes_cluster1
#> 40                77             98.23    genes_cluster1
#> 41                79             98.23    genes_cluster1
#> 42                80             98.23    genes_cluster1
#> 43                82             98.23    genes_cluster1
#> 44                83             98.23    genes_cluster1
#> 45                84             98.23    genes_cluster1
#> 46                85             98.23    genes_cluster1
#> 47 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] 65.62500 66.07143Here 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.2.2 (2022-10-31)
#>  os       Ubuntu 20.04.5 LTS
#>  system   x86_64, linux-gnu
#>  ui       X11
#>  language (EN)
#>  collate  C
#>  ctype    en_US.UTF-8
#>  tz       America/New_York
#>  date     2023-01-19
#>  pandoc   2.5 @ /usr/bin/ (via rmarkdown)
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────
#>  package              * version    date (UTC) lib source
#>  annotate               1.76.0     2023-01-19 [2] Bioconductor
#>  AnnotationDbi        * 1.60.0     2023-01-19 [2] Bioconductor
#>  ape                    5.6-2      2022-03-02 [2] CRAN (R 4.2.2)
#>  aplot                  0.1.9      2022-11-24 [2] CRAN (R 4.2.2)
#>  assertthat             0.2.1      2019-03-21 [2] CRAN (R 4.2.2)
#>  Biobase              * 2.58.0     2023-01-19 [2] Bioconductor
#>  BiocGenerics         * 0.44.0     2023-01-19 [2] Bioconductor
#>  BiocParallel           1.32.5     2023-01-19 [2] Bioconductor
#>  Biostrings             2.66.0     2023-01-19 [2] Bioconductor
#>  bit                    4.0.5      2022-11-15 [2] CRAN (R 4.2.2)
#>  bit64                  4.0.5      2020-08-30 [2] CRAN (R 4.2.2)
#>  bitops                 1.0-7      2021-04-24 [2] CRAN (R 4.2.2)
#>  blob                   1.2.3      2022-04-10 [2] CRAN (R 4.2.2)
#>  bslib                  0.4.2      2022-12-16 [2] CRAN (R 4.2.2)
#>  cachem                 1.0.6      2021-08-19 [2] CRAN (R 4.2.2)
#>  callr                  3.7.3      2022-11-02 [2] CRAN (R 4.2.2)
#>  caret                * 6.0-93     2022-08-09 [2] CRAN (R 4.2.2)
#>  class                  7.3-20     2022-01-16 [2] CRAN (R 4.2.2)
#>  cli                    3.6.0      2023-01-09 [2] CRAN (R 4.2.2)
#>  clusterProfiler      * 4.6.0      2023-01-19 [2] Bioconductor
#>  codetools              0.2-18     2020-11-04 [2] CRAN (R 4.2.2)
#>  colorspace             2.0-3      2022-02-21 [2] CRAN (R 4.2.2)
#>  cowplot                1.1.1      2020-12-30 [2] CRAN (R 4.2.2)
#>  crayon                 1.5.2      2022-09-29 [2] CRAN (R 4.2.2)
#>  data.table             1.14.6     2022-11-16 [2] CRAN (R 4.2.2)
#>  DBI                    1.1.3      2022-06-18 [2] CRAN (R 4.2.2)
#>  DelayedArray           0.24.0     2023-01-19 [2] Bioconductor
#>  dendextend             1.16.0     2022-07-04 [2] CRAN (R 4.2.2)
#>  DESeq2                 1.38.3     2023-01-19 [2] Bioconductor
#>  devtools               2.4.5      2022-10-11 [2] CRAN (R 4.2.2)
#>  digest                 0.6.31     2022-12-11 [2] CRAN (R 4.2.2)
#>  DOSE                 * 3.24.2     2023-01-19 [2] Bioconductor
#>  downloader             0.4        2015-07-09 [2] CRAN (R 4.2.2)
#>  dplyr                * 1.0.10     2022-09-01 [2] CRAN (R 4.2.2)
#>  dynamicTreeCut       * 1.63-1     2016-03-11 [2] CRAN (R 4.2.2)
#>  e1071                  1.7-12     2022-10-24 [2] CRAN (R 4.2.2)
#>  ellipsis               0.3.2      2021-04-29 [2] CRAN (R 4.2.2)
#>  enrichplot             1.18.3     2023-01-19 [2] Bioconductor
#>  evaluate               0.20       2023-01-17 [2] CRAN (R 4.2.2)
#>  fansi                  1.0.3      2022-03-24 [2] CRAN (R 4.2.2)
#>  farver                 2.1.1      2022-07-06 [2] CRAN (R 4.2.2)
#>  fastcluster            1.2.3      2021-05-24 [2] CRAN (R 4.2.2)
#>  fastmap                1.1.0      2021-01-25 [2] CRAN (R 4.2.2)
#>  fastmatch              1.1-3      2021-07-23 [2] CRAN (R 4.2.2)
#>  fgsea                  1.24.0     2023-01-19 [2] Bioconductor
#>  foreach                1.5.2      2022-02-02 [2] CRAN (R 4.2.2)
#>  fs                     1.5.2      2021-12-08 [2] CRAN (R 4.2.2)
#>  future                 1.30.0     2022-12-16 [2] CRAN (R 4.2.2)
#>  future.apply           1.10.0     2022-11-05 [2] CRAN (R 4.2.2)
#>  geneplotter            1.76.0     2023-01-19 [2] Bioconductor
#>  generics               0.1.3      2022-07-05 [2] CRAN (R 4.2.2)
#>  GenomeInfoDb         * 1.34.7     2023-01-19 [2] Bioconductor
#>  GenomeInfoDbData       1.2.9      2022-11-08 [2] Bioconductor
#>  GenomicRanges        * 1.50.2     2023-01-19 [2] Bioconductor
#>  ggforce                0.4.1      2022-10-04 [2] CRAN (R 4.2.2)
#>  ggfun                  0.0.9      2022-11-21 [2] CRAN (R 4.2.2)
#>  ggplot2              * 3.4.0      2022-11-04 [2] CRAN (R 4.2.2)
#>  ggplotify              0.1.0      2021-09-02 [2] CRAN (R 4.2.2)
#>  ggraph                 2.1.0      2022-10-09 [2] CRAN (R 4.2.2)
#>  ggrepel                0.9.2      2022-11-06 [2] CRAN (R 4.2.2)
#>  ggtree                 3.6.2      2023-01-19 [2] Bioconductor
#>  glmnet                 4.1-6      2022-11-27 [2] CRAN (R 4.2.2)
#>  globals                0.16.2     2022-11-21 [2] CRAN (R 4.2.2)
#>  glue                   1.6.2      2022-02-24 [2] CRAN (R 4.2.2)
#>  GO.db                  3.16.0     2022-11-08 [2] Bioconductor
#>  GOSemSim               2.24.0     2023-01-19 [2] Bioconductor
#>  gower                  1.0.1      2022-12-22 [2] CRAN (R 4.2.2)
#>  graph                  1.76.0     2023-01-19 [2] Bioconductor
#>  graphite               1.44.0     2023-01-19 [2] Bioconductor
#>  graphlayouts           0.8.4      2022-11-24 [2] CRAN (R 4.2.2)
#>  gridExtra              2.3        2017-09-09 [2] CRAN (R 4.2.2)
#>  gridGraphics           0.5-1      2020-12-13 [2] CRAN (R 4.2.2)
#>  gson                   0.0.9      2022-09-06 [2] CRAN (R 4.2.2)
#>  gtable                 0.3.1      2022-09-01 [2] CRAN (R 4.2.2)
#>  hardhat                1.2.0      2022-06-30 [2] CRAN (R 4.2.2)
#>  HDO.db                 0.99.1     2022-11-08 [2] Bioconductor
#>  highr                  0.10       2022-12-22 [2] CRAN (R 4.2.2)
#>  htmltools              0.5.4      2022-12-07 [2] CRAN (R 4.2.2)
#>  htmlwidgets            1.6.1      2023-01-07 [2] CRAN (R 4.2.2)
#>  httpuv                 1.6.8      2023-01-12 [2] CRAN (R 4.2.2)
#>  httr                   1.4.4      2022-08-17 [2] CRAN (R 4.2.2)
#>  igraph                 1.3.5      2022-09-22 [2] CRAN (R 4.2.2)
#>  ipred                  0.9-13     2022-06-02 [2] CRAN (R 4.2.2)
#>  IRanges              * 2.32.0     2023-01-19 [2] Bioconductor
#>  iterators              1.0.14     2022-02-05 [2] CRAN (R 4.2.2)
#>  jquerylib              0.1.4      2021-04-26 [2] CRAN (R 4.2.2)
#>  jsonlite               1.8.4      2022-12-06 [2] CRAN (R 4.2.2)
#>  KEGGREST               1.38.0     2023-01-19 [2] Bioconductor
#>  knitr                  1.41       2022-11-18 [2] CRAN (R 4.2.2)
#>  labeling               0.4.2      2020-10-20 [2] CRAN (R 4.2.2)
#>  later                  1.3.0      2021-08-18 [2] CRAN (R 4.2.2)
#>  lattice              * 0.20-45    2021-09-22 [2] CRAN (R 4.2.2)
#>  lava                   1.7.1      2023-01-06 [2] CRAN (R 4.2.2)
#>  lazyeval               0.2.2      2019-03-15 [2] CRAN (R 4.2.2)
#>  lifecycle              1.0.3      2022-10-07 [2] CRAN (R 4.2.2)
#>  listenv                0.9.0      2022-12-16 [2] CRAN (R 4.2.2)
#>  locfit               * 1.5-9.7    2023-01-02 [2] CRAN (R 4.2.2)
#>  lubridate              1.9.0      2022-11-06 [2] CRAN (R 4.2.2)
#>  magrittr               2.0.3      2022-03-30 [2] CRAN (R 4.2.2)
#>  MASS                   7.3-58.1   2022-08-03 [2] CRAN (R 4.2.2)
#>  Matrix                 1.5-3      2022-11-11 [2] CRAN (R 4.2.2)
#>  MatrixGenerics       * 1.10.0     2023-01-19 [2] Bioconductor
#>  matrixStats          * 0.63.0     2022-11-18 [2] CRAN (R 4.2.2)
#>  memoise                2.0.1      2021-11-26 [2] CRAN (R 4.2.2)
#>  mime                   0.12       2021-09-28 [2] CRAN (R 4.2.2)
#>  miniUI                 0.1.1.1    2018-05-18 [2] CRAN (R 4.2.2)
#>  ModelMetrics           1.2.2.2    2020-03-17 [2] CRAN (R 4.2.2)
#>  munsell                0.5.0      2018-06-12 [2] CRAN (R 4.2.2)
#>  networkD3            * 0.4        2017-03-18 [2] CRAN (R 4.2.2)
#>  nlme                   3.1-161    2022-12-15 [2] CRAN (R 4.2.2)
#>  nnet                   7.3-18     2022-09-28 [2] CRAN (R 4.2.2)
#>  org.Hs.eg.db         * 3.16.0     2022-11-08 [2] Bioconductor
#>  parallelly             1.34.0     2023-01-13 [2] CRAN (R 4.2.2)
#>  patchwork              1.1.2      2022-08-19 [2] CRAN (R 4.2.2)
#>  pillar                 1.8.1      2022-08-19 [2] CRAN (R 4.2.2)
#>  pkgbuild               1.4.0      2022-11-27 [2] CRAN (R 4.2.2)
#>  pkgconfig              2.0.3      2019-09-22 [2] CRAN (R 4.2.2)
#>  pkgload                1.3.2      2022-11-16 [2] CRAN (R 4.2.2)
#>  plyr                   1.8.8      2022-11-11 [2] CRAN (R 4.2.2)
#>  png                    0.1-8      2022-11-29 [2] CRAN (R 4.2.2)
#>  polyclip               1.10-4     2022-10-20 [2] CRAN (R 4.2.2)
#>  prettyunits            1.1.1      2020-01-24 [2] CRAN (R 4.2.2)
#>  pROC                   1.18.0     2021-09-03 [2] CRAN (R 4.2.2)
#>  processx               3.8.0      2022-10-26 [2] CRAN (R 4.2.2)
#>  prodlim                2019.11.13 2019-11-17 [2] CRAN (R 4.2.2)
#>  profvis                0.3.7      2020-11-02 [2] CRAN (R 4.2.2)
#>  promises               1.2.0.1    2021-02-11 [2] CRAN (R 4.2.2)
#>  proxy                  0.4-27     2022-06-09 [2] CRAN (R 4.2.2)
#>  ps                     1.7.2      2022-10-26 [2] CRAN (R 4.2.2)
#>  purrr                  1.0.1      2023-01-10 [2] CRAN (R 4.2.2)
#>  qvalue                 2.30.0     2023-01-19 [2] Bioconductor
#>  R6                     2.5.1      2021-08-19 [2] CRAN (R 4.2.2)
#>  rappdirs               0.3.3      2021-01-31 [2] CRAN (R 4.2.2)
#>  RColorBrewer           1.1-3      2022-04-03 [2] CRAN (R 4.2.2)
#>  Rcpp                   1.0.9      2022-07-08 [2] CRAN (R 4.2.2)
#>  RcppArmadillo          0.11.4.3.1 2023-01-15 [2] CRAN (R 4.2.2)
#>  RcppParallel           5.1.6      2023-01-09 [2] CRAN (R 4.2.2)
#>  RCurl                  1.98-1.9   2022-10-03 [2] CRAN (R 4.2.2)
#>  reactome.db            1.82.0     2022-11-08 [2] Bioconductor
#>  ReactomePA           * 1.42.0     2023-01-19 [2] Bioconductor
#>  recipes                1.0.4      2023-01-11 [2] CRAN (R 4.2.2)
#>  remotes                2.4.2      2021-11-30 [2] CRAN (R 4.2.2)
#>  reshape2               1.4.4      2020-04-09 [2] CRAN (R 4.2.2)
#>  rlang                  1.0.6      2022-09-24 [2] CRAN (R 4.2.2)
#>  rmarkdown              2.19       2022-12-15 [2] CRAN (R 4.2.2)
#>  rpart                  4.1.19     2022-10-21 [2] CRAN (R 4.2.2)
#>  RSQLite                2.2.20     2022-12-22 [2] CRAN (R 4.2.2)
#>  Rtsne                  0.16       2022-04-17 [2] CRAN (R 4.2.2)
#>  S4Vectors            * 0.36.1     2023-01-19 [2] Bioconductor
#>  sass                   0.4.4      2022-11-24 [2] CRAN (R 4.2.2)
#>  scales                 1.2.1      2022-08-20 [2] CRAN (R 4.2.2)
#>  scatterpie             0.1.8      2022-09-03 [2] CRAN (R 4.2.2)
#>  scGPS                * 1.12.2     2023-01-19 [1] Bioconductor
#>  sessioninfo            1.2.2      2021-12-06 [2] CRAN (R 4.2.2)
#>  shadowtext             0.1.2      2022-04-22 [2] CRAN (R 4.2.2)
#>  shape                  1.4.6      2021-05-19 [2] CRAN (R 4.2.2)
#>  shiny                  1.7.4      2022-12-15 [2] CRAN (R 4.2.2)
#>  SingleCellExperiment * 1.20.0     2023-01-19 [2] Bioconductor
#>  stringi                1.7.12     2023-01-11 [2] CRAN (R 4.2.2)
#>  stringr                1.5.0      2022-12-02 [2] CRAN (R 4.2.2)
#>  SummarizedExperiment * 1.28.0     2023-01-19 [2] Bioconductor
#>  survival               3.5-0      2023-01-09 [2] CRAN (R 4.2.2)
#>  tibble                 3.1.8      2022-07-22 [2] CRAN (R 4.2.2)
#>  tidygraph              1.2.2      2022-08-22 [2] CRAN (R 4.2.2)
#>  tidyr                  1.2.1      2022-09-08 [2] CRAN (R 4.2.2)
#>  tidyselect             1.2.0      2022-10-10 [2] CRAN (R 4.2.2)
#>  tidytree               0.4.2      2022-12-18 [2] CRAN (R 4.2.2)
#>  timechange             0.2.0      2023-01-11 [2] CRAN (R 4.2.2)
#>  timeDate               4022.108   2023-01-07 [2] CRAN (R 4.2.2)
#>  treeio                 1.22.0     2023-01-19 [2] Bioconductor
#>  tweenr                 2.0.2      2022-09-06 [2] CRAN (R 4.2.2)
#>  urlchecker             1.0.1      2021-11-30 [2] CRAN (R 4.2.2)
#>  usethis                2.1.6      2022-05-25 [2] CRAN (R 4.2.2)
#>  utf8                   1.2.2      2021-07-24 [2] CRAN (R 4.2.2)
#>  vctrs                  0.5.1      2022-11-16 [2] CRAN (R 4.2.2)
#>  viridis                0.6.2      2021-10-13 [2] CRAN (R 4.2.2)
#>  viridisLite            0.4.1      2022-08-22 [2] CRAN (R 4.2.2)
#>  withr                  2.5.0      2022-03-03 [2] CRAN (R 4.2.2)
#>  xfun                   0.36       2022-12-21 [2] CRAN (R 4.2.2)
#>  XML                    3.99-0.13  2022-12-04 [2] CRAN (R 4.2.2)
#>  xtable                 1.8-4      2019-04-21 [2] CRAN (R 4.2.2)
#>  XVector                0.38.0     2023-01-19 [2] Bioconductor
#>  yaml                   2.3.6      2022-10-18 [2] CRAN (R 4.2.2)
#>  yulab.utils            0.0.6      2022-12-20 [2] CRAN (R 4.2.2)
#>  zlibbioc               1.44.0     2023-01-19 [2] Bioconductor
#> 
#>  [1] /tmp/Rtmpfm8KkJ/Rinst26133f400b4ec3
#>  [2] /home/biocbuild/bbs-3.16-bioc/R/library
#> 
#> ──────────────────────────────────────────────────────────────────────────────