1 Introduction

This vignette describes how to use the iSEEfier package to configure various initial states of iSEE instances, in order to simplify the task of visualizing single-cell RNA-seq, bulk RNA-seq data, or even your proteomics data in iSEE. In the remainder of this vignette, we will illustrate the main features of r BiocStyle::Biocpkg("iSEEfier") on a publicly available dataset from Baron et al. “A Single-Cell Transcriptomic Map of the Human and Mouse Pancreas Reveals Inter- and Intra-cell Population Structure”, published in Cell Systems in 2016. doi:10.1016/j.cels.2016.08.011. The data is made available via the scRNAseq Bioconductor package. We’ll simply use the mouse dataset, consisting of islets isolated from five C57BL/6 and ICR mice. # Getting started {#gettingstarted} To install iSEEfier package, we start R and enter:

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

Once installed, the package can be loaded and attached to the current workspace as follows:

library("iSEEfier")

2 Create an initial state for gene expression visualization using iSEEinit()

When we have all input elements ready, we can create an iSEE initial state by running:

iSEEinit(sce = sce_obj,
         features = feature_list,
         reddim.type = reduced_dim,
         clusters = cluster,
         groups = group,
         add_markdown_panel = FALSE,
         add_dynamicTable_panel = FALSE)

To configure the initial state of our iSEE instance using iSEEinit(), we need five parameters:

  1. sce : A SingleCellExperiment object. This object stores information of different quantifications (counts, log-expression…), dimensionality reduction coordinates (t-SNE, UMAP…), as well as some metadata related to the samples and features. We’ll start by loading the sce object:
library("scRNAseq")
sce <- BaronPancreasData('mouse')
sce
#> class: SingleCellExperiment 
#> dim: 14878 1886 
#> metadata(0):
#> assays(1): counts
#> rownames(14878): X0610007P14Rik X0610009B22Rik ... Zzz3 l7Rn6
#> rowData names(0):
#> colnames(1886): mouse1_lib1.final_cell_0001 mouse1_lib1.final_cell_0002
#>   ... mouse2_lib3.final_cell_0394 mouse2_lib3.final_cell_0395
#> colData names(2): strain label
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):

Let’s add the normalized counts

library("scuttle")
sce <- logNormCounts(sce)

Now we can add different dimensionality reduction coordinates

library("scater")
sce <- runPCA(sce)
sce <- runTSNE(sce)
sce <- runUMAP(sce)

Now our sce is ready, we can move on to the next argument.

  1. features : which is a list or a vector of genes/features of interest. Let’s say we would like to visualize the expression of some genes that were identified as marker genes for different cell population.
gene_list <- c("Gcg", # alpha
               "Ins1") # beta
  1. reddim_type : In this example we decided to plot our data as a t-SNE plot.
reddim_type <- "TSNE"
  1. clusters : Now we specify what clusters/cell-types/states/samples we would like to color/split our data with
# cell populations
cluster <- "label" #the name should match what's in the colData names
  1. groups : Here we can add the groups/conditions/cell-types
# ICR vs C57BL/6
group <- "strain" #the name should match what's in the colData names

We can choose to include in this initial step a MarkdownBoard and a DynamicMarkerTable, along with its linked panels by setting the arguments add_markdown_panel and add_dynamicTable_panel to TRUE. At this point, all the elements are ready to be transferred into iSEEinit()

initial1 <- iSEEinit(sce = sce,
                    features = gene_list,
                    clusters = cluster,
                    groups = group,
                    add_markdown_panel = TRUE,
                    add_dynamicTable_panel = TRUE)

Now we are one step away from visualizing our list of genes of interest. All that’s left to do is to run iSEE with the initial state created with iSEEinit()

library("iSEE")
iSEE(sce, initial= initial1)

This instance, generated with iSEEinit(), returns a combination of panels, linked to each other, with the goal of visualizing the expression of certain marker genes in each cell population/group:

  • A ReducedDimensionPlot, FeatureAssayPlot and RowDataTable for each single gene in features.
  • A ComplexHeatmapPlot with all genes in features
  • A DynamicMarkerTable that identifies marker genes from a sample selection.
  • A ColumnDataPlot panel
  • A MarkdownBoard panel

3 Create an initial state for feature sets exploration using iSEEnrich()

Sometimes it is interesting to look at some specific feature sets and the associated genes. That’s when the utility of iSEEnrich becomes apparent. We will need 4 elements to explore feature sets of interest:

  • sce: A SingleCellExperiment object
  • collection: A character vector specifying the gene set collections of interest (it is possible to use GO or KEGG terms)
  • gene_identifier: A character string specifying the identifier to use to extract gene IDs for the organism package. This can be “ENS” for ENSEMBL ids, “SYMBOL” for gene names…
  • organism: A character string of the org.*.eg.db package to use to extract mappings of gene sets to gene IDs.
