iSEEu 1.6.0
The iSEE package (Rue-Albrecht et al. 2018) provides a general and flexible framework for interactively exploring SummarizedExperiment objects.
However, in many cases, more specialized panels are required for effective visualization of specific data types.
The iSEEu package implements a collection of such dedicated panel classes that work directly in the iSEE application and can smoothly interact with other panels.
This allows users to quickly parametrize bespoke apps for their data to address scientific questions of interest.
We first load in the package:
library(iSEEu)All the panels described in this document can be deployed by simply passing them into the iSEE() function via the initial= argument, as shown in the following examples.
To demonstrate the use of these panels,
we will perform a differential expression analysis on the airway dataset with the edgeR package.
We store the resulting statistics in the rowData of the SummarizedExperiment so that it can be accessed by iSEE panels.
library(airway)
data(airway)
library(edgeR)
y <- DGEList(assay(airway), samples=colData(airway))
y <- y[filterByExpr(y, group=y$samples$dex),]
y <- calcNormFactors(y)
design <- model.matrix(~dex, y$samples)
y <- estimateDisp(y, design)
fit <- glmQLFit(y, design)
res <- glmQLFTest(fit, coef=2)
tab <- topTags(res, n=Inf)$table
rowData(airway) <- cbind(rowData(airway), tab[rownames(airway),])The MAPlot class creates a MA plot, i.e., with the log-fold change on the y-axis and the average expression on the x-axis.
Features with significant differences in each direction are highlighted and counted on the legend.
Users can vary the significance threshold and apply ad hoc filters on the log-fold change.
This is a subclass of the RowDataPlot so points can be transmitted to other panels as multiple row selections.
Instances of this class are created like:
ma.panel <- MAPlot(PanelWidth=6L)
app <- iSEE(airway, initial=list(ma.panel))The VolcanoPlot class creates a volcano plot with the log-fold change on the x-axis and the negative log-p-value on the y-axis.
Features with significant differences in each direction are highlighted and counted on the legend.
Users can vary the significance threshold and apply ad hoc filters on the log-fold change.
This is a subclass of the RowDataPlot so points can be transmitted to other panels as multiple row selections.
Instances of this class are created like:
vol.panel <- VolcanoPlot(PanelWidth=6L)
app <- iSEE(airway, initial=list(vol.panel))The LogFCLogFCPlot class creates a scatter plot of two log-fold changes from different DE comparisons.
This allows us to compare DE results on the same dataset - or even from different datasets, as long as the row names are shared.
Users can vary the significant threshold used to identify DE genes in either or both comparisons.
This is a subclass of the RowDataPlot so points can be transmitted to other panels as multiple row selections.
Instances of this class are created like:
# Creating another comparison, this time by blocking on the cell line
design.alt <- model.matrix(~cell + dex, y$samples)
y.alt <- estimateDisp(y, design.alt)
fit.alt <- glmQLFit(y.alt, design.alt)
res.alt <- glmQLFTest(fit.alt, coef=2)
tab.alt <- topTags(res.alt, n=Inf)$table
rowData(airway) <- cbind(rowData(airway), alt=tab.alt[rownames(airway),])
lfc.panel <- LogFCLogFCPlot(PanelWidth=6L, YAxis="alt.logFC", 
    YPValueField="alt.PValue")
app <- iSEE(airway, initial=list(lfc.panel))To demonstrate, we will perform a quick analysis of a small dataset from the scRNAseq package. This involves computing normalized expression values and low-dimensional results using the scater package.
library(scRNAseq)
sce <- ReprocessedAllenData(assays="tophat_counts")
library(scater)
sce <- logNormCounts(sce, exprs_values="tophat_counts")
sce <- runPCA(sce, ncomponents=4)
sce <- runTSNE(sce)The DynamicReducedDimensionPlot class creates a scatter plot with a dimensionality reduction result, namely principal components analysis (PCA), \(t\)-stochastic neighbor embedding (\(t\)-SNE) or uniform manifold and approximate projection (UMAP).
It does so dynamically on the subset of points that are selected in a transmitting panel,
allowing users to focus on finer structure when dealing with a heterogeneous population.
Calculations are performed using relevant functions from the scater package.
# Receives a selection from a reduced dimension plot.
dyn.panel <- DynamicReducedDimensionPlot(Type="UMAP", Assay="logcounts",
    ColumnSelectionSource="ReducedDimensionPlot1", PanelWidth=6L)
