Analysing single-cell RNA-sequencing Data

Johannes Griss

2024-02-02

Introduction

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.

Citation

To cite this package, use

Griss J. ReactomeGSA, https://github.com/reactome/ReactomeGSA (2019)

Installation

The ReactomeGSA package can be directly installed from Bioconductor:

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")

if (!require(ReactomeGSA))
  BiocManager::install("ReactomeGSA")

# install the ReactomeGSA.data package for the example data
if (!require(ReactomeGSA.data))
  BiocManager::install("ReactomeGSA.data")

For more information, see https://bioconductor.org/install/.

Example data

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.

library(ReactomeGSA.data)
#> Loading required package: limma
#> Loading required package: edgeR
#> Loading required package: ReactomeGSA
#> Loading required package: Seurat
#> Loading required package: SeuratObject
#> Loading required package: sp
#> 
#> Attaching package: 'SeuratObject'
#> The following object is masked from 'package:base':
#> 
#>     intersect
data(jerby_b_cells)

jerby_b_cells
#> An object of class Seurat 
#> 23686 features across 920 samples within 1 assay 
#> Active assay: RNA (23686 features, 0 variable features)
#>  2 layers present: counts, data

Pathway analysis of cell clusters

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: 86
#>   Results:
#>   - Seurat:
#>     1773 pathways
#>     11122 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.1092170 0.09778349  0.1470607
#> R-HSA-109606  Intrinsic Pathway for Apoptosis 0.1175318 0.11313998  0.1158441
#> R-HSA-109703              PKB-mediated events 0.1273845 0.05246229  0.1064824
#>               Cluster_12 Cluster_13  Cluster_2  Cluster_3  Cluster_4  Cluster_5
#> R-HSA-1059683 0.10983399 0.10497527 0.11718678 0.11532617 0.11482856 0.10741936
#> R-HSA-109606  0.11991881 0.13197679 0.11289277 0.11524909 0.11788849 0.10862934
#> R-HSA-109703  0.09546617 0.07332126 0.08331065 0.08430835 0.05562504 0.04623421
#>                Cluster_6  Cluster_7  Cluster_8  Cluster_9
#> R-HSA-1059683 0.09962466 0.12065203 0.13868658 0.10552683
#> R-HSA-109606  0.11353425 0.12080057 0.12515946 0.11837913
#> R-HSA-109703  0.12392905 0.07718021 0.07833784 0.01423891

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.4875863
#> R-HSA-8964540                              Alanine metabolism -0.5061732
#> R-HSA-190374  FGFR1c and Klotho ligand binding and activation -0.3432648
#> R-HSA-140180                                    COX reactions -0.3450870
#> R-HSA-9024909           BDNF activates NTRK2 (TRKB) signaling -0.3755033
#> R-HSA-9025046           NTF3 activates NTRK2 (TRKB) signaling -0.3963853
#>                     max      diff
#> R-HSA-350864  0.3756096 0.8631959
#> R-HSA-8964540 0.2560705 0.7622437
#> R-HSA-190374  0.4160633 0.7593281
#> R-HSA-140180  0.3723295 0.7174166
#> R-HSA-9024909 0.3235067 0.6990100
#> R-HSA-9025046 0.2993049 0.6956902

Plotting the results

The ReactomeGSA package contains two functions to visualize these pathway results. The first simply plots the expression for a selected pathway:

plot_gsva_pathway(gsva_result, pathway_id = rownames(max_difference)[1])

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.

Pathway-level PCA

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:

plot_gsva_pca(gsva_result)

In this analysis, cluster 11 is a clear outlier from the other B cell subtypes and therefore might be prioritised for further evaluation.

Session Info

