The iSEEfier User’s Guide
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")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)To configure the initial state of our iSEE instance using iSEEinit(), we
need five parameters:
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.
features : which is a vector or a dataframe containing the 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") # betareddim_type : In this example we decided to plot our data as a t-SNE plot.reddim_type <- "TSNE"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 namesgroups : Here we can add the groups/conditions/cell-types# ICR vs C57BL/6
group <- "strain" #the name should match what's in the colData namesWe can choose to include in this initial step a MarkdownBoard by setting the arguments
add_markdown_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)In case our features parameter was a data.frame, we could assign the name of the column containing the features to the gene_id parameter.
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:
ReducedDimensionPlot, FeatureAssayPlot and RowDataTable for each single gene in features.ComplexHeatmapPlot with all genes in featuresColumnDataPlot panelMarkdownBoard paneliSEEnrich()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 objectcollection: 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.reddim_type: A string vector containing the dimensionality reduction typeclusters: A character string containing the name of the clusters/cell-type/state…(as listed in the colData of the sce)groups: A character string of the groups/conditions…(as it appears in the colData of the sce)GO_collection <- "GO"
Mm_organism <- "org.Mm.eg.db"
gene_id <- "SYMBOL"
cluster <- "label"
group <- "strain"
reddim_type <- "PCA"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,
                     clusters = cluster,
                     groups = group,
                     reddim_type = reddim_type)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)iSEEmarker()In many cases, we are interested in determining the identity of our clusters, or further subset our cells types. That’s where iSEEmarker() comes in handy. Similar to iSEEinit(), we need the following parameters:
sce: a SingleCellExperiment objectclusters: the name of the clusters/cell-type/stategroups: the groups/conditionsselection_plot_format: the class of the panel that we will be using to select the clusters of interest.initial3 <- iSEEmarker(
  sce = sce,
  clusters = cluster,
  groups = group,
  selection_plot_format = "ColumnDataPlot")This function returns a list of panels, with the goal of visualizing the expression of