GO_collection <- "GO"
Mm_organism <- "org.Mm.eg.db"
gene_id <- "SYMBOL"

Now let’s create this initial setup for iSEE using iSEEnrich()

results <- iSEEnrich(sce = sce,
                     collection = GO_collection,
                     organism = Mm_organism,
                     gene_identifier = gene_id)

iSEEnrich will specifically return a list with the updated sce object and its associated initial configuration. To start the iSEE instance we run:

iSEE(results$sce, initial = results$initial)

4 Visualize a preview of the initial configurations with view_initial_tiles()

Previously, we successfully generated two distinct initial configurations for iSEE. However, understanding the expected content of our iSEE instances is not always straightforward. That’s when we can use view_initial_tiles(). We only need as an input the initial configuration to obtain a graphical visualization of the expected the corresponding iSEE instance:

library(ggplot2)
view_initial_tiles(initial = initial1)

view_initial_tiles(initial = results$initial)

5 Visualize network connections between panels with view_initial_network()

As some of these panels are linked to each other, we can visualize these networks with view_initial_network(). Similar to iSEEconfigviewer(), this function takes the initial setup as input: This function always returns the igraph object underlying the visualizations that can be displayed as a side effect.

library("igraph")
library("visNetwork")
g1 <- view_initial_network(initial1, plot_format = "igraph")

g1
#> IGRAPH 04925e3 DN-- 14 5 -- 
#> + attr: name (v/c), color (v/c)
#> + edges from 04925e3 (vertex names):
#> [1] ReducedDimensionPlot1->ColumnDataPlot1      
#> [2] ReducedDimensionPlot2->ColumnDataPlot1      
#> [3] ReducedDimensionPlot3->FeatureAssayPlot3    
#> [4] ReducedDimensionPlot4->FeatureAssayPlot4    
#> [5] DynamicMarkerTable1  ->ReducedDimensionPlot4
initial2 <- results$initial
g2 <- view_initial_network(initial2, plot_format = "visNetwork")

6 Merge different initial configurations with glue_initials()

Sometimes, it would be interesting to merge different iSEE initial configurations to visualize all different panel in the same iSEE instance.

merged_config <- glue_initials(initial1,initial2)

We can then preview the content of this initial configuration

view_initial_tiles(merged_config)

Session info

