The ReactomeGSA package is a client to the web-based Reactome Analysis System. Essentially, it performs a gene set analysis using the latest version of the Reactome pathway database as a backend.
This vignette shows how the ReactomeGSA package can be used to perform a pathway analysis of cell clusters in single-cell RNA-sequencing data.
To cite this package, use
Griss J. ReactomeGSA, https://github.com/reactome/ReactomeGSA (2019)
The ReactomeGSA
package can be directly installed from Bioconductor:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
if (!require(ReactomeGSA))
BiocManager::install("ReactomeGSA")
#> Loading required package: ReactomeGSA
# install the ReactomeGSA.data package for the example data
if (!require(ReactomeGSA.data))
BiocManager::install("ReactomeGSA.data")
#> Loading required package: ReactomeGSA.data
#> Loading required package: limma
#> Loading required package: edgeR
#> Loading required package: Seurat
#> Attaching SeuratObject
For more information, see https://bioconductor.org/install/.
As an example we load single-cell RNA-sequencing data of B cells extracted from the dataset published by Jerby-Arnon et al. (Cell, 2018).
Note: This is not a complete Seurat object. To decrease the size, the object only contains gene expression values and cluster annotations.
The pathway analysis is at the very end of a scRNA-seq workflow. This means, that any Q/C was already performed, the data was normalized and cells were already clustered.
The ReactomeGSA package can now be used to get pathway-level expression values for every cell cluster. This is achieved by calculating the mean gene expression for every cluster and then submitting this data to a gene set variation analysis.
All of this is wrapped in the single analyse_sc_clusters
function.
library(ReactomeGSA)
gsva_result <- analyse_sc_clusters(jerby_b_cells, verbose = TRUE)
#> Calculating average cluster expression...
#> Converting expression data to string... (This may take a moment)
#> Conversion complete
#> Submitting request to Reactome API...
#> Compressing request data...
#> Reactome Analysis submitted succesfully
#> Converting dataset Seurat...
#> Mapping identifiers...
#> Performing gene set analysis using ssGSEA
#> Analysing dataset 'Seurat' using ssGSEA
#> Retrieving result...
The resulting object is a standard ReactomeAnalysisResult
object.
gsva_result
#> ReactomeAnalysisResult object
#> Reactome Release: 77
#> Results:
#> - Seurat:
#> 1728 pathways
#> 10899 fold changes for genes
#> No Reactome visualizations available
#> ReactomeAnalysisResult
pathways
returns the pathway-level expression values per cell cluster:
pathway_expression <- pathways(gsva_result)
# simplify the column names by removing the default dataset identifier
colnames(pathway_expression) <- gsub("\\.Seurat", "", colnames(pathway_expression))
pathway_expression[1:3,]
#> Name Cluster_1 Cluster_10 Cluster_11
#> R-HSA-1059683 Interleukin-6 signaling 0.1065684 0.10341343 0.1509753
#> R-HSA-109606 Intrinsic Pathway for Apoptosis 0.1085210 0.10553758 0.1128828
#> R-HSA-109703 PKB-mediated events 0.1272747 0.05257374 0.1060544
#> Cluster_12 Cluster_13 Cluster_2 Cluster_3 Cluster_4 Cluster_5
#> R-HSA-1059683 0.10846242 0.10910438 0.11966042 0.11517448 0.11644377 0.11112119
#> R-HSA-109606 0.11381796 0.12695957 0.10476205 0.10735390 0.11049311 0.10289929
#> R-HSA-109703 0.09543865 0.07311941 0.08310135 0.08400314 0.05539675 0.04613623
#> Cluster_6 Cluster_7 Cluster_8 Cluster_9
#> R-HSA-1059683 0.09762711 0.12081193 0.13870308 0.10850241
#> R-HSA-109606 0.10827093 0.11286957 0.11733161 0.11316059
#> R-HSA-109703 0.12388218 0.07717628 0.07794358 0.01419815
A simple approach to find the most relevant pathways is to assess the maximum difference in expression for every pathway:
# find the maximum differently expressed pathway
max_difference <- do.call(rbind, apply(pathway_expression, 1, function(row) {
values <- as.numeric(row[2:length(row)])
return(data.frame(name = row[1], min = min(values), max = max(values)))
}))
max_difference$diff <- max_difference$max - max_difference$min
# sort based on the difference
max_difference <- max_difference[order(max_difference$diff, decreasing = T), ]
head(max_difference)
#> name min
#> R-HSA-350864 Regulation of thyroid hormone activity -0.4872382
#> R-HSA-8964540 Alanine metabolism -0.5062071
#> R-HSA-190374 FGFR1c and Klotho ligand binding and activation -0.3435773
#> R-HSA-140180 COX reactions -0.3451581
#> R-HSA-9024909 BDNF activates NTRK2 (TRKB) signaling -0.3743071
#> R-HSA-9025046 NTF3 activates NTRK2 (TRKB) signaling -0.3953189
#> max diff
#> R-HSA-350864 0.3744482 0.8616864
#> R-HSA-8964540 0.2553352 0.7615423
#> R-HSA-190374 0.4149894 0.7585667
#> R-HSA-140180 0.3722165 0.7173746
#> R-HSA-9024909 0.3229013 0.6972084
#> R-HSA-9025046 0.2987376 0.6940564
The ReactomeGSA package contains two functions to visualize these pathway results. The first simply plots the expression for a selected pathway:
For a better overview, the expression of multiple pathways can be shown as a heatmap using gplots
heatmap.2
function:
# Additional parameters are directly passed to gplots heatmap.2 function
plot_gsva_heatmap(gsva_result, max_pathways = 15, margins = c(6,20))
The plot_gsva_heatmap
function can also be used to only display specific pahtways:
# limit to selected B cell related pathways
relevant_pathways <- c("R-HSA-983170", "R-HSA-388841", "R-HSA-2132295", "R-HSA-983705", "R-HSA-5690714")
plot_gsva_heatmap(gsva_result,
pathway_ids = relevant_pathways, # limit to these pathways
margins = c(6,30), # adapt the figure margins in heatmap.2
dendrogram = "col", # only plot column dendrogram
scale = "row", # scale for each pathway
key = FALSE, # don't display the color key
lwid=c(0.1,4)) # remove the white space on the left
This analysis shows us that cluster 8 has a marked up-regulation of B Cell receptor signalling, which is linked to a co-stimulation of the CD28 family. Additionally, there is a gradient among the cluster with respect to genes releated to antigen presentation.
Therefore, we are able to further classify the observed B cell subtypes based on their pathway activity.
The pathway-level expression analysis can also be used to run a Principal Component Analysis on the samples. This is simplified through the function plot_gsva_pca
:
In this analysis, cluster 11 is a clear outlier from the other B cell subtypes and therefore might be prioritised for further evaluation.
sessionInfo()
#> R version 4.1.1 (2021-08-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.13-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] ReactomeGSA.data_1.6.0 SeuratObject_4.0.2 Seurat_4.0.4
#> [4] edgeR_3.34.1 limma_3.48.3 ReactomeGSA_1.6.1
#>
#> loaded via a namespace (and not attached):
#> [1] Rtsne_0.15 colorspace_2.0-2 deldir_0.2-10
#> [4] ellipsis_0.3.2 ggridges_0.5.3 spatstat.data_2.1-0
#> [7] farver_2.1.0 leiden_0.3.9 listenv_0.8.0
#> [10] ggrepel_0.9.1 fansi_0.5.0 codetools_0.2-18
#> [13] splines_4.1.1 knitr_1.34 polyclip_1.10-0
#> [16] jsonlite_1.7.2 ica_1.0-2 cluster_2.1.2
#> [19] png_0.1-7 uwot_0.1.10 shiny_1.7.0
#> [22] sctransform_0.3.2 spatstat.sparse_2.0-0 BiocManager_1.30.16
#> [25] compiler_4.1.1 httr_1.4.2 assertthat_0.2.1
#> [28] Matrix_1.3-4 fastmap_1.1.0 lazyeval_0.2.2
#> [31] later_1.3.0 prettyunits_1.1.1 htmltools_0.5.2
#> [34] tools_4.1.1 igraph_1.2.6 gtable_0.3.0
#> [37] glue_1.4.2 RANN_2.6.1 reshape2_1.4.4
#> [40] dplyr_1.0.7 Rcpp_1.0.7 scattermore_0.7
#> [43] jquerylib_0.1.4 vctrs_0.3.8 nlme_3.1-153
#> [46] lmtest_0.9-38 xfun_0.26 stringr_1.4.0
#> [49] globals_0.14.0 mime_0.11 miniUI_0.1.1.1
#> [52] lifecycle_1.0.0 irlba_2.3.3 gtools_3.9.2
#> [55] goftest_1.2-2 future_1.22.1 MASS_7.3-54
#> [58] zoo_1.8-9 scales_1.1.1 spatstat.core_2.3-0
#> [61] hms_1.1.0 promises_1.2.0.1 spatstat.utils_2.2-0
#> [64] parallel_4.1.1 RColorBrewer_1.1-2 curl_4.3.2
#> [67] yaml_2.2.1 reticulate_1.22 pbapply_1.5-0
#> [70] gridExtra_2.3 ggplot2_3.3.5 sass_0.4.0
#> [73] rpart_4.1-15 stringi_1.7.4 highr_0.9
#> [76] caTools_1.18.2 rlang_0.4.11 pkgconfig_2.0.3
#> [79] matrixStats_0.61.0 bitops_1.0-7 evaluate_0.14
#> [82] lattice_0.20-45 ROCR_1.0-11 purrr_0.3.4
#> [85] tensor_1.5 labeling_0.4.2 patchwork_1.1.1
#> [88] htmlwidgets_1.5.4 cowplot_1.1.1 tidyselect_1.1.1
#> [91] parallelly_1.28.1 RcppAnnoy_0.0.19 plyr_1.8.6
#> [94] magrittr_2.0.1 R6_2.5.1 gplots_3.1.1
#> [97] generics_0.1.0 DBI_1.1.1 mgcv_1.8-36
#> [100] pillar_1.6.2 fitdistrplus_1.1-5 survival_3.2-13
#> [103] abind_1.4-5 tibble_3.1.4 future.apply_1.8.1
#> [106] crayon_1.4.1 KernSmooth_2.23-20 utf8_1.2.2
#> [109] spatstat.geom_2.2-2 plotly_4.9.4.1 rmarkdown_2.11
#> [112] progress_1.2.2 locfit_1.5-9.4 grid_4.1.1
#> [115] data.table_1.14.0 digest_0.6.27 xtable_1.8-4
#> [118] tidyr_1.1.3 httpuv_1.6.3 munsell_0.5.0
#> [121] viridisLite_0.4.0 bslib_0.3.0