Analysing single-cell RNA-sequencing Data

Johannes Griss

2022-11-01

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")
#> 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
#> Attaching sp

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

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
#> 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: 82
#>   Results:
#>   - Seurat:
#>     1748 pathways
#>     11011 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.1083024 0.09933853  0.1419717
#> R-HSA-109606  Intrinsic Pathway for Apoptosis 0.1120826 0.10765166  0.1129437
#> R-HSA-109703              PKB-mediated events 0.1269551 0.05260470  0.1062727
#>               Cluster_12 Cluster_13  Cluster_2  Cluster_3  Cluster_4  Cluster_5
#> R-HSA-1059683 0.11064768 0.10186327 0.11589537 0.11294177 0.11265002 0.10465541
#> R-HSA-109606  0.11687392 0.12568517 0.10688977 0.10872405 0.11113868 0.10329639
#> R-HSA-109703  0.09538291 0.07349972 0.08298411 0.08385108 0.05546712 0.04589538
#>                Cluster_6  Cluster_7  Cluster_8  Cluster_9
#> R-HSA-1059683 0.09748855 0.12269029 0.13386683 0.10123911
#> R-HSA-109606  0.10807664 0.11563910 0.11937031 0.11260225
#> R-HSA-109703  0.12357026 0.07697748 0.07796547 0.01425426

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.4876458
#> R-HSA-8964540                              Alanine metabolism -0.5061431
#> R-HSA-190374  FGFR1c and Klotho ligand binding and activation -0.3434411
#> R-HSA-140180                                    COX reactions -0.3450976
#> R-HSA-9024909           BDNF activates NTRK2 (TRKB) signaling -0.3746819
#> R-HSA-5263617       Metabolism of ingested MeSeO2H into MeSeH -0.1941002
#>                     max      diff
#> R-HSA-350864  0.3752821 0.8629279
#> R-HSA-8964540 0.2551881 0.7613312
#> R-HSA-190374  0.4157735 0.7592146
#> R-HSA-140180  0.3721532 0.7172508
#> R-HSA-9024909 0.3232355 0.6979174
#> R-HSA-5263617 0.4938569 0.6879571

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 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] ReactomeGSA.data_1.11.0 sp_1.5-0                SeuratObject_4.1.2     
#> [4] Seurat_4.2.0            edgeR_3.40.0            limma_3.54.0           
#> [7] ReactomeGSA_1.12.0     
#> 
#> loaded via a namespace (and not attached):
#>   [1] Rtsne_0.16            colorspace_2.0-3      deldir_1.0-6         
#>   [4] ellipsis_0.3.2        ggridges_0.5.4        spatstat.data_3.0-0  
#>   [7] farver_2.1.1          leiden_0.4.3          listenv_0.8.0        
#>  [10] ggrepel_0.9.1         fansi_1.0.3           codetools_0.2-18     
#>  [13] splines_4.2.1         cachem_1.0.6          knitr_1.40           
#>  [16] polyclip_1.10-4       jsonlite_1.8.3        ica_1.0-3            
#>  [19] cluster_2.1.4         png_0.1-7             rgeos_0.5-9          
#>  [22] uwot_0.1.14           shiny_1.7.3           sctransform_0.3.5    
#>  [25] spatstat.sparse_3.0-0 BiocManager_1.30.19   compiler_4.2.1       
#>  [28] httr_1.4.4            assertthat_0.2.1      Matrix_1.5-1         
#>  [31] fastmap_1.1.0         lazyeval_0.2.2        cli_3.4.1            
#>  [34] later_1.3.0           prettyunits_1.1.1     htmltools_0.5.3      
#>  [37] tools_4.2.1           igraph_1.3.5          gtable_0.3.1         
#>  [40] glue_1.6.2            RANN_2.6.1            reshape2_1.4.4       
#>  [43] dplyr_1.0.10          Rcpp_1.0.9            scattermore_0.8      
#>  [46] jquerylib_0.1.4       vctrs_0.5.0           nlme_3.1-160         
#>  [49] progressr_0.11.0      lmtest_0.9-40         spatstat.random_3.0-0
#>  [52] xfun_0.34             stringr_1.4.1         globals_0.16.1       
#>  [55] mime_0.12             miniUI_0.1.1.1        lifecycle_1.0.3      
#>  [58] irlba_2.3.5.1         gtools_3.9.3          goftest_1.2-3        
#>  [61] future_1.28.0         MASS_7.3-58.1         zoo_1.8-11           
#>  [64] scales_1.2.1          spatstat.core_2.4-4   hms_1.1.2            
#>  [67] promises_1.2.0.1      spatstat.utils_3.0-1  parallel_4.2.1       
#>  [70] RColorBrewer_1.1-3    curl_4.3.3            yaml_2.3.6           
#>  [73] reticulate_1.26       pbapply_1.5-0         gridExtra_2.3        
#>  [76] ggplot2_3.3.6         sass_0.4.2            rpart_4.1.19         
#>  [79] stringi_1.7.8         highr_0.9             caTools_1.18.2       
#>  [82] rlang_1.0.6           pkgconfig_2.0.3       matrixStats_0.62.0   
#>  [85] bitops_1.0-7          evaluate_0.17         lattice_0.20-45      
#>  [88] tensor_1.5            ROCR_1.0-11           purrr_0.3.5          
#>  [91] labeling_0.4.2        patchwork_1.1.2       htmlwidgets_1.5.4    
#>  [94] cowplot_1.1.1         tidyselect_1.2.0      parallelly_1.32.1    
#>  [97] RcppAnnoy_0.0.20      plyr_1.8.7            magrittr_2.0.3       
#> [100] R6_2.5.1              gplots_3.1.3          generics_0.1.3       
#> [103] DBI_1.1.3             mgcv_1.8-41           pillar_1.8.1         
#> [106] fitdistrplus_1.1-8    abind_1.4-5           survival_3.4-0       
#> [109] tibble_3.1.8          future.apply_1.9.1    crayon_1.5.2         
#> [112] KernSmooth_2.23-20    utf8_1.2.2            spatstat.geom_3.0-3  
#> [115] plotly_4.10.0         rmarkdown_2.17        progress_1.2.2       
#> [118] locfit_1.5-9.6        grid_4.2.1            data.table_1.14.4    
#> [121] digest_0.6.30         xtable_1.8-4          tidyr_1.2.1          
#> [124] httpuv_1.6.6          munsell_0.5.0         viridisLite_0.4.1    
#> [127] bslib_0.4.0