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.
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:
We can also visualise the cell type tree of the reference data.
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:
sessionInfo()
#> R version 4.4.1 (2024-06-14)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] scClassify_1.17.0 BiocStyle_2.33.1
#>
#> loaded via a namespace (and not attached):
#> [1] gridExtra_2.3 rlang_1.1.4
#> [3] magrittr_2.0.3 matrixStats_1.3.0
#> [5] compiler_4.4.1 mgcv_1.9-1
#> [7] DelayedMatrixStats_1.27.3 vctrs_0.6.5
#> [9] reshape2_1.4.4 stringr_1.5.1
#> [11] pkgconfig_2.0.3 crayon_1.5.3
#> [13] fastmap_1.2.0 XVector_0.45.0
#> [15] labeling_0.4.3 ggraph_2.2.1
#> [17] utf8_1.2.4 rmarkdown_2.28
#> [19] UCSC.utils_1.1.0 purrr_1.0.2
#> [21] xfun_0.47 zlibbioc_1.51.1
#> [23] cachem_1.1.0 GenomeInfoDb_1.41.1
#> [25] jsonlite_1.8.8 highr_0.11
#> [27] rhdf5filters_1.17.0 DelayedArray_0.31.11
#> [29] Rhdf5lib_1.27.0 BiocParallel_1.39.0
#> [31] tweenr_2.0.3 parallel_4.4.1
#> [33] cluster_2.1.6 R6_2.5.1
#> [35] bslib_0.8.0 stringi_1.8.4
#> [37] limma_3.61.9 diptest_0.77-1
#> [39] GenomicRanges_1.57.1 jquerylib_0.1.4
#> [41] Rcpp_1.0.13 SummarizedExperiment_1.35.1
#> [43] knitr_1.48 mixtools_2.0.0
#> [45] IRanges_2.39.2 Matrix_1.7-0
#> [47] splines_4.4.1 igraph_2.0.3
#> [49] tidyselect_1.2.1 abind_1.4-5
#> [51] yaml_2.3.10 hopach_2.65.0
#> [53] viridis_0.6.5 minpack.lm_1.2-4
#> [55] codetools_0.2-20 Cepo_1.11.2
#> [57] lattice_0.22-6 tibble_3.2.1
#> [59] plyr_1.8.9 Biobase_2.65.1
#> [61] withr_3.0.1 evaluate_0.24.0
#> [63] survival_3.7-0 proxy_0.4-27
#> [65] polyclip_1.10-7 kernlab_0.9-33
#> [67] pillar_1.9.0 BiocManager_1.30.25
#> [69] MatrixGenerics_1.17.0 stats4_4.4.1
#> [71] plotly_4.10.4 generics_0.1.3
#> [73] S4Vectors_0.43.2 ggplot2_3.5.1
#> [75] sparseMatrixStats_1.17.2 munsell_0.5.1
#> [77] scales_1.3.0 glue_1.7.0
#> [79] lazyeval_0.2.2 proxyC_0.4.1
#> [81] maketools_1.3.0 tools_4.4.1
#> [83] data.table_1.16.0 sys_3.4.2
#> [85] buildtools_1.0.0 graphlayouts_1.1.1
#> [87] tidygraph_1.3.1 rhdf5_2.49.0
#> [89] grid_4.4.1 tidyr_1.3.1
#> [91] colorspace_2.1-1 SingleCellExperiment_1.27.2
#> [93] nlme_3.1-166 GenomeInfoDbData_1.2.12
#> [95] patchwork_1.2.0 ggforce_0.4.2
#> [97] HDF5Array_1.33.6 cli_3.6.3
#> [99] fansi_1.0.6 segmented_2.1-1
#> [101] S4Arrays_1.5.7 viridisLite_0.4.2
#> [103] dplyr_1.1.4 gtable_0.3.5
#> [105] sass_0.4.9 digest_0.6.37
#> [107] BiocGenerics_0.51.0 SparseArray_1.5.31
#> [109] ggrepel_0.9.5 htmlwidgets_1.6.4
#> [111] farver_2.1.2 memoise_2.0.1
#> [113] htmltools_0.5.8.1 lifecycle_1.0.4
#> [115] httr_1.4.7 statmod_1.5.0
#> [117] MASS_7.3-61