marker genes selected from the DynamicMarkerTable in each cell cell type. Unlike iSEEinit(), which requires us to specify a list of genes, iSEEmarker() utilizes the DynamicMarkerTable that performs statistical testing through the findMarkers() function from the scran package.
To start exploring the marker genes of each cell type with iSEE, we run:
iSEE(sce, initial = initial3)view_initial_tiles()Previously, we successfully generated three 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)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 7ac5d95 DN-- 11 3 -- 
#> + attr: name (v/c), color (v/c)
#> + edges from 7ac5d95 (vertex names):
#> [1] ReducedDimensionPlot1->ColumnDataPlot1  
#> [2] ReducedDimensionPlot2->ColumnDataPlot1  
#> [3] ReducedDimensionPlot3->FeatureAssayPlot3
initial2 <- results$initial
g2 <- view_initial_network(initial2, plot_format = "visNetwork")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)sessionInfo()
#> R version 4.5.0 RC (2025-04-04 r88126)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.2 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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.1.4               
#>  [3] scater_1.36.0               ggplot2_3.5.2              
#>  [5] scuttle_1.18.0              scRNAseq_2.21.1            
#>  [7] SingleCellExperiment_1.30.0 SummarizedExperiment_1.38.0
#>  [9] Biobase_2.68.0              GenomicRanges_1.60.0       
#> [11] GenomeInfoDb_1.44.0         IRanges_2.42.0             
#> [13] S4Vectors_0.46.0            BiocGenerics_0.54.0        
#> [15] generics_0.1.3              MatrixGenerics_1.20.0      
#> [17] matrixStats_1.5.0           iSEEfier_1.4.0             
#> [19] BiocStyle_2.36.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] splines_4.5.0            later_1.4.2              BiocIO_1.18.0           
#>   [4] bitops_1.0-9             filelock_1.0.3           tibble_3.2.1            
#>   [7] XML_3.99-0.18            lifecycle_1.0.4          httr2_1.1.2             
#>  [10] doParallel_1.0.17        lattice_0.22-7           ensembldb_2.32.0        
#>  [13] alabaster.base_1.8.0     magrittr_2.0.3           sass_0.4.10             
#>  [16] rmarkdown_2.29           jquerylib_0.1.4          yaml_2.3.10             
#>  [19] httpuv_1.6.15            DBI_1.2.3                RColorBrewer_1.1-3      
#>  [22] abind_1.4-8              Rtsne_0.17               AnnotationFilter_1.32.0 
#>  [25] RCurl_1.98-1.17          rappdirs_0.3.3           circlize_0.4.16         
#>  [28] GenomeInfoDbData_1.2.14  ggrepel_0.9.6            irlba_2.3.5.1           
#>  [31] alabaster.sce_1.8.0      iSEEhex_1.10.0           codetools_0.2-20        
#>  [34] DelayedArray_0.34.0      DT_0.33                  tidyselect_1.2.1        
#>  [37] shape_1.4.6.1            farver_2.1.2             UCSC.utils_1.4.0        
#>  [40] viridis_0.6.5            ScaledMatrix_1.16.0      shinyWidgets_0.9.0      
#>  [43] BiocFileCache_2.16.0     GenomicAlignments_1.44.0 jsonlite_2.0.0          
#>  [46] GetoptLong_1.0.5         BiocNeighbors_2.2.0      iterators_1.0.14        
#>  [49] foreach_1.5.2            tools_4.5.0              Rcpp_1.0.14             
#>  [52] glue_1.8.0               gridExtra_2.3            SparseArray_1.8.0       
#>  [55] BiocBaseUtils_1.10.0     xfun_0.52                mgcv_1.9-3              
#>  [58] dplyr_1.1.4              HDF5Array_1.36.0         gypsum_1.4.0            
#>  [61] shinydashboard_0.7.2     withr_3.0.2              BiocManager_1.30.25     
#>  [64] fastmap_1.2.0            rhdf5filters_1.20.0      shinyjs_2.1.0           
#>  [67] digest_0.6.37            rsvd_1.0.5               R6_2.6.1                
#>  [70] mime_0.13                colorspace_2.1-1         listviewer_4.0.0        
#>  [73] RSQLite_2.3.9            h5mread_1.0.0            hexbin_1.28.5           
#>  [76] FNN_1.1.4.1              rtracklayer_1.68.0       httr_1.4.7              
#>  [79] htmlwidgets_1.6.4        S4Arrays_1.8.0           org.Mm.eg.db_3.21.0     
#>  [82] uwot_0.2.3               iSEE_2.20.0              pkgconfig_2.0.3         
#>  [85] gtable_0.3.6             blob_1.2.4               ComplexHeatmap_2.24.0   
#>  [88] XVector_0.48.0           htmltools_0.5.8.1        bookdown_0.43           
#>  [91] ProtGenerics_1.40.0      rintrojs_0.3.4           clue_0.3-66             
#>  [94] scales_1.3.0             alabaster.matrix_1.8.0   png_0.1-8               
#>  [97] knitr_1.50               rjson_0.2.23             nlme_3.1-168            
#> [100] curl_6.2.2               shinyAce_0.4.4           cachem_1.1.0            
#> [103] rhdf5_2.52.0             GlobalOptions_0.1.2      BiocVersion_3.21.1      
#> [106] parallel_4.5.0           miniUI_0.1.1.1           vipor_0.4.7             
#> [109] AnnotationDbi_1.70.0     restfulr_0.0.15          pillar_1.10.2           
#> [112] grid_4.5.0               alabaster.schemas_1.8.0  vctrs_0.6.5             
#> [115] promises_1.3.2           BiocSingular_1.24.0      dbplyr_2.5.0            
#> [118] iSEEu_1.20.0             beachmat_2.24.0          xtable_1.8-4            
#> [121] cluster_2.1.8.1          beeswarm_0.4.0           evaluate_1.0.3          
#> [124] magick_2.8.6             tinytex_0.57             GenomicFeatures_1.60.0  
#> [127] cli_3.6.4                compiler_4.5.0           Rsamtools_2.24.0        
#> [130] rlang_1.1.6              crayon_1.5.3             ggbeeswarm_0.7.2        
#> [133] viridisLite_0.4.2        alabaster.se_1.8.0       BiocParallel_1.42.0     
#> [136] munsell_0.5.1            Biostrings_2.76.0        lazyeval_0.2.2          
#> [139] colourpicker_1.3.0       Matrix_1.7-3             ExperimentHub_2.16.0    
#> [142] bit64_4.6.0-1            Rhdf5lib_1.30.0          KEGGREST_1.48.0         
#> [145] shiny_1.10.0             alabaster.ranges_1.8.0   AnnotationHub_3.16.0    
#> [148] memoise_2.0.1            bslib_0.9.0              bit_4.6.0