The EWCE
R package is designed to facilitate expression weighted cell type
enrichment analysis as described in our Frontiers in Neuroscience
paper (???). EWCE
can be applied to any gene list.
Using EWCE
essentially involves two steps:
ewceData
(which comes with EWCE
).bootstrap_enrichment_test
function.NOTE: This documentation is for the development version of EWCE
. See Bioconductor for documentation on the current release version.
library(EWCE)
## Loading required package: RNOmni
set.seed(1234)
#### Package name ####
pkg <- tolower("EWCE")
#### Username of DockerHub account ####
docker_user <- "neurogenomicslab"
Load a CTD previously generated from mouse cortex and hypothalamus single-cell RNA-seq data (from the Karolinska Institute).
Each level of a CTD corresponds to increasingly refined cell-type/-subtype annotations. For example, in the CTD ewceData::ctd()
level 1 includes the cell-type “interneurons”, while level 2 breaks these this group into 16 different interneuron subtypes (“Int…”).
ctd <- ewceData::ctd()
## snapshotDate(): 2022-10-24
## see ?ewceData and browseVignettes('ewceData') for documentation
## loading from cache
Plot the expression of four markers genes across all cell types in the CTD.
plt_exp <- EWCE::plot_ctd(ctd = ctd,
level = 1,
genes = c("Apoe","Gfap","Gapdh"),
metric = "mean_exp")
plt_spec <- EWCE::plot_ctd(ctd = ctd,
level = 2,
genes = c("Apoe","Gfap","Gapdh"),
metric = "specificity")
Gene lists input into EWCE can comes from any source (e.g. GWAS, candidate genes, pathways).
Here, we provide an example gene list of Alzheimer’s disease-related nominated from a GWAS.
hits <- ewceData::example_genelist()
## see ?ewceData and browseVignettes('ewceData') for documentation
## loading from cache
print(hits)
## [1] "APOE" "BIN1" "CLU" "ABCA7" "CR1" "PICALM"
## [7] "MS4A6A" "CD33" "MS4A4E" "CD2AP" "EOGA1" "INPP5D"
## [13] "MEF2C" "HLA-DRB5" "ZCWPW1" "NME8" "PTK2B" "CELF1"
## [19] "SORL1" "FERMT2" "SLC24A4" "CASS4"
We now run the cell type enrichment tests on the gene list. Since the CTD is from mouse data (and is annotated using mouse genes) we specify the argument sctSpecies="mouse"
. bootstrap_enrichment_test
will automaticlaly convert the mouse genes to human genes.
Since the gene list came from GWAS in humans, we set genelistSpecies="human"
.
Note: We set the seed at the top of this vignette to ensure reproducibility in the bootstrap sampling function.
Note: We use 100 repetitions here for the purposes of a quick example, but in practice you would want to use reps=10000
for publishable results.
You can now speed up the bootstrapping process by parallelising across
multiple cores with the parameter no_cores
(=1
by default).
reps <- 100
annotLevel <- 1
full_results <- EWCE::bootstrap_enrichment_test(sct_data = ctd,
sctSpecies = "mouse",
genelistSpecies = "human",
hits = hits,
reps = reps,
annotLevel = annotLevel)
## 1 core(s) assigned as workers (71 reserved).
## Generating gene background for mouse x human ==> human
## Retrieving all genes using: homologene.
## Retrieving all organisms available in homologene.
## Mapping species name: mouse
## Common name mapping found for mouse
## 1 organism identified from search: 10090
## Gene table with 21,207 rows retrieved.
## Returning all 21,207 genes from mouse.
## --
## Retrieving all genes using: homologene.
## Retrieving all organisms available in homologene.
## Mapping species name: human
## Common name mapping found for human
## 1 organism identified from search: 9606
## Gene table with 19,129 rows retrieved.
## Returning all 19,129 genes from human.
## --
## Preparing gene_df.
## data.frame format detected.
## Extracting genes from Gene.Symbol.
## 21,207 genes extracted.
## Converting mouse ==> human orthologs using: homologene
## Retrieving all organisms available in homologene.
## Mapping species name: mouse
## Common name mapping found for mouse
## 1 organism identified from search: 10090
## Retrieving all organisms available in homologene.
## Mapping species name: human
## Common name mapping found for human
## 1 organism identified from search: 9606
## Checking for genes without orthologs in human.
## Extracting genes from input_gene.
## 17,355 genes extracted.
## Extracting genes from ortholog_gene.
## 17,355 genes extracted.
## Checking for genes without 1:1 orthologs.
## Dropping 131 genes that have multiple input_gene per ortholog_gene (many:1).
## Dropping 498 genes that have multiple ortholog_gene per input_gene (1:many).
## Filtering gene_df with gene_map
## Adding input_gene col to gene_df.
## Adding ortholog_gene col to gene_df.
##
## =========== REPORT SUMMARY ===========
## Total genes dropped after convert_orthologs :
## 4,725 / 21,207 (22%)
## Total genes remaining after convert_orthologs :
## 16,482 / 21,207 (78%)
## --
##
## =========== REPORT SUMMARY ===========
## 16,482 / 21,207 (77.72%) target_species genes remain after ortholog conversion.
## 16,482 / 19,129 (86.16%) reference_species genes remain after ortholog conversion.
## Retrieving all genes using: homologene.
## Retrieving all organisms available in homologene.
## Mapping species name: human
## Common name mapping found for human
## 1 organism identified from search: 9606
## Gene table with 19,129 rows retrieved.
## Returning all 19,129 genes from human.
## --
##
## =========== REPORT SUMMARY ===========
## 19,129 / 19,129 (100%) target_species genes remain after ortholog conversion.
## 19,129 / 19,129 (100%) reference_species genes remain after ortholog conversion.
## 16,482 intersect background genes used.
## Standardising CellTypeDataset
## Converting to sparse matrix.
## Converting to sparse matrix.
## Checking gene list inputs.
## Returning 16,482 unique genes from the user-supplied bg.
## Standardising sct_data.
## Converting gene list input to standardised human genes.
## Running without gene size control.
## 17 hit genes remain after filtering.
## Testing for enrichment in 7 cell types...
## Sorting results by p-value.
## Computing BH-corrected q-values.
## 1 significant cell type enrichment results @ q<0.05 :
## CellType annotLevel p fold_change sd_from_mean q
## 1 microglia 1 0 1.965915 3.938119 0
The main table of results is stored in full_results$results
.
In this case, microglia were the only cell type that was significantly enriched in the Alzheimer’s disease gene list.
knitr::kable(full_results$results)
CellType | annotLevel | p | fold_change | sd_from_mean | q | |
---|---|---|---|---|---|---|
microglia | microglia | 1 | 0.00 | 1.9659148 | 3.9381188 | 0.000 |
astrocytes_ependymal | astrocytes_ependymal | 1 | 0.13 | 1.2624889 | 1.1553910 | 0.455 |
pyramidal_SS | pyramidal_SS | 1 | 0.80 | 0.8699242 | -0.8226268 | 1.000 |
oligodendrocytes | oligodendrocytes | 1 | 0.87 | 0.7631149 | -1.0861761 | 1.000 |
pyramidal_CA1 | pyramidal_CA1 | 1 | 0.89 | 0.8202496 | -1.1738063 | 1.000 |
endothelial_mural | endothelial_mural | 1 | 0.90 | 0.7674534 | -1.1811797 | 1.000 |
interneurons | interneurons | 1 | 1.00 | 0.4012954 | -3.4703413 | 1.000 |
The results can be visualised using another function, which shows for each cell type, the number of standard deviations from the mean the level of expression was found to be in the target gene list, relative to the bootstrapped mean.
The dendrogram at the top shows how the cell types are hierarchically clustered by transcriptional similarity.
plot_list <- EWCE::ewce_plot(total_res = full_results$results,
mtc_method = "BH",
ctd = ctd)
## Loading required namespace: cowplot
## Loading required namespace: gridExtra
## Scale for 'x' is already present. Adding another scale for 'x', which will
## replace the existing scale.
print(plot_list$withDendro)
ewce is now available via DockerHub as a containerised environment with Rstudio and all necessary dependencies pre-installed.
First, install Docker if you have not already.
Create an image of the Docker container in command line:
docker pull neurogenomicslab/ewce
Once the image has been created, you can launch it with:
docker run \
-d \
-e ROOT=true \
-e PASSWORD=bioc \
-v ~/Desktop:/Desktop \
-v /Volumes:/Volumes \
-p 8787:8787 \
neurogenomicslab/ewce
-d
ensures the container will run in “detached” mode,
which means it will persist even after you’ve closed your command line session.-e PASSWORD=...
flag.If you are using a system that does not allow Docker (as is the case for many institutional computing clusters), you can instead install Docker images via Singularity.
singularity pull docker://neurogenomicslab/ewce
Finally, launch the containerised Rstudio by entering the following URL in any web browser: http://localhost:8787/
Login using the credentials set during the Installation steps.
utils::sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-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] ewceData_1.5.0 ExperimentHub_2.6.0 AnnotationHub_3.6.0
## [4] BiocFileCache_2.6.0 dbplyr_2.2.1 BiocGenerics_0.44.0
## [7] EWCE_1.6.0 RNOmni_1.0.1 BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] colorspace_2.0-3 ggtree_3.6.0
## [3] grr_0.9.5 ggsignif_0.6.4
## [5] ellipsis_0.3.2 XVector_0.38.0
## [7] GenomicRanges_1.50.0 aplot_0.1.8
## [9] farver_2.1.1 ggpubr_0.4.0
## [11] bit64_4.0.5 interactiveDisplayBase_1.36.0
## [13] AnnotationDbi_1.60.0 fansi_1.0.3
## [15] orthogene_1.4.0 codetools_0.2-18
## [17] cachem_1.0.6 knitr_1.40
## [19] jsonlite_1.8.3 broom_1.0.1
## [21] png_0.1-7 shiny_1.7.3
## [23] BiocManager_1.30.19 compiler_4.2.1
## [25] httr_1.4.4 backports_1.4.1
## [27] lazyeval_0.2.2 assertthat_0.2.1
## [29] Matrix_1.5-1 fastmap_1.1.0
## [31] limma_3.54.0 cli_3.4.1
## [33] later_1.3.0 htmltools_0.5.3
## [35] tools_4.2.1 gtable_0.3.1
## [37] glue_1.6.2 GenomeInfoDbData_1.2.9
## [39] reshape2_1.4.4 dplyr_1.0.10
## [41] rappdirs_0.3.3 Rcpp_1.0.9
## [43] carData_3.0-5 Biobase_2.58.0
## [45] jquerylib_0.1.4 vctrs_0.5.0
## [47] Biostrings_2.66.0 babelgene_22.9
## [49] ape_5.6-2 nlme_3.1-160
## [51] xfun_0.34 stringr_1.4.1
## [53] mime_0.12 lifecycle_1.0.3
## [55] rstatix_0.7.0 zlibbioc_1.44.0
## [57] scales_1.2.1 promises_1.2.0.1
## [59] MatrixGenerics_1.10.0 parallel_4.2.1
## [61] SummarizedExperiment_1.28.0 gprofiler2_0.2.1
## [63] SingleCellExperiment_1.20.0 yaml_2.3.6
## [65] curl_4.3.3 gridExtra_2.3
## [67] memoise_2.0.1 ggplot2_3.3.6
## [69] ggfun_0.0.7 yulab.utils_0.0.5
## [71] sass_0.4.2 stringi_1.7.8
## [73] RSQLite_2.2.18 highr_0.9
## [75] BiocVersion_3.16.0 S4Vectors_0.36.0
## [77] tidytree_0.4.1 filelock_1.0.2
## [79] BiocParallel_1.32.0 GenomeInfoDb_1.34.0
## [81] rlang_1.0.6 pkgconfig_2.0.3
## [83] matrixStats_0.62.0 bitops_1.0-7
## [85] evaluate_0.17 lattice_0.20-45
## [87] purrr_0.3.5 labeling_0.4.2
## [89] htmlwidgets_1.5.4 treeio_1.22.0
## [91] patchwork_1.1.2 cowplot_1.1.1
## [93] bit_4.0.4 tidyselect_1.2.0
## [95] plyr_1.8.7 magrittr_2.0.3
## [97] bookdown_0.29 R6_2.5.1
## [99] magick_2.7.3 IRanges_2.32.0
## [101] generics_0.1.3 DelayedArray_0.24.0
## [103] DBI_1.1.3 withr_2.5.0
## [105] pillar_1.8.1 KEGGREST_1.38.0
## [107] abind_1.4-5 RCurl_1.98-1.9
## [109] tibble_3.1.8 homologene_1.4.68.19.3.27
## [111] crayon_1.5.2 car_3.1-1
## [113] utf8_1.2.2 plotly_4.10.0
## [115] rmarkdown_2.17 grid_4.2.1
## [117] data.table_1.14.4 blob_1.2.3
## [119] digest_0.6.30 xtable_1.8-4
## [121] HGNChelper_0.8.1 tidyr_1.2.1
## [123] httpuv_1.6.6 gridGraphics_0.5-1
## [125] stats4_4.2.1 munsell_0.5.0
## [127] viridisLite_0.4.1 ggplotify_0.1.0
## [129] bslib_0.4.0