Introduction

A common application of single-cell RNA sequencing (RNA-seq) data is to identify discrete cell types. To take advantage of the large collection of well-annotated scRNA-seq datasets, scClassify package implements a set of methods to perform accurate cell type classification based on ensemble learning and sample size calculation. This vignette demonstrates the usage of scClassify, providing a pithy description of each method with workable examples.

First, install scClassify via BiocManager.

# installation of scClassify
if (!requireNamespace("BiocManager", quietly = TRUE)) {
    install.packages("BiocManager")
}
BiocManager::install("scClassify")

Setting up the data

We assume that you have log-transformed (size-factor normalized) matrices where each row is a gene and each column a cell for a reference dataset and a query dataset. For demonstration purposes, we will take a subset of single-cell pancreas datasets from two independent studies (Wang et al., and Xin et al.).

library("scClassify")
data("scClassify_example")
xin_cellTypes <- scClassify_example$xin_cellTypes
exprsMat_xin_subset <- scClassify_example$exprsMat_xin_subset
wang_cellTypes <- scClassify_example$wang_cellTypes
exprsMat_wang_subset <- scClassify_example$exprsMat_wang_subset
exprsMat_xin_subset <- as(exprsMat_xin_subset, "dgCMatrix")
exprsMat_wang_subset <- as(exprsMat_wang_subset, "dgCMatrix")

The original cell type annotations and compositions of the example datasets can be easily accessed as shown below.

table(xin_cellTypes)
#> xin_cellTypes
#> alpha  beta delta gamma 
#>   285   261    49    79
table(wang_cellTypes)
#> wang_cellTypes
#>   acinar    alpha     beta    delta   ductal    gamma stellate 
#>        5      206      118       10       96       21       45

We can see that Xin et al. data only have 4 cell types, while Wang et al. has 7 cell types.

scClassify

Non-ensemble scClassify

We first perform non-ensemble scClassify by using Xin et al. as our reference dataset and Wang et al. data as ur query data. We use WKNN as the KNN algorithm, DE (differential expression genes) as the gene selection method, and lastly pearson as the similarity calculation method.

scClassify_res <- scClassify(exprsMat_train = exprsMat_xin_subset,
                             cellTypes_train = xin_cellTypes,
                             exprsMat_test = list(wang = exprsMat_wang_subset),
                             cellTypes_test = list(wang = wang_cellTypes),
                             tree = "HOPACH",
                             algorithm = "WKNN",
                             selectFeatures = c("limma"),
                             similarity = c("pearson"),
                             returnList = FALSE,
                             verbose = FALSE)

We can check the cell type tree generated by the reference data:

scClassify_res$trainRes
#> Class: scClassifyTrainModel 
#> Model name: training 
#> Feature selection methods: limma 
#> Number of cells in the training data: 674 
#> Number of cell types in the training data: 4
plotCellTypeTree(cellTypeTree(scClassify_res$trainRes))

Noted that scClassify_res$trainRes is a scClassifyTrainModel class.

Check the prediction results.

table(scClassify_res$testRes$wang$pearson_WKNN_limma$predRes, wang_cellTypes)
#>                   wang_cellTypes
#>                    acinar alpha beta delta ductal gamma stellate
#>   alpha                 0   206    0     0      0     2        0
#>   beta                  0     0  118     0      1     0        0
#>   beta_delta_gamma      0     0    0     0     25     0        0
#>   delta                 0     0    0    10      0     0        0
#>   gamma                 0     0    0     0      0    19        0
#>   unassigned            5     0    0     0     70     0       45

Ensemble Classify

We next perform ensemble scClassify by using Xin et al. as our reference dataset and Wang et al. data as our query data. We use WKNN as the KNN algorithm, DE as the gene selection method, and pearson and spearman as the similarity calculation methods. Thus, we will generate two combinations of gene selection models and similarity metrics as training classifiers:

  1. WKNN + DE + pearson
  2. WKNN + DE + spearman

Here, we will weight these two classifiers equally by setting weighted_ensemble = FALSE. By default this is set as TRUE, so each base classifier will be weighted by the accuracy rates trained in the reference data.

scClassify_res_ensemble <- scClassify(exprsMat_train = exprsMat_xin_subset,
                                      cellTypes_train = xin_cellTypes,
                                      exprsMat_test = list(wang = exprsMat_wang_subset),
                                      cellTypes_test = list(wang = wang_cellTypes),
                                      tree = "HOPACH",
                                      algorithm = "WKNN",
                                      selectFeatures = c("limma"),
                                      similarity = c("pearson", "cosine"),
                                      weighted_ensemble = FALSE,
                                      returnList = FALSE,
                                      verbose = FALSE)

We can compare the two base classifiers predictions as below.

table(scClassify_res_ensemble$testRes$wang$pearson_WKNN_limma$predRes,
      scClassify_res_ensemble$testRes$wang$cosine_WKNN_limma$predRes)
