The EasyCellType
package was designed to examine an input marker list using
the databases and provide annotation recommendations in graphical outcomes.
The package refers to 3 public available marker gene data bases,
and provides two approaches to conduct the annotation anaysis:
gene set enrichment analysis(GSEA) and a modified Fisher’s exact test.
The package has been submitted to bioconductor
to achieve an easy access for researchers.
This vignette shows a simple workflow illustrating how EasyCellType package works. The data set that will be used throughout the example is freely available from 10X Genomics.
The package can be installed using BiocManager
by the following commands
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("EasyCellType")
Alternatively, the package can also be installed using devtools
and launched by
library(devtools)
install_github("rx-li/EasyCellType")
After the installation, the package can be loaded with
library(EasyCellType)
We use the Peripheral Blood Mononuclear Cells (PBMC) data freely available from 10X Genomics.
The data can be downladed from
https://cf.10xgenomics.com/samples/cell/pbmc3k/pbmc3k_filtered_gene_bc_matrices.tar.gz.
After downloading the data, it can be read using function Read10X
.
We have included the data in our package, which can be loaded with
data(pbmc_data)
We followed the standard workflow provided by Seurat
package(Hao et al. 2021) to process the PBMC data set.
The detailed technical explanations can be found in
https://satijalab.org/seurat/articles/pbmc3k_tutorial.html.
library(Seurat)
# Initialize the Seurat object
pbmc <- CreateSeuratObject(counts = pbmc_data, project = "pbmc3k", min.cells = 3, min.features = 200)
# QC and select samples
pbmc[["percent.mt"]] <- PercentageFeatureSet(pbmc, pattern = "^MT-")
pbmc <- subset(pbmc, subset = nFeature_RNA > 200 & nFeature_RNA < 2500 & percent.mt < 5)
# Normalize the data
pbmc <- NormalizeData(pbmc)
# Identify highly variable features
pbmc <- FindVariableFeatures(pbmc, selection.method = "vst", nfeatures = 2000)
# Scale the data
all.genes <- rownames(pbmc)
pbmc <- ScaleData(pbmc, features = all.genes)
# Perfom linear dimensional reduction
pbmc <- RunPCA(pbmc, features = VariableFeatures(object = pbmc))
# Cluster the cells
pbmc <- FindNeighbors(pbmc, dims = 1:10)
pbmc <- FindClusters(pbmc, resolution = 0.5)
# Find differentially expressed features
markers <- FindAllMarkers(pbmc, only.pos = TRUE, min.pct = 0.25, logfc.threshold = 0.25)
Now we get the expressed markers for each cluster. We then convert the gene symbols to Entrez IDs.
library(org.Hs.eg.db)
library(AnnotationDbi)
markers$entrezid <- mapIds(org.Hs.eg.db,
keys=markers$gene, #Column containing Ensembl gene ids
column="ENTREZID",
keytype="SYMBOL",
multiVals="first")
markers <- na.omit(markers)
In case the data is measured in mouse, we would replace the package org.Hs.eg.db
with org.Mm.eg.db
and do the above analysis.
The input for EasyCellType
package should be a data frame containing Entrez IDs,
clusters and expression scores. The order of columns should follow this rule.
In each cluster, the gene should be sorted by the expression score.
library(dplyr)
markers_sort <- data.frame(gene=markers$entrezid, cluster=markers$cluster,
score=markers$avg_log2FC) %>%
group_by(cluster) %>%
mutate(rank = rank(score), ties.method = "random") %>%
arrange(desc(rank))
input.d <- as.data.frame(markers_sort[, 1:3])
We have include the processed data in the package. It can be loaded with
data("gene_pbmc")
input.d <- gene_pbmc
Now we can call the annot
function to run annotation analysis.
annot.GSEA <- easyct(input.d, db="cellmarker", species="Human",
tissue=c("Blood", "Peripheral blood", "Blood vessel",
"Umbilical cord blood", "Venous blood"), p_cut=0.3,
test="GSEA")
We used the GSEA approach to do the annotation. In our package, we use GSEA
function in clusterProfiler
package(Wu et al. 2021) to conduct the enrichment analysis.
You can replace ‘GSEA’ with ‘fisher’ if you would like to use Fisher exact test
to do the annotation.
The candidate tissues can be seen using data(cellmarker_tissue)
,
data(clustermole_tissue)
and data(panglao_tissue)
.
The dot plot showing the overall annotation results can be created by
plot_dot(test="GSEA", annot.GSEA)
Bar plot can be created by
plot_bar(test="GSEA", annot.GSEA)
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: /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
#>
#> 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
#>
#> time zone: America/New_York
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] dplyr_1.1.4 org.Hs.eg.db_3.20.0 AnnotationDbi_1.68.0
#> [4] IRanges_2.40.0 S4Vectors_0.44.0 Biobase_2.66.0
#> [7] BiocGenerics_0.52.0 Seurat_5.1.0 SeuratObject_5.0.2
#> [10] sp_2.1-4 EasyCellType_1.5.4 devtools_2.4.5
#> [13] usethis_3.0.0 BiocStyle_2.34.0
#>
#> loaded via a namespace (and not attached):
#> [1] fs_1.6.4 matrixStats_1.4.1 spatstat.sparse_3.1-0
#> [4] enrichplot_1.26.0 httr_1.4.7 RColorBrewer_1.1-3
#> [7] profvis_0.4.0 tools_4.4.1 sctransform_0.4.1
#> [10] utf8_1.2.4 R6_2.5.1 lazyeval_0.2.2
#> [13] uwot_0.2.2 urlchecker_1.0.1 withr_3.0.2
#> [16] gridExtra_2.3 progressr_0.15.0 cli_3.6.3
#> [19] spatstat.explore_3.3-3 fastDummies_1.7.4 sass_0.4.9
#> [22] spatstat.data_3.1-2 ggridges_0.5.6 pbapply_1.7-2
#> [25] yulab.utils_0.1.7 gson_0.1.0 DOSE_4.0.0
#> [28] R.utils_2.12.3 parallelly_1.38.0 sessioninfo_1.2.2
#> [31] limma_3.62.0 RSQLite_2.3.7 generics_0.1.3
#> [34] gridGraphics_0.5-1 ica_1.0-3 spatstat.random_3.3-2
#> [37] GO.db_3.20.0 Matrix_1.7-1 fansi_1.0.6
#> [40] abind_1.4-8 R.methodsS3_1.8.2 lifecycle_1.0.4
#> [43] yaml_2.3.10 qvalue_2.38.0 Rtsne_0.17
#> [46] grid_4.4.1 blob_1.2.4 promises_1.3.0
#> [49] crayon_1.5.3 miniUI_0.1.1.1 ggtangle_0.0.3
#> [52] lattice_0.22-6 cowplot_1.1.3 KEGGREST_1.46.0
#> [55] magick_2.8.5 pillar_1.9.0 knitr_1.48
#> [58] fgsea_1.32.0 future.apply_1.11.3 codetools_0.2-20
#> [61] fastmatch_1.1-4 leiden_0.4.3.1 glue_1.8.0
#> [64] ggfun_0.1.7 spatstat.univar_3.0-1 data.table_1.16.2
#> [67] remotes_2.5.0 vctrs_0.6.5 png_0.1-8
#> [70] treeio_1.30.0 spam_2.11-0 org.Mm.eg.db_3.20.0
#> [73] gtable_0.3.6 cachem_1.1.0 xfun_0.48
#> [76] mime_0.12 survival_3.7-0 tinytex_0.53
#> [79] statmod_1.5.0 ellipsis_0.3.2 fitdistrplus_1.2-1
#> [82] ROCR_1.0-11 nlme_3.1-166 ggtree_3.14.0
#> [85] bit64_4.5.2 RcppAnnoy_0.0.22 GenomeInfoDb_1.42.0
#> [88] bslib_0.8.0 irlba_2.3.5.1 KernSmooth_2.23-24
#> [91] colorspace_2.1-1 DBI_1.2.3 tidyselect_1.2.1
#> [94] processx_3.8.4 bit_4.5.0 compiler_4.4.1
#> [97] curl_5.2.3 desc_1.4.3 plotly_4.10.4
#> [100] bookdown_0.41 scales_1.3.0 lmtest_0.9-40
#> [103] callr_3.7.6 stringr_1.5.1 digest_0.6.37
#> [106] goftest_1.2-3 spatstat.utils_3.1-0 rmarkdown_2.28
#> [109] XVector_0.46.0 htmltools_0.5.8.1 pkgconfig_2.0.3
#> [112] highr_0.11 fastmap_1.2.0 rlang_1.1.4
#> [115] htmlwidgets_1.6.4 UCSC.utils_1.2.0 shiny_1.9.1
#> [118] farver_2.1.2 jquerylib_0.1.4 zoo_1.8-12
#> [121] jsonlite_1.8.9 BiocParallel_1.40.0 GOSemSim_2.32.0
#> [124] R.oo_1.26.0 magrittr_2.0.3 GenomeInfoDbData_1.2.13
#> [127] ggplotify_0.1.2 dotCall64_1.2 patchwork_1.3.0
#> [130] munsell_0.5.1 Rcpp_1.0.13 ape_5.8
#> [133] reticulate_1.39.0 stringi_1.8.4 zlibbioc_1.52.0
#> [136] MASS_7.3-61 plyr_1.8.9 pkgbuild_1.4.5
#> [139] parallel_4.4.1 listenv_0.9.1 ggrepel_0.9.6
#> [142] forcats_1.0.0 deldir_2.0-4 Biostrings_2.74.0
#> [145] splines_4.4.1 tensor_1.5 ps_1.8.1
#> [148] igraph_2.1.1 spatstat.geom_3.3-3 RcppHNSW_0.6.0
#> [151] reshape2_1.4.4 pkgload_1.4.0 evaluate_1.0.1
#> [154] BiocManager_1.30.25 httpuv_1.6.15 RANN_2.6.2
#> [157] tidyr_1.3.1 purrr_1.0.2 polyclip_1.10-7
#> [160] future_1.34.0 scattermore_1.2 ggplot2_3.5.1
#> [163] xtable_1.8-4 RSpectra_0.16-2 tidytree_0.4.6
#> [166] later_1.3.2 viridisLite_0.4.2 tibble_3.2.1
#> [169] clusterProfiler_4.14.0 aplot_0.2.3 memoise_2.0.1
#> [172] cluster_2.1.6 globals_0.16.3
Hao, Yuhan, Stephanie Hao, Erica Andersen-Nissen, William M. Mauck III, Shiwei Zheng, Andrew Butler, Maddie J. Lee, et al. 2021. “Integrated Analysis of Multimodal Single-Cell Data.” Cell. https://doi.org/10.1016/j.cell.2021.04.048.
Wu, Tianzhi, Erqiang Hu, Shuangbin Xu, Meijun Chen, Pingfan Guo, Zehan Dai, Tingze Feng, et al. 2021. “ClusterProfiler 4.0: A Universal Enrichment Tool for Interpreting Omics Data.” The Innovation. https://doi.org/10.1016/j.xinn.2021.100141.