# NOTE: users do not have to manually create this, just 
# copy it from the "Panel Settings" of an already open app.
red.panel <- ReducedDimensionPlot(PanelId=1L, PanelWidth=6L,
    BrushData = list(
        xmin = -45.943, xmax = -15.399, ymin = -58.560, 
        ymax = 49.701, coords_css = list(xmin = 51.009, 
            xmax = 165.009, ymin = 39.009, 
            ymax = 422.009), coords_img = list(xmin = 66.313, 
            xmax = 214.514, ymin = 50.712, 
            ymax = 548.612), img_css_ratio = list(x = 1.300, 
            y = 1.299), mapping = list(x = "X", y = "Y"), 
        domain = list(left = -49.101, right = 57.228, 
            bottom = -70.389, top = 53.519), 
        range = list(left = 50.986, right = 566.922, 
            bottom = 603.013, top = 33.155), 
        log = list(x = NULL, y = NULL), direction = "xy", 
        brushId = "ReducedDimensionPlot1_Brush", 
        outputId = "ReducedDimensionPlot1"
    )
)
app <- iSEE(sce, initial=list(red.panel, dyn.panel))The DynamicMarkerTable class dynamically computes basic differential statistics comparing assay values across groups of multiple selections in a transmitting panel.
If only the active selection exists in the transmitting panel, a comparison is performed between the points in that selection and all unselected points.
If saved selections are present, pairwise comparisons between the active selection and each saved selection is performed and the results are combined into a single table using the findMarkers() function from scran.
diff.panel <- DynamicMarkerTable(PanelWidth=8L, Assay="logcounts",
    ColumnSelectionSource="ReducedDimensionPlot1",)
# Recycling the reduced dimension panel above, adding a saved selection to
# compare to the active selection.
red.panel[["SelectionHistory"]] <- list(
    BrushData = list(
        xmin = 15.143, xmax = 57.228, ymin = -40.752, 
        ymax = 25.674, coords_css = list(xmin = 279.009, 
            xmax = 436.089, ymin = 124.009, 
            ymax = 359.009), coords_img = list(xmin = 362.716, 
            xmax = 566.922, ymin = 161.212, 
            ymax = 466.712), img_css_ratio = list(x = 1.300, 
            y = 1.299), mapping = list(x = "X", y = "Y"), 
        domain = list(left = -49.101, right = 57.228, 
            bottom = -70.389, top = 53.519), 
        range = list(left = 50.986, right = 566.922, 
            bottom = 603.013, top = 33.155), 
        log = list(x = NULL, y = NULL), direction = "xy", 
        brushId = "ReducedDimensionPlot1_Brush", 
        outputId = "ReducedDimensionPlot1"
    )
)
red.panel[["PanelWidth"]] <- 4L # To fit onto one line.
app <- iSEE(sce, initial=list(red.panel, diff.panel))The FeatureSetTable() class is a bit unusual in that its rows do not correspond to any dimension of the SummarizedExperiment.
Rather, each row is a feature set (e.g., from GO or KEGG) that, upon click, transmits a multiple row selection to other panels.
The multiple selection consists of all rows in the chosen feature set,
allowing users to identify the positions of all genes in a pathway of interest on, say, a volcano plot.
This is also a rare example of a panel that only transmits and does not receive any selections from other panels.
setFeatureSetCommands(createGeneSetCommands(identifier="ENSEMBL"))
gset.tab <- FeatureSetTable(Selected="GO:0002576", 
    Search="platelet", PanelWidth=6L)
# This volcano plot will highlight the genes in the selected gene set.
vol.panel <- VolcanoPlot(RowSelectionSource="FeatureSetTable1",
    ColorBy="Row selection", PanelWidth=6L)