sessionInfo()
#> R version 4.4.0 RC (2024-04-16 r86468)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.20-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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] visNetwork_2.1.2            igraph_2.0.3               
#>  [3] scater_1.33.0               ggplot2_3.5.1              
#>  [5] scuttle_1.15.0              scRNAseq_2.17.10           
#>  [7] SingleCellExperiment_1.27.0 SummarizedExperiment_1.35.0
#>  [9] Biobase_2.65.0              GenomicRanges_1.57.0       
#> [11] GenomeInfoDb_1.41.0         IRanges_2.39.0             
#> [13] S4Vectors_0.43.0            BiocGenerics_0.51.0        
#> [15] MatrixGenerics_1.17.0       matrixStats_1.3.0          
#> [17] iSEEfier_1.1.0              BiocStyle_2.33.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] splines_4.4.0             later_1.3.2              
#>   [3] BiocIO_1.15.0             bitops_1.0-7             
#>   [5] filelock_1.0.3            tibble_3.2.1             
#>   [7] XML_3.99-0.16.1           lifecycle_1.0.4          
#>   [9] httr2_1.0.1               doParallel_1.0.17        
#>  [11] lattice_0.22-6            ensembldb_2.29.0         
#>  [13] alabaster.base_1.5.0      magrittr_2.0.3           
#>  [15] sass_0.4.9                rmarkdown_2.26           
#>  [17] jquerylib_0.1.4           yaml_2.3.8               
#>  [19] httpuv_1.6.15             DBI_1.2.2                
#>  [21] RColorBrewer_1.1-3        abind_1.4-5              
#>  [23] zlibbioc_1.51.0           Rtsne_0.17               
#>  [25] AnnotationFilter_1.29.0   RCurl_1.98-1.14          
#>  [27] rappdirs_0.3.3            circlize_0.4.16          
#>  [29] GenomeInfoDbData_1.2.12   ggrepel_0.9.5            
#>  [31] irlba_2.3.5.1             alabaster.sce_1.5.0      
#>  [33] iSEEhex_1.7.0             DelayedMatrixStats_1.27.0
#>  [35] codetools_0.2-20          DelayedArray_0.31.0      
#>  [37] DT_0.33                   tidyselect_1.2.1         
#>  [39] shape_1.4.6.1             farver_2.1.1             
#>  [41] UCSC.utils_1.1.0          viridis_0.6.5            
#>  [43] ScaledMatrix_1.13.0       shinyWidgets_0.8.6       
#>  [45] BiocFileCache_2.13.0      GenomicAlignments_1.41.0 
#>  [47] jsonlite_1.8.8            BiocNeighbors_1.23.0     
#>  [49] GetoptLong_1.0.5          iterators_1.0.14         
#>  [51] foreach_1.5.2             tools_4.4.0              
#>  [53] Rcpp_1.0.12               glue_1.7.0               
#>  [55] gridExtra_2.3             SparseArray_1.5.0        
#>  [57] BiocBaseUtils_1.7.0       xfun_0.43                
#>  [59] mgcv_1.9-1                dplyr_1.1.4              
#>  [61] HDF5Array_1.33.0          gypsum_1.1.0             
#>  [63] shinydashboard_0.7.2      withr_3.0.0              
#>  [65] BiocManager_1.30.22       fastmap_1.1.1            
#>  [67] rhdf5filters_1.17.0       fansi_1.0.6              
#>  [69] shinyjs_2.1.0             rsvd_1.0.5               
#>  [71] digest_0.6.35             R6_2.5.1                 
#>  [73] mime_0.12                 colorspace_2.1-0         
#>  [75] listviewer_4.0.0          RSQLite_2.3.6            
#>  [77] paws.storage_0.5.0        utf8_1.2.4               
#>  [79] generics_0.1.3            hexbin_1.28.3            
#>  [81] FNN_1.1.4                 rtracklayer_1.65.0       
#>  [83] httr_1.4.7                htmlwidgets_1.6.4        
#>  [85] S4Arrays_1.5.0            org.Mm.eg.db_3.19.1      
#>  [87] uwot_0.2.2                iSEE_2.17.0              
#>  [89] pkgconfig_2.0.3           gtable_0.3.5             
#>  [91] blob_1.2.4                ComplexHeatmap_2.21.0    
#>  [93] XVector_0.45.0            htmltools_0.5.8.1        
#>  [95] bookdown_0.39             ProtGenerics_1.37.0      
#>  [97] rintrojs_0.3.4            clue_0.3-65              
#>  [99] scales_1.3.0              alabaster.matrix_1.5.0   
#> [101] png_0.1-8                 knitr_1.46               
#> [103] rjson_0.2.21              nlme_3.1-164             
#> [105] curl_5.2.1                shinyAce_0.4.2           
#> [107] cachem_1.0.8              rhdf5_2.49.0             
#> [109] GlobalOptions_0.1.2       BiocVersion_3.20.0       
#> [111] parallel_4.4.0            miniUI_0.1.1.1           
#> [113] vipor_0.4.7               AnnotationDbi_1.67.0     
#> [115] restfulr_0.0.15           pillar_1.9.0             
#> [117] grid_4.4.0                alabaster.schemas_1.5.0  
#> [119] vctrs_0.6.5               promises_1.3.0           
#> [121] BiocSingular_1.21.0       dbplyr_2.5.0             
#> [123] iSEEu_1.17.0              beachmat_2.21.0          
#> [125] xtable_1.8-4              cluster_2.1.6            
#> [127] beeswarm_0.4.0            evaluate_0.23            
#> [129] magick_2.8.3              tinytex_0.50             
#> [131] GenomicFeatures_1.57.0    cli_3.6.2                
#> [133] compiler_4.4.0            Rsamtools_2.21.0         
#> [135] rlang_1.1.3               crayon_1.5.2             
#> [137] paws.common_0.7.2         ggbeeswarm_0.7.2         
#> [139] viridisLite_0.4.2         alabaster.se_1.5.0       
#> [141] BiocParallel_1.39.0       munsell_0.5.1            
#> [143] Biostrings_2.73.0         lazyeval_0.2.2           
#> [145] colourpicker_1.3.0        Matrix_1.7-0             
#> [147] ExperimentHub_2.13.0      sparseMatrixStats_1.17.0 
#> [149] bit64_4.0.5               Rhdf5lib_1.27.0          
#> [151] KEGGREST_1.45.0           shiny_1.8.1.1            
#> [153] highr_0.10                alabaster.ranges_1.5.0   
#> [155] AnnotationHub_3.13.0      memoise_2.0.1            
#> [157] bslib_0.7.0               bit_4.0.5