Contents

1 Introduction

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:

  1. Prepare a single-cell reference; i.e. CellTypeDataset (CTD). Alternatively, you can use one of the pre-generated CTDs we provide via the package ewceData (which comes with EWCE).
  2. Run cell type enrichment on a gene list using the bootstrap_enrichment_test function.

NOTE: This documentation is for the development version of EWCE. See Bioconductor for documentation on the current release version.

2 Setup

library(EWCE) 
## Loading required package: RNOmni
set.seed(1234)

#### Package name ####
pkg <- tolower("EWCE")
#### Username of DockerHub account ####
docker_user <- "neurogenomicslab"

3 Run cell-type enrichment tests

3.1 1. Prepare input data

3.1.1 CellTypeDataset

Load a CTD previously generated from mouse cortex and hypothalamus single-cell RNA-seq data (from the Karolinska Institute).

3.1.1.1 CTD levels

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()
## see ?ewceData and browseVignettes('ewceData') for documentation
## loading from cache

3.1.1.2 Plot CTD mean_exp

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")

3.1.2 Gene list

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"

3.2 2. Run cell type enrichment tests

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.

3.2.0.1 Hyperparameters

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.

3.2.0.2 Parallelisation

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
## Gathering ortholog reports.
## 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.
## 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.
## --
## --
## 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.
## Gathering ortholog reports.
## 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.
## 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
## Checking gene list inputs.
## Running without gene size control.
## 17 hit gene(s) remain after filtering.
## Computing gene scores.
## Using previously sampled genes.
## Computing gene counts.
## 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)
print(plot_list$withDendro)
## NULL

4 Docker

ewce is now available via DockerHub as a containerised environment with Rstudio and all necessary dependencies pre-installed.

4.1 Installation

4.2 Method 1: via Docker

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
  • The -d ensures the container will run in “detached” mode, which means it will persist even after you’ve closed your command line session.
  • Optionally, you can also install the Docker Desktop to easily manage your containers.
  • You can set the password to whatever you like by changing the -e PASSWORD=... flag.
  • The username will be “rstudio” by default.

4.3 Method 2: via Singularity

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

4.4 Usage

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.

5 Session Info

utils::sessionInfo()
## R Under development (unstable) (2023-11-11 r85510)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] ewceData_1.11.0      ExperimentHub_2.11.0 AnnotationHub_3.11.0
## [4] BiocFileCache_2.11.1 dbplyr_2.4.0         BiocGenerics_0.49.1 
## [7] EWCE_1.11.3          RNOmni_1.0.1.2       BiocStyle_2.31.0    
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.7                magrittr_2.0.3               
##   [3] magick_2.8.1                  farver_2.1.1                 
##   [5] rmarkdown_2.25                fs_1.6.3                     
##   [7] zlibbioc_1.49.0               vctrs_0.6.4                  
##   [9] memoise_2.0.1                 RCurl_1.98-1.13              
##  [11] ggtree_3.11.0                 rstatix_0.7.2                
##  [13] htmltools_0.5.7               S4Arrays_1.3.0               
##  [15] curl_5.1.0                    broom_1.0.5                  
##  [17] SparseArray_1.3.1             gridGraphics_0.5-1           
##  [19] sass_0.4.7                    bslib_0.5.1                  
##  [21] htmlwidgets_1.6.2             HGNChelper_0.8.1             
##  [23] plyr_1.8.9                    plotly_4.10.3                
##  [25] cachem_1.0.8                  mime_0.12                    
##  [27] lifecycle_1.0.4               pkgconfig_2.0.3              
##  [29] Matrix_1.6-3                  R6_2.5.1                     
##  [31] fastmap_1.1.1                 GenomeInfoDbData_1.2.11      
##  [33] MatrixGenerics_1.15.0         shiny_1.8.0                  
##  [35] digest_0.6.33                 aplot_0.2.2                  
##  [37] colorspace_2.1-0              patchwork_1.1.3              
##  [39] AnnotationDbi_1.65.2          S4Vectors_0.41.1             
##  [41] grr_0.9.5                     GenomicRanges_1.55.1         
##  [43] RSQLite_2.3.3                 ggpubr_0.6.0                 
##  [45] labeling_0.4.3                filelock_1.0.2               
##  [47] fansi_1.0.5                   httr_1.4.7                   
##  [49] abind_1.4-5                   compiler_4.4.0               
##  [51] withr_2.5.2                   bit64_4.0.5                  
##  [53] backports_1.4.1               BiocParallel_1.37.0          
##  [55] orthogene_1.9.0               carData_3.0-5                
##  [57] DBI_1.1.3                     homologene_1.4.68.19.3.27    
##  [59] highr_0.10                    ggsignif_0.6.4               
##  [61] rappdirs_0.3.3                DelayedArray_0.29.0          
##  [63] tools_4.4.0                   ape_5.7-1                    
##  [65] interactiveDisplayBase_1.41.0 httpuv_1.6.12                
##  [67] glue_1.6.2                    nlme_3.1-163                 
##  [69] promises_1.2.1                grid_4.4.0                   
##  [71] reshape2_1.4.4                generics_0.1.3               
##  [73] gtable_0.3.4                  tidyr_1.3.0                  
##  [75] data.table_1.14.8             car_3.1-2                    
##  [77] utf8_1.2.4                    XVector_0.43.0               
##  [79] BiocVersion_3.19.1            pillar_1.9.0                 
##  [81] stringr_1.5.1                 yulab.utils_0.1.0            
##  [83] babelgene_22.9                limma_3.59.1                 
##  [85] later_1.3.1                   dplyr_1.1.4                  
##  [87] treeio_1.27.0                 lattice_0.22-5               
##  [89] bit_4.0.5                     tidyselect_1.2.0             
##  [91] SingleCellExperiment_1.25.0   Biostrings_2.71.1            
##  [93] knitr_1.45                    bookdown_0.36                
##  [95] IRanges_2.37.0                SummarizedExperiment_1.33.0  
##  [97] stats4_4.4.0                  xfun_0.41                    
##  [99] Biobase_2.63.0                statmod_1.5.0                
## [101] matrixStats_1.1.0             stringi_1.8.1                
## [103] lazyeval_0.2.2                ggfun_0.1.3                  
## [105] yaml_2.3.7                    codetools_0.2-19             
## [107] evaluate_0.23                 tibble_3.2.1                 
## [109] BiocManager_1.30.22           ggplotify_0.1.2              
## [111] cli_3.6.1                     xtable_1.8-4                 
## [113] munsell_0.5.0                 jquerylib_0.1.4              
## [115] Rcpp_1.0.11                   GenomeInfoDb_1.39.1          
## [117] gprofiler2_0.2.2              png_0.1-8                    
## [119] parallel_4.4.0                ellipsis_0.3.2               
## [121] ggplot2_3.4.4                 blob_1.2.4                   
## [123] bitops_1.0-7                  viridisLite_0.4.2            
## [125] tidytree_0.4.5                scales_1.2.1                 
## [127] purrr_1.0.2                   crayon_1.5.2                 
## [129] rlang_1.1.2                   KEGGREST_1.43.0

6 References