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 will provide an example showing how users can use a pretrained
model of scClassify to predict cell types. A pretrained model is a
scClassifyTrainModel
object returned by train_scClassify()
.
A list of pretrained model can be found in
https://sydneybiox.github.io/scClassify/index.html.
First, install scClassify
, install BiocManager
and use
BiocManager::install
to install scClassify
package.
# installation of scClassify
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("scClassify")
We assume that you have log-transformed (size-factor normalized) matrices as query datasets, where each row refers to a gene and each column a cell. For demonstration purposes, we will take a subset of single-cell pancreas datasets from one independent study (Wang et al.).
library(scClassify)
data("scClassify_example")
wang_cellTypes <- scClassify_example$wang_cellTypes
exprsMat_wang_subset <- scClassify_example$exprsMat_wang_subset
exprsMat_wang_subset <- as(exprsMat_wang_subset, "dgCMatrix")
Here, we load our pretrained model using a subset of the Xin et al. human pancreas dataset as our reference data.
First, let us check basic information relating to our pretrained model.
data("trainClassExample_xin")
trainClassExample_xin
#> 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
In this pretrained model, we have selected the genes based on Differential Expression using limma. To check the genes that are available in the pretrained model:
features(trainClassExample_xin)
#> [1] "limma"
We can also visualise the cell type tree of the reference data.
plotCellTypeTree(cellTypeTree(trainClassExample_xin))
Next, we perform predict_scClassify
with our pretrained model
trainRes = trainClassExample
to predict the cell types of our
query data matrix exprsMat_wang_subset_sparse
. Here,
we used pearson
and spearman
as similarity metrics.
pred_res <- predict_scClassify(exprsMat_test = exprsMat_wang_subset,
trainRes = trainClassExample_xin,
cellTypes_test = wang_cellTypes,
algorithm = "WKNN",
features = c("limma"),
similarity = c("pearson", "spearman"),
prob_threshold = 0.7,
verbose = TRUE)
#> Performing unweighted ensemble learning...
#> Using parameters:
#> similarity algorithm features
#> "pearson" "WKNN" "limma"
#> [1] "Using dynamic correlation cutoff..."
#> [1] "Using dynamic correlation cutoff..."
#> classify_res
#> correct correctly unassigned intermediate
#> 0.704590818 0.239520958 0.000000000
#> incorrectly unassigned error assigned misclassified
#> 0.000000000 0.051896208 0.003992016
#> Using parameters:
#> similarity algorithm features
#> "spearman" "WKNN" "limma"
#> [1] "Using dynamic correlation cutoff..."
#> [1] "Using dynamic correlation cutoff..."
#> classify_res
#> correct correctly unassigned intermediate
#> 0.702594810 0.013972056 0.000000000
#> incorrectly unassigned error assigned misclassified
#> 0.001996008 0.277445110 0.003992016
#> weights for each base method:
#> [1] NA NA
Noted that the cellType_test
is not a required input.
For datasets with unknown labels, users can simply leave it
as cellType_test = NULL
.
Prediction results for pearson as the similarity metric:
table(pred_res$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
Prediction results for spearman as the similarity metric:
table(pred_res$spearman_WKNN_limma$predRes, wang_cellTypes)
#> wang_cellTypes
#> acinar alpha beta delta ductal gamma stellate
#> alpha 0 206 0 0 0 2 2
#> beta 2 0 118 0 29 0 6
#> beta_delta_gamma 1 0 0 0 66 0 31
#> delta 0 0 0 10 0 0 2
#> gamma 0 0 0 0 0 18 0
#> unassigned 2 0 0 0 1 1 4
sessionInfo()
#> R version 4.2.0 RC (2022-04-19 r82224)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.4 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_GB LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] scClassify_1.8.0 BiocStyle_2.24.0
#>
#> loaded via a namespace (and not attached):
#> [1] segmented_1.5-0 nlme_3.1-157
#> [3] bitops_1.0-7 matrixStats_0.62.0
#> [5] hopach_2.56.0 GenomeInfoDb_1.32.0
#> [7] tools_4.2.0 bslib_0.3.1
#> [9] utf8_1.2.2 R6_2.5.1
#> [11] HDF5Array_1.24.0 mgcv_1.8-40
#> [13] DBI_1.1.2 BiocGenerics_0.42.0
#> [15] colorspace_2.0-3 rhdf5filters_1.8.0
#> [17] gridExtra_2.3 tidyselect_1.1.2
#> [19] proxyC_0.2.4 compiler_4.2.0
#> [21] cli_3.3.0 Biobase_2.56.0
#> [23] DelayedArray_0.22.0 labeling_0.4.2
#> [25] bookdown_0.26 sass_0.4.1
#> [27] diptest_0.76-0 scales_1.2.0
#> [29] proxy_0.4-26 stringr_1.4.0
#> [31] digest_0.6.29 mixtools_1.2.0
#> [33] rmarkdown_2.14 XVector_0.36.0
#> [35] pkgconfig_2.0.3 htmltools_0.5.2
#> [37] sparseMatrixStats_1.8.0 Cepo_1.2.0
#> [39] MatrixGenerics_1.8.0 highr_0.9
#> [41] fastmap_1.1.0 limma_3.52.0
#> [43] rlang_1.0.2 DelayedMatrixStats_1.18.0
#> [45] jquerylib_0.1.4 generics_0.1.2
#> [47] farver_2.1.0 jsonlite_1.8.0
#> [49] BiocParallel_1.30.0 dplyr_1.0.8
#> [51] RCurl_1.98-1.6 magrittr_2.0.3
#> [53] GenomeInfoDbData_1.2.8 patchwork_1.1.1
#> [55] Matrix_1.4-1 Rcpp_1.0.8.3
#> [57] munsell_0.5.0 S4Vectors_0.34.0
#> [59] Rhdf5lib_1.18.0 fansi_1.0.3
#> [61] viridis_0.6.2 lifecycle_1.0.1
#> [63] stringi_1.7.6 yaml_2.3.5
#> [65] ggraph_2.0.5 MASS_7.3-57
#> [67] SummarizedExperiment_1.26.0 zlibbioc_1.42.0
#> [69] rhdf5_2.40.0 plyr_1.8.7
#> [71] grid_4.2.0 parallel_4.2.0
#> [73] ggrepel_0.9.1 crayon_1.5.1
#> [75] lattice_0.20-45 splines_4.2.0
#> [77] graphlayouts_0.8.0 magick_2.7.3
#> [79] knitr_1.38 pillar_1.7.0
#> [81] igraph_1.3.1 GenomicRanges_1.48.0
#> [83] reshape2_1.4.4 stats4_4.2.0
#> [85] glue_1.6.2 evaluate_0.15
#> [87] RcppParallel_5.1.5 BiocManager_1.30.17
#> [89] vctrs_0.4.1 tweenr_1.0.2
#> [91] gtable_0.3.0 purrr_0.3.4
#> [93] polyclip_1.10-0 tidyr_1.2.0
#> [95] kernlab_0.9-30 assertthat_0.2.1
#> [97] ggplot2_3.3.5 xfun_0.30
#> [99] ggforce_0.3.3 tidygraph_1.2.1
#> [101] survival_3.3-1 viridisLite_0.4.0
#> [103] minpack.lm_1.2-2 SingleCellExperiment_1.18.0
#> [105] tibble_3.1.6 IRanges_2.30.0
#> [107] cluster_2.1.3 statmod_1.4.36
#> [109] ellipsis_0.3.2