sessionInfo()
#> R Under development (unstable) (2024-01-16 r85808)
#> 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] ReactomeGSA.data_1.17.1 Seurat_5.0.1            SeuratObject_5.0.1     
#> [4] sp_2.1-3                ReactomeGSA_1.17.3      edgeR_4.1.15           
#> [7] limma_3.59.1           
#> 
#> loaded via a namespace (and not attached):
#>   [1] RColorBrewer_1.1-3     jsonlite_1.8.8         magrittr_2.0.3        
#>   [4] spatstat.utils_3.0-4   farver_2.1.1           rmarkdown_2.25        
#>   [7] vctrs_0.6.5            ROCR_1.0-11            spatstat.explore_3.2-6
#>  [10] progress_1.2.3         htmltools_0.5.7        curl_5.2.0            
#>  [13] sass_0.4.8             sctransform_0.4.1      parallelly_1.36.0     
#>  [16] KernSmooth_2.23-22     bslib_0.6.1            htmlwidgets_1.6.4     
#>  [19] ica_1.0-3              plyr_1.8.9             plotly_4.10.4         
#>  [22] zoo_1.8-12             cachem_1.0.8           igraph_2.0.1.1        
#>  [25] mime_0.12              lifecycle_1.0.4        pkgconfig_2.0.3       
#>  [28] Matrix_1.6-5           R6_2.5.1               fastmap_1.1.1         
#>  [31] fitdistrplus_1.1-11    future_1.33.1          shiny_1.8.0           
#>  [34] digest_0.6.34          colorspace_2.1-0       patchwork_1.2.0       
#>  [37] tensor_1.5             RSpectra_0.16-1        irlba_2.3.5.1         
#>  [40] labeling_0.4.3         progressr_0.14.0       fansi_1.0.6           
#>  [43] spatstat.sparse_3.0-3  httr_1.4.7             polyclip_1.10-6       
#>  [46] abind_1.4-5            compiler_4.4.0         withr_3.0.0           
#>  [49] fastDummies_1.7.3      highr_0.10             gplots_3.1.3.1        
#>  [52] MASS_7.3-60.2          gtools_3.9.5           caTools_1.18.2        
#>  [55] tools_4.4.0            lmtest_0.9-40          httpuv_1.6.14         
#>  [58] future.apply_1.11.1    goftest_1.2-3          glue_1.7.0            
#>  [61] nlme_3.1-164           promises_1.2.1         grid_4.4.0            
#>  [64] Rtsne_0.17             cluster_2.1.6          reshape2_1.4.4        
#>  [67] generics_0.1.3         gtable_0.3.4           spatstat.data_3.0-4   
#>  [70] tidyr_1.3.1            hms_1.1.3              data.table_1.15.0     
#>  [73] utf8_1.2.4             BiocGenerics_0.49.1    spatstat.geom_3.2-8   
#>  [76] RcppAnnoy_0.0.22       ggrepel_0.9.5          RANN_2.6.1            
#>  [79] pillar_1.9.0           stringr_1.5.1          spam_2.10-0           
#>  [82] RcppHNSW_0.5.0         later_1.3.2            splines_4.4.0         
#>  [85] dplyr_1.1.4            lattice_0.22-5         survival_3.5-7        
#>  [88] deldir_2.0-2           tidyselect_1.2.0       locfit_1.5-9.8        
#>  [91] miniUI_0.1.1.1         pbapply_1.7-2          knitr_1.45            
#>  [94] gridExtra_2.3          scattermore_1.2        xfun_0.41             
#>  [97] Biobase_2.63.0         statmod_1.5.0          matrixStats_1.2.0     
#> [100] stringi_1.8.3          lazyeval_0.2.2         yaml_2.3.8            
#> [103] evaluate_0.23          codetools_0.2-19       tibble_3.2.1          
#> [106] cli_3.6.2              uwot_0.1.16            xtable_1.8-4          
#> [109] reticulate_1.35.0      munsell_0.5.0          jquerylib_0.1.4       
#> [112] Rcpp_1.0.12            globals_0.16.2         spatstat.random_3.2-2 
#> [115] png_0.1-8              parallel_4.4.0         ellipsis_0.3.2        
#> [118] ggplot2_3.4.4          prettyunits_1.2.0      dotCall64_1.1-1       
#> [121] bitops_1.0-7           listenv_0.9.1          viridisLite_0.4.2     
#> [124] scales_1.3.0           ggridges_0.5.6         crayon_1.5.2          
#> [127] leiden_0.4.3.1         purrr_1.0.2            rlang_1.1.3           
#> [130] cowplot_1.1.3