#>                   
#>                    alpha beta beta_delta_gamma delta gamma unassigned
#>   alpha              208    0                0     0     0          0
#>   beta                 0  119                0     0     0          0
#>   beta_delta_gamma     0    0               22     0     0          3
#>   delta                0    0                0    10     0          0
#>   gamma                0    0                0     0    19          0
#>   unassigned           1    3               15     0     0        101

Now, check the final ensemble results:

table(scClassify_res_ensemble$testRes$wang$ensembleRes$cellTypes,
      wang_cellTypes)
#>                   wang_cellTypes
#>                    acinar alpha beta delta ductal gamma stellate
#>   alpha                 0   206    0     0      0     2        0
#>   beta                  0     0  118     0      1     0        0
#>   beta_delta_gamma      0     0    0     0     25     0        0
#>   delta                 0     0    0    10      0     0        0
#>   gamma                 0     0    0     0      0    19        0
#>   unassigned            5     0    0     0     70     0       45

Train your own model

You can also train your own model scClassifyTrainModel using train_scClassify(). Note that by setting weightsCal = TRUE, we will calculate the training error of the reference data as the weights for the individual classifiers.

Here, we illustrate the training function with gene selection methods based on differential expression (“limma”) and biomodal distribution (“BI”).

trainClass <- train_scClassify(exprsMat_train = exprsMat_xin_subset,
                               cellTypes_train = xin_cellTypes,
                               selectFeatures = c("limma", "BI"),
                               returnList = FALSE
)
#> after filtering not expressed genes 
#> [1] 980 674
#> [1] "Feature Selection..."
#> [1] "Number of genes selected to construct HOPACH tree 160"
#> [1] "Constructing tree ..."
#> [1] "Training...."
#> [1] "=== selecting features by: limma ===="
#> [1] "=== selecting features by: BI ===="
trainClass
#> Class: scClassifyTrainModel 
#> Model name: training 
#> Feature selection methods: limma BI 
#> Number of cells in the training data: 674 
#> Number of cell types in the training data: 4

Session Info

sessionInfo()
#> R Under development (unstable) (2020-03-25 r78063)
#> Platform: x86_64-apple-darwin15.6.0 (64-bit)
#> Running under: macOS Catalina 10.15.4
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.0/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.0/Resources/lib/libRlapack.dylib
#> 
#> locale:
#> [1] en_AU.UTF-8/en_AU.UTF-8/en_AU.UTF-8/C/en_AU.UTF-8/en_AU.UTF-8
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] scClassify_0.99.2 BiocStyle_2.15.6 
#> 
#> loaded via a namespace (and not attached):
#>  [1] Biobase_2.47.3      viridis_0.5.1       mixtools_1.2.0     
#>  [4] tidyr_1.0.2         tidygraph_1.1.2     viridisLite_0.3.0  
#>  [7] splines_4.0.0       ggraph_2.0.2        RcppParallel_5.0.0 
#> [10] assertthat_0.2.1    statmod_1.4.34      BiocManager_1.30.10
#> [13] stats4_4.0.0        yaml_2.2.1          ggrepel_0.8.2      
#> [16] pillar_1.4.3        backports_1.1.6     lattice_0.20-41    
#> [19] glue_1.4.0          limma_3.43.5        digest_0.6.25      
#> [22] polyclip_1.10-0     colorspace_1.4-1    htmltools_0.4.0    
#> [25] Matrix_1.2-18       pkgconfig_2.0.3     bookdown_0.18      
#> [28] purrr_0.3.3         scales_1.1.0        tweenr_1.0.1       
#> [31] hopach_2.47.0       BiocParallel_1.21.2 ggforce_0.3.1      
#> [34] proxy_0.4-23        tibble_3.0.0        mgcv_1.8-31        
#> [37] farver_2.0.3        ggplot2_3.3.0       ellipsis_0.3.0     
#> [40] BiocGenerics_0.33.3 cli_2.0.2           survival_3.1-12    
#> [43] magrittr_1.5        crayon_1.3.4        memoise_1.1.0      
#> [46] evaluate_0.14       fs_1.4.1            fansi_0.4.1        
#> [49] nlme_3.1-145        MASS_7.3-51.5       segmented_1.1-0    
#> [52] tools_4.0.0         minpack.lm_1.2-1    lifecycle_0.2.0    
#> [55] stringr_1.4.0       kernlab_0.9-29      S4Vectors_0.25.15  
#> [58] munsell_0.5.0       cluster_2.1.0       compiler_4.0.0     
#> [61] pkgdown_1.5.1.9000  proxyC_0.1.5        rlang_0.4.5        
#> [64] grid_4.0.0          rstudioapi_0.11     igraph_1.2.5       
#> [67] labeling_0.3        rmarkdown_2.1       gtable_0.3.0       
#> [70] graphlayouts_0.6.0  R6_2.4.1            gridExtra_2.3      
#> [73] knitr_1.28          dplyr_0.8.5         rprojroot_1.3-2    
#> [76] desc_1.2.0          stringi_1.4.6       parallel_4.0.0     
#> [79] Rcpp_1.0.4.6        vctrs_0.2.4         diptest_0.75-7     
#> [82] tidyselect_1.0.0    xfun_0.12