scGPS introduction

Quan Nguyen and Michael Thompson

2020-04-28

1. Installation instruction

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

2. A simple workflow of the scGPS:

The purpose of this workflow is to solve the following task:

2.1 Create scGPS objects

2.2 Run prediction

2.3 Summarise results

3. A complete workflow of the scGPS:

The purpose of this workflow is to solve the following task:

3.1 Identify clusters in a dataset using CORE

(skip this step if clusters are known)

3.2 Identify clusters in a dataset using SCORE (Stable Clustering at Optimal REsolution)

(skip this step if clusters are known)

(SCORE aims to get stable subpopulation results by introducing bagging aggregation and bootstrapping to the CORE algorithm)

3.3 Visualise all cluster results in all iterations

3.4 Compare clustering results with other dimensional reduction methods (e.g., tSNE)

3.5 Find gene markers and annotate clusters

#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))
suppressMessages(library(DESeq))

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] NA                      NA                      NA                     
#> [4] NA                      "DE_Subpop1vsRemaining" "DE_Subpop2vsRemaining"
#> [7] "DE_Subpop3vsRemaining" "DE_Subpop4vsRemaining"

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

4. Relationship between clusters within one sample or between two samples

The purpose of this workflow is to solve the following task:

4.1 Start the scGPS prediction to find relationship between clusters

4.2 Display summary results for the prediction

#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 94.1176470588235 1.06951871657754 LDA for subpop 1 in target mixedpop2
#> 2 98.5714285714286 71.4285714285714 LDA for subpop 2 in target mixedpop2
#> 3 90.2255639097744 3.00751879699248 LDA for subpop 3 in target mixedpop2
#> 4               95               20 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 90.3743315508021 43.8502673796791 ElasticNet for subpop1 in target mixedpop2
#> 2 97.8571428571428 97.1428571428571 ElasticNet for subpop2 in target mixedpop2
#> 3 92.4812030075188 71.4285714285714 ElasticNet for subpop3 in target mixedpop2
#> 4              100             72.5 ElasticNet for subpop4 in target mixedpop2

# summary maximum deviance explained by the model during the model training
summary_deviance(object = LSOLDA_dat)
#> $allDeviance
#> [1] "0.709"  "0.5711"
#> 
#> $DeviMax
#>           dat_DE$Dfd          Deviance           DEgenes
#> 1                  0             0.709    genes_cluster1
#> 2                  1             0.709    genes_cluster1
#> 3                  2             0.709    genes_cluster1
#> 4                  3             0.709    genes_cluster1
#> 5                  6             0.709    genes_cluster1
#> 6                  7             0.709    genes_cluster1
#> 7                 11             0.709    genes_cluster1
#> 8                 13             0.709    genes_cluster1
#> 9                 16             0.709    genes_cluster1
#> 10                18             0.709    genes_cluster1
#> 11                20             0.709    genes_cluster1
#> 12                22             0.709    genes_cluster1
#> 13                24             0.709    genes_cluster1
#> 14                25             0.709    genes_cluster1
#> 15                29             0.709    genes_cluster1
#> 16                32             0.709    genes_cluster1
#> 17                35             0.709    genes_cluster1
#> 18                39             0.709    genes_cluster1
#> 19                40             0.709    genes_cluster1
#> 20                43             0.709    genes_cluster1
#> 21                44             0.709    genes_cluster1
#> 22                48             0.709    genes_cluster1
#> 23 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] 70.98214 66.51786

4.3 Plot the relationship between clusters in one sample

Here 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$DE_Subpop1vsRemaining$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$DE_Subpop2vsRemaining$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$DE_Subpop3vsRemaining$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$DE_Subpop4vsRemaining$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 )
p1
C2_MP2 → C2_MP2
99 
C4_MP2 → C4_MP2
90 
C1_MP2 → C1_MP2
78 
C4_MP2 → C2_MP2
57 
C3_MP2 → C1_MP2
52 
C4_MP2 → C1_MP2
48 
C4_MP2 → C3_MP2
45 
C2_MP2 → C4_MP2
35 
C3_MP2 → C3_MP2
33 
C1_MP2 → C4_MP2
16 
C2_MP2 → C3_MP2
11 
C3_MP2 → C4_MP2
6 
C2_MP2 → C1_MP2
6 
C1_MP2 → C3_MP2
4 
C1_MP2 → C2_MP2
0 
C3_MP2 → C2_MP2
0 
C1_MP2
184
C1_MP2
C2_MP2
156
C2_MP2
C3_MP2
92
C3_MP2
C4_MP2
240
C4_MP2