app <- iSEE(airway, initial=list(gset.tab, vol.panel))iSEEu contains a number of “modes” that allow users to conveniently load an iSEE instance in one of several common configurations:
modeEmpty() will launch an empty app, i.e., with no panels.
This is occasionally useful to jump to the landing page where a user can then upload a SummarizedExperiment object.modeGating() will launch an app with multiple feature assay panels that are linked to each other.
This is useful for applying sequential restrictions on the data, equivalent to gating in a flow cytometry experiment.modeReducedDim() will launch an app with multiple reduced dimension plots.
This is useful for examining different views of large high-dimensional datasets (e.g., single-cell studies).If you want to contribute to the development of the iSEEu package, here is a quick step-by-step guide:
iSEEu repository from GitHub (https://github.com/iSEE/iSEEu) and clone it locally.git clone https://github.com/[your_github_username]/iSEEu.gitAdd the desired new files - start from the R folder, then document via roxygen2 - and push to your fork.
As an example you can check out to understand how things are supposed to work, there are several modes already defined in the R/ directory.
A typical contribution could include e.g. a function defining an iSEE mode, named modeXXX, where XXX provides a clear representation of the purpose of the mode.
Please place each mode in a file of its own, with the same name as the function.
The function should be documented (including an example), and any required packages should be added to the DESCRIPTION file.
Once your mode/function is done, consider adding some information in the package.
Some examples might be a screenshot of the mode in action (to be placed in the folder inst/modes_img), and a well-documented example use case (maybe an entry in the vignettes folder). Also add yourself as a contributor (ctb) to the DESCRIPTION file.
Make a pull request to the original repo - the GitHub site offers a practical framework to do so, enabling comments, code reviews, and other goodies. The iSEE core team will evaluate the contribution and get back to you!
That’s pretty much it!
Example data sets can often be obtained from an ExperimentHub package (e.g. from the scRNAseq package for single-cell RNA-sequencing data), and should not be added to the iSEEu package.
testthat frameworkWe do follow some guidelines regarding the names given to variables, please abide to these for consistency with the rest of the codebase. Here are a few pointers:
git diff operations easier to check.camelCase for modes and other functions.function_name for internalsPanelClassName for panels.genericFunction for the API.scope1.scope2.name for variable names in the cached infoIf you intend to understand more in depth the internals of the iSEE framework, consider checking out the bookdown resource we put together at https://isee.github.io/iSEE-book/
Many of the “global” variables that are used in several places in iSEE are defined in the constants.R script in iSEE. We suggest to refer to those constants by their actual value rather than their internal variable name in downstream panel code. Both constant variable names and values may change at any time, but we will only announce changes to the constant value.
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.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-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] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] scater_1.22.0               ggplot2_3.3.5              
##  [3] scuttle_1.4.0               scRNAseq_2.7.2             
##  [5] edgeR_3.36.0                limma_3.50.0               
##  [7] airway_1.13.0               iSEEu_1.6.0                
##  [9] iSEE_2.6.0                  SingleCellExperiment_1.16.0
## [11] SummarizedExperiment_1.24.0 Biobase_2.54.0             
## [13] GenomicRanges_1.46.0        GenomeInfoDb_1.30.0        
## [15] IRanges_2.28.0              S4Vectors_0.32.0           
## [17] BiocGenerics_0.40.0         MatrixGenerics_1.6.0       
## [19] matrixStats_0.61.0          BiocStyle_2.22.0           
## 
## loaded via a namespace (and not attached):
##   [1] circlize_0.4.13               AnnotationHub_3.2.0          
##   [3] BiocFileCache_2.2.0           igraph_1.2.7                 
##   [5] lazyeval_0.2.2                shinydashboard_0.7.2         
##   [7] splines_4.1.1                 BiocParallel_1.28.0          
##   [9] digest_0.6.28                 ensembldb_2.18.0             
##  [11] foreach_1.5.1                 htmltools_0.5.2              
##  [13] viridis_0.6.2                 fansi_0.5.0                  
##  [15] magrittr_2.0.1                memoise_2.0.0                
##  [17] ScaledMatrix_1.2.0            cluster_2.1.2                
##  [19] doParallel_1.0.16             ComplexHeatmap_2.10.0        
##  [21] Biostrings_2.62.0             prettyunits_1.1.1            
##  [23] colorspace_2.0-2              blob_1.2.2                   
##  [25] rappdirs_0.3.3                ggrepel_0.9.1                
##  [27] xfun_0.27                     dplyr_1.0.7                  
##  [29] crayon_1.4.1                  RCurl_1.98-1.5               
##  [31] jsonlite_1.7.2                iterators_1.0.13             
##  [33] glue_1.4.2                    gtable_0.3.0                 
##  [35] zlibbioc_1.40.0               XVector_0.34.0               
##  [37] GetoptLong_1.0.5              DelayedArray_0.20.0          
##  [39] BiocSingular_1.10.0           shape_1.4.6                  
##  [41] scales_1.1.1                  DBI_1.1.1                    
##  [43] miniUI_0.1.1.1                Rcpp_1.0.7                   
##  [45] progress_1.2.2                viridisLite_0.4.0            
##  [47] xtable_1.8-4                  clue_0.3-60                  
##  [49] rsvd_1.0.5                    bit_4.0.4                    
##  [51] DT_0.19                       htmlwidgets_1.5.4            
##  [53] httr_1.4.2                    RColorBrewer_1.1-2           
##  [55] shinyAce_0.4.1                ellipsis_0.3.2               
##  [57] XML_3.99-0.8                  pkgconfig_2.0.3              
##  [59] sass_0.4.0                    dbplyr_2.1.1                 
##  [61] locfit_1.5-9.4                utf8_1.2.2                   
##  [63] tidyselect_1.1.1              rlang_0.4.12                 
##  [65] later_1.3.0                   AnnotationDbi_1.56.0         
##  [67] munsell_0.5.0                 BiocVersion_3.14.0           
##  [69] tools_4.1.1                   cachem_1.0.6                 
##  [71] generics_0.1.1                RSQLite_2.2.8                
##  [73] ExperimentHub_2.2.0           rintrojs_0.3.0               
##  [75] evaluate_0.14                 stringr_1.4.0                
##  [77] fastmap_1.1.0                 yaml_2.2.1                   
##  [79] knitr_1.36                    bit64_4.0.5                  
##  [81] purrr_0.3.4                   AnnotationFilter_1.18.0      
##  [83] KEGGREST_1.34.0               sparseMatrixStats_1.6.0      
##  [85] nlme_3.1-153                  mime_0.12                    
##  [87] xml2_1.3.2                    biomaRt_2.50.0               
##  [89] compiler_4.1.1                beeswarm_0.4.0               
##  [91] filelock_1.0.2                curl_4.3.2                   
##  [93] png_0.1-7                     interactiveDisplayBase_1.32.0
##  [95] tibble_3.1.5                  bslib_0.3.1                  
##  [97] stringi_1.7.5                 highr_0.9                    
##  [99] GenomicFeatures_1.46.0        lattice_0.20-45              
## [101] ProtGenerics_1.26.0           Matrix_1.3-4                 
## [103] shinyjs_2.0.0                 vctrs_0.3.8                  
## [105] pillar_1.6.4                  lifecycle_1.0.1              
## [107] BiocManager_1.30.16           jquerylib_0.1.4              
## [109] GlobalOptions_0.1.2           BiocNeighbors_1.12.0         
## [111] irlba_2.3.3                   bitops_1.0-7                 
## [113] rtracklayer_1.54.0            httpuv_1.6.3                 
## [115] BiocIO_1.4.0                  R6_2.5.1                     
## [117] bookdown_0.24                 promises_1.2.0.1             
## [119] gridExtra_2.3                 vipor_0.4.5                  
## [121] codetools_0.2-18              colourpicker_1.1.1           
## [123] assertthat_0.2.1              fontawesome_0.2.2            
## [125] rjson_0.2.20                  withr_2.4.2                  
## [127] shinyWidgets_0.6.2            GenomicAlignments_1.30.0     
## [129] Rsamtools_2.10.0              GenomeInfoDbData_1.2.7       
## [131] hms_1.1.1                     mgcv_1.8-38                  
## [133] parallel_4.1.1                beachmat_2.10.0              
## [135] grid_4.1.1                    DelayedMatrixStats_1.16.0    
## [137] rmarkdown_2.11                Rtsne_0.15                   
## [139] shiny_1.7.1                   ggbeeswarm_0.6.0             
## [141] restfulr_0.0.13Rue-Albrecht, Kevin, Federico Marini, Charlotte Soneson, and Aaron T. L. Lun. 2018. “ISEE: Interactive Summarizedexperiment Explorer.” F1000Research 7 (June): 741. https://doi.org/10.12688/f1000research.14966.1.