4.3 Plot the relationship between clusters in two samples

Here 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$DE_Subpop1vsRemaining$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$DE_Subpop2vsRemaining$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$DE_Subpop3vsRemaining$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)
p1
C3_MP1 → C1_MP2
57 
C1_MP1 → C2_MP2
56 
C1_MP1 → C4_MP2
50 
C2_MP1 → C3_MP2
47 
C1_MP1 → C1_MP2
47 
C2_MP1 → C4_MP2
41 
C2_MP1 → C1_MP2
35 
C3_MP1 → C3_MP2
33 
C1_MP1 → C3_MP2
26 
C3_MP1 → C4_MP2
21 
C3_MP1 → C2_MP2
17 
C2_MP1 → C2_MP2
14 
C1_MP1
179
C1_MP1
C2_MP1
137
C2_MP1
C3_MP1
128
C3_MP1
C1_MP2
139
C1_MP2
C2_MP2
86
C2_MP2
C3_MP2
106
C3_MP2
C4_MP2
113
C4_MP2
devtools::session_info()
#> ─ Session info ───────────────────────────────────────────────────────────────
#>  setting  value                       
#>  version  R version 4.0.0 (2020-04-24)
#>  os       Ubuntu 18.04.4 LTS          
#>  system   x86_64, linux-gnu           
#>  ui       X11                         
#>  language (EN)                        
#>  collate  C                           
#>  ctype    en_US.UTF-8                 
#>  tz       America/New_York            
#>  date     2020-04-28                  
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────
#>  package              * version     date       lib source        
#>  annotate               1.66.0      2020-04-27 [2] Bioconductor  
#>  AnnotationDbi        * 1.50.0      2020-04-27 [2] Bioconductor  
#>  assertthat             0.2.1       2019-03-21 [2] CRAN (R 4.0.0)
#>  backports              1.1.6       2020-04-05 [2] CRAN (R 4.0.0)
#>  Biobase              * 2.48.0      2020-04-27 [2] Bioconductor  
#>  BiocGenerics         * 0.34.0      2020-04-27 [2] Bioconductor  
#>  BiocManager            1.30.10     2019-11-16 [2] CRAN (R 4.0.0)
#>  BiocParallel           1.22.0      2020-04-27 [2] Bioconductor  
#>  bit                    1.1-15.2    2020-02-10 [2] CRAN (R 4.0.0)
#>  bit64                  0.9-7       2017-05-08 [2] CRAN (R 4.0.0)
#>  bitops                 1.0-6       2013-08-17 [2] CRAN (R 4.0.0)
#>  blob                   1.2.1       2020-01-20 [2] CRAN (R 4.0.0)
#>  callr                  3.4.3       2020-03-28 [2] CRAN (R 4.0.0)
#>  caret                * 6.0-86      2020-03-20 [2] CRAN (R 4.0.0)
#>  checkmate              2.0.0       2020-02-06 [2] CRAN (R 4.0.0)
#>  class                  7.3-17      2020-04-26 [2] CRAN (R 4.0.0)
#>  cli                    2.0.2       2020-02-28 [2] CRAN (R 4.0.0)
#>  clusterProfiler      * 3.16.0      2020-04-27 [2] Bioconductor  
#>  codetools              0.2-16      2018-12-24 [2] CRAN (R 4.0.0)
#>  colorspace             1.4-1       2019-03-18 [2] CRAN (R 4.0.0)
#>  cowplot                1.0.0       2019-07-11 [2] CRAN (R 4.0.0)
#>  crayon                 1.3.4       2017-09-16 [2] CRAN (R 4.0.0)
#>  data.table             1.12.8      2019-12-09 [2] CRAN (R 4.0.0)
#>  DBI                    1.1.0       2019-12-15 [2] CRAN (R 4.0.0)
#>  DelayedArray         * 0.14.0      2020-04-27 [2] Bioconductor  
#>  dendextend             1.13.4      2020-02-28 [2] CRAN (R 4.0.0)
#>  desc                   1.2.0       2018-05-01 [2] CRAN (R 4.0.0)
#>  DESeq                * 1.40.0      2020-04-27 [2] Bioconductor  
#>  devtools               2.3.0       2020-04-10 [2] CRAN (R 4.0.0)
#>  digest                 0.6.25      2020-02-23 [2] CRAN (R 4.0.0)
#>  DO.db                  2.9         2020-04-24 [2] Bioconductor  
#>  DOSE                 * 3.14.0      2020-04-27 [2] Bioconductor  
#>  downloader             0.4         2015-07-09 [2] CRAN (R 4.0.0)
#>  dplyr                * 0.8.5       2020-03-07 [2] CRAN (R 4.0.0)
#>  dynamicTreeCut       * 1.63-1      2016-03-11 [2] CRAN (R 4.0.0)
#>  e1071                  1.7-3       2019-11-26 [2] CRAN (R 4.0.0)
#>  ellipsis               0.3.0       2019-09-20 [2] CRAN (R 4.0.0)
#>  enrichplot             1.8.0       2020-04-27 [2] Bioconductor  
#>  europepmc              0.3         2018-04-20 [2] CRAN (R 4.0.0)
#>  evaluate               0.14        2019-05-28 [2] CRAN (R 4.0.0)
#>  fansi                  0.4.1       2020-01-08 [2] CRAN (R 4.0.0)
#>  farver                 2.0.3       2020-01-16 [2] CRAN (R 4.0.0)
#>  fastcluster            1.1.25      2018-06-07 [2] CRAN (R 4.0.0)
#>  fastmatch              1.1-0       2017-01-28 [2] CRAN (R 4.0.0)
#>  fgsea                  1.14.0      2020-04-27 [2] Bioconductor  
#>  foreach                1.5.0       2020-03-30 [2] CRAN (R 4.0.0)
#>  fs                     1.4.1       2020-04-04 [2] CRAN (R 4.0.0)
#>  genefilter             1.70.0      2020-04-27 [2] Bioconductor  
#>  geneplotter            1.66.0      2020-04-27 [2] Bioconductor  
#>  generics               0.0.2       2018-11-29 [2] CRAN (R 4.0.0)
#>  GenomeInfoDb         * 1.24.0      2020-04-27 [2] Bioconductor  
#>  GenomeInfoDbData       1.2.3       2020-04-24 [2] Bioconductor  
#>  GenomicRanges        * 1.40.0      2020-04-27 [2] Bioconductor  
#>  ggforce                0.3.1       2019-08-20 [2] CRAN (R 4.0.0)
#>  ggplot2              * 3.3.0       2020-03-05 [2] CRAN (R 4.0.0)
#>  ggplotify              0.0.5       2020-03-12 [2] CRAN (R 4.0.0)
#>  ggraph                 2.0.2       2020-03-17 [2] CRAN (R 4.0.0)
#>  ggrepel                0.8.2       2020-03-08 [2] CRAN (R 4.0.0)
#>  ggridges               0.5.2       2020-01-12 [2] CRAN (R 4.0.0)
#>  glmnet                 3.0-2       2019-12-11 [2] CRAN (R 4.0.0)
#>  glue                   1.4.0       2020-04-03 [2] CRAN (R 4.0.0)
#>  GO.db                  3.10.0      2020-04-24 [2] Bioconductor  
#>  GOSemSim               2.14.0      2020-04-27 [2] Bioconductor  
#>  gower                  0.2.1       2019-05-14 [2] CRAN (R 4.0.0)
#>  graph                  1.66.0      2020-04-27 [2] Bioconductor  
#>  graphite               1.34.0      2020-04-27 [2] Bioconductor  
#>  graphlayouts           0.7.0       2020-04-25 [2] CRAN (R 4.0.0)
#>  gridExtra              2.3         2017-09-09 [2] CRAN (R 4.0.0)
#>  gridGraphics           0.5-0       2020-02-25 [2] CRAN (R 4.0.0)
#>  gtable                 0.3.0       2019-03-25 [2] CRAN (R 4.0.0)
#>  hms                    0.5.3       2020-01-08 [2] CRAN (R 4.0.0)
#>  htmltools              0.4.0       2019-10-04 [2] CRAN (R 4.0.0)
#>  htmlwidgets            1.5.1       2019-10-08 [2] CRAN (R 4.0.0)
#>  httr                   1.4.1       2019-08-05 [2] CRAN (R 4.0.0)
#>  igraph                 1.2.5       2020-03-19 [2] CRAN (R 4.0.0)
#>  ipred                  0.9-9       2019-04-28 [2] CRAN (R 4.0.0)
#>  IRanges              * 2.22.0      2020-04-27 [2] Bioconductor  
#>  iterators              1.0.12      2019-07-26 [2] CRAN (R 4.0.0)
#>  jsonlite               1.6.1       2020-02-02 [2] CRAN (R 4.0.0)
#>  knitr                  1.28        2020-02-06 [2] CRAN (R 4.0.0)
#>  labeling               0.3         2014-08-23 [2] CRAN (R 4.0.0)
#>  lattice              * 0.20-41     2020-04-02 [2] CRAN (R 4.0.0)
#>  lava                   1.6.7       2020-03-05 [2] CRAN (R 4.0.0)
#>  lifecycle              0.2.0       2020-03-06 [2] CRAN (R 4.0.0)
#>  locfit               * 1.5-9.4     2020-03-25 [2] CRAN (R 4.0.0)
#>  lubridate              1.7.8       2020-04-06 [2] CRAN (R 4.0.0)
#>  magrittr               1.5         2014-11-22 [2] CRAN (R 4.0.0)
#>  MASS                   7.3-51.6    2020-04-26 [2] CRAN (R 4.0.0)
#>  Matrix                 1.2-18      2019-11-27 [2] CRAN (R 4.0.0)
#>  matrixStats          * 0.56.0      2020-03-13 [2] CRAN (R 4.0.0)
#>  memoise                1.1.0       2017-04-21 [2] CRAN (R 4.0.0)
#>  ModelMetrics           1.2.2.2     2020-03-17 [2] CRAN (R 4.0.0)
#>  munsell                0.5.0       2018-06-12 [2] CRAN (R 4.0.0)
#>  networkD3            * 0.4         2017-03-18 [2] CRAN (R 4.0.0)
#>  nlme                   3.1-147     2020-04-13 [2] CRAN (R 4.0.0)
#>  nnet                   7.3-14      2020-04-26 [2] CRAN (R 4.0.0)
#>  org.Hs.eg.db         * 3.10.0      2020-04-24 [2] Bioconductor  
#>  pillar                 1.4.3       2019-12-20 [2] CRAN (R 4.0.0)
#>  pkgbuild               1.0.7       2020-04-25 [2] CRAN (R 4.0.0)
#>  pkgconfig              2.0.3       2019-09-22 [2] CRAN (R 4.0.0)
#>  pkgload                1.0.2       2018-10-29 [2] CRAN (R 4.0.0)
#>  plyr                   1.8.6       2020-03-03 [2] CRAN (R 4.0.0)
#>  polyclip               1.10-0      2019-03-14 [2] CRAN (R 4.0.0)
#>  prettyunits            1.1.1       2020-01-24 [2] CRAN (R 4.0.0)
#>  pROC                   1.16.2      2020-03-19 [2] CRAN (R 4.0.0)
#>  processx               3.4.2       2020-02-09 [2] CRAN (R 4.0.0)
#>  prodlim                2019.11.13  2019-11-17 [2] CRAN (R 4.0.0)
#>  progress               1.2.2       2019-05-16 [2] CRAN (R 4.0.0)
#>  ps                     1.3.2       2020-02-13 [2] CRAN (R 4.0.0)
#>  purrr                  0.3.4       2020-04-17 [2] CRAN (R 4.0.0)
#>  qvalue                 2.20.0      2020-04-27 [2] Bioconductor  
#>  R6                     2.4.1       2019-11-12 [2] CRAN (R 4.0.0)
#>  rappdirs               0.3.1       2016-03-28 [2] CRAN (R 4.0.0)
#>  RColorBrewer           1.1-2       2014-12-07 [2] CRAN (R 4.0.0)
#>  Rcpp                   1.0.4.6     2020-04-09 [2] CRAN (R 4.0.0)
#>  RcppArmadillo          0.9.870.2.0 2020-04-27 [2] CRAN (R 4.0.0)
#>  RcppParallel           5.0.0       2020-03-11 [2] CRAN (R 4.0.0)
#>  RCurl                  1.98-1.2    2020-04-18 [2] CRAN (R 4.0.0)
#>  reactome.db            1.70.0      2020-04-24 [2] Bioconductor  
#>  ReactomePA           * 1.32.0      2020-04-27 [2] Bioconductor  
#>  recipes                0.1.10      2020-03-18 [2] CRAN (R 4.0.0)
#>  remotes                2.1.1       2020-02-15 [2] CRAN (R 4.0.0)
#>  reshape2               1.4.4       2020-04-09 [2] CRAN (R 4.0.0)
#>  rlang                  0.4.5       2020-03-01 [2] CRAN (R 4.0.0)
#>  rmarkdown              2.1         2020-01-20 [2] CRAN (R 4.0.0)
#>  rpart                  4.1-15      2019-04-12 [2] CRAN (R 4.0.0)
#>  rprojroot              1.3-2       2018-01-03 [2] CRAN (R 4.0.0)
#>  RSQLite                2.2.0       2020-01-07 [2] CRAN (R 4.0.0)
#>  Rtsne                  0.15        2018-11-10 [2] CRAN (R 4.0.0)
#>  rvcheck                0.1.8       2020-03-01 [2] CRAN (R 4.0.0)
#>  S4Vectors            * 0.26.0      2020-04-27 [2] Bioconductor  
#>  scales                 1.1.0       2019-11-18 [2] CRAN (R 4.0.0)
#>  scatterpie             0.1.4       2019-11-08 [2] CRAN (R 4.0.0)
#>  scGPS                * 1.2.0       2020-04-27 [1] Bioconductor  
#>  sessioninfo            1.1.1       2018-11-05 [2] CRAN (R 4.0.0)
#>  shape                  1.4.4       2018-02-07 [2] CRAN (R 4.0.0)
#>  SingleCellExperiment * 1.10.0      2020-04-27 [2] Bioconductor  
#>  stringi                1.4.6       2020-02-17 [2] CRAN (R 4.0.0)
#>  stringr                1.4.0       2019-02-10 [2] CRAN (R 4.0.0)
#>  SummarizedExperiment * 1.18.0      2020-04-27 [2] Bioconductor  
#>  survival               3.1-12      2020-04-10 [2] CRAN (R 4.0.0)
#>  testthat               2.3.2       2020-03-02 [2] CRAN (R 4.0.0)
#>  tibble                 3.0.1       2020-04-20 [2] CRAN (R 4.0.0)
#>  tidygraph              1.1.2       2019-02-18 [2] CRAN (R 4.0.0)
#>  tidyr                  1.0.2       2020-01-24 [2] CRAN (R 4.0.0)
#>  tidyselect             1.0.0       2020-01-27 [2] CRAN (R 4.0.0)
#>  timeDate               3043.102    2018-02-21 [2] CRAN (R 4.0.0)
#>  triebeard              0.3.0       2016-08-04 [2] CRAN (R 4.0.0)
#>  tweenr                 1.0.1       2018-12-14 [2] CRAN (R 4.0.0)
#>  urltools               1.7.3       2019-04-14 [2] CRAN (R 4.0.0)
#>  usethis                1.6.0       2020-04-09 [2] CRAN (R 4.0.0)
#>  vctrs                  0.2.4       2020-03-10 [2] CRAN (R 4.0.0)
#>  viridis                0.5.1       2018-03-29 [2] CRAN (R 4.0.0)
#>  viridisLite            0.3.0       2018-02-01 [2] CRAN (R 4.0.0)
#>  withr                  2.2.0       2020-04-20 [2] CRAN (R 4.0.0)
#>  xfun                   0.13        2020-04-13 [2] CRAN (R 4.0.0)
#>  XML                    3.99-0.3    2020-01-20 [2] CRAN (R 4.0.0)
#>  xml2                   1.3.2       2020-04-23 [2] CRAN (R 4.0.0)
#>  xtable                 1.8-4       2019-04-21 [2] CRAN (R 4.0.0)
#>  XVector                0.28.0      2020-04-27 [2] Bioconductor  
#>  yaml                   2.2.1       2020-02-01 [2] CRAN (R 4.0.0)
#>  zlibbioc               1.34.0      2020-04-27 [2] Bioconductor  
#> 
#> [1] /tmp/Rtmp0myVxh/Rinst2d0c1948fb50
#> [2] /home/biocbuild/bbs-3.11-bioc/R/library