Compiled date: 2023-10-24

Last edited: 2018-03-08

License: MIT + file LICENSE

1 Overview

Some tweaks can be performed to enable iSEE to run efficiently on large datasets. This includes datasets with many features (methylation, SNPs) or many columns (cytometry, single-cell RNA-seq). To demonstate some of this functionality, we will use a dataset from the TENxPBMCData dataset:

library(TENxPBMCData)
sce.pbmc <- TENxPBMCData("pbmc68k")
sce.pbmc$Library <- factor(sce.pbmc$Library)
sce.pbmc
#> class: SingleCellExperiment 
#> dim: 32738 68579 
#> metadata(0):
#> assays(1): counts
#> rownames(32738): ENSG00000243485 ENSG00000237613 ... ENSG00000215616
#>   ENSG00000215611
#> rowData names(3): ENSEMBL_ID Symbol_TENx Symbol
#> colnames: NULL
#> colData names(11): Sample Barcode ... Individual Date_published
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):

2 Using out-of-memory matrices

Many SummarizedExperiment objects store assay matrices as in-memory matrix-like objects, be they ordinary matrices or alternative representations such as sparse matrices from the Matrix package. For example, if we looked at the Allen data, we would see that the counts are stored as an ordinary matrix.

library(scRNAseq)
sce.allen <- ReprocessedAllenData("tophat_counts")
class(assay(sce.allen, "tophat_counts"))
#> [1] "matrix" "array"

In situations involving large datasets and limited computational resources, storing the entire assay in memory may not be feasible. Rather, we can represent the data as a file-backed matrix where contents are stored on disk and retrieved on demand. Within the Bioconductor ecosystem, the easiest way of doing this is to create a HDF5Matrix, which uses a HDF5 file to store all the assay data. We see that this has already been done for use in the 68K PBMC dataset:

counts(sce.pbmc, withDimnames=FALSE)
#> <32738 x 68579> sparse HDF5Matrix object of type "integer":
#>              [,1]     [,2]     [,3]     [,4] ... [,68576] [,68577] [,68578]
#>     [1,]        0        0        0        0   .        0        0        0
#>     [2,]        0        0        0        0   .        0        0        0
#>     [3,]        0        0        0        0   .        0        0        0
#>     [4,]        0        0        0        0   .        0        0        0
#>     [5,]        0        0        0        0   .        0        0        0
#>      ...        .        .        .        .   .        .        .        .
#> [32734,]        0        0        0        0   .        0        0        0
#> [32735,]        0        0        0        0   .        0        0        0
#> [32736,]        0        0        0        0   .        0        0        0
#> [32737,]        0        0        0        0   .        0        0        0
#> [32738,]        0        0        0        0   .        0        0        0
#>          [,68579]
#>     [1,]        0
#>     [2,]        0
#>     [3,]        0
#>     [4,]        0
#>     [5,]        0
#>      ...        .
#> [32734,]        0
#> [32735,]        0
#> [32736,]        0
#> [32737,]        0
#> [32738,]        0

Despite the dimensions of this matrix, the HDF5Matrix object occupies very little space in memory.

object.size(counts(sce.pbmc, withDimnames=FALSE))
#> 2496 bytes

However, parts of the data can still be read in on demand. For all intents and purposes, the HDF5Matrix appears to be an ordinary matrix to downstream applications and can be used as such.

first.gene <- counts(sce.pbmc)[1,]
head(first.gene)
#> [1] 0 0 0 0 0 0

This means that we can use the 68K PBMC SingleCellExperiment object in iSEE() without any extra work. The app below shows the distribution of counts for everyone’s favorite gene MALAT1 across libraries. Here, iSEE() is simply retrieving data on demand from the HDF5Matrix without ever loading the entire assay matrix into memory. This enables it to run efficiently on arbitrary large datasets with limited resources.

library(iSEE)
app <- iSEE(sce.pbmc, initial=
    list(RowDataTable(Selected="ENSG00000251562", Search="MALAT1"), 
        FeatureAssayPlot(XAxis="Column data", XAxisColumnData="Library",
            YAxisFeatureSource="RowDataTable1")
    )
)

Generally speaking, these HDF5 files are written once by a process with sufficient computational resources (i.e., memory and time). We typically create HDF5Matrix objects using the writeHDF5Array() function from the HDF5Array package. After the file is created, the objects can be read many times in more deprived environments.

sce.h5 <- sce.allen
library(HDF5Array)
assay(sce.h5, "tophat_counts", withDimnames=FALSE)  <- 
    writeHDF5Array(assay(sce.h5, "tophat_counts"), file="assay.h5", name="counts")
class(assay(sce.h5, "tophat_counts", withDimnames=FALSE))
#> [1] "HDF5Matrix"
#> attr(,"package")
#> [1] "HDF5Array"
list.files("assay.h5")
#> character(0)

It is worth noting that iSEE() does not know or care that the data is stored in a HDF5 file. The app is fully compatible with any matrix-like representation of the assay data that supports dim() and [,. As such, iSEE() can be directly used with other memory-efficient objects like the DeferredMatrix and LowRankMatrix from the BiocSingular package, or perhaps the ResidualMatrix from the batchelor package.

3 Downsampling points

It is also possible to downsample points to reduce the time required to generate the plot. This involves subsetting the dataset so that only the most recently plotted point for an overlapping set of points is shown. In this manner, we avoid wasting time in plotting many points that would not be visible anyway. To demonstrate, we will re-use the 68K PBMC example and perform downsampling on the feature assay plot; we can see that its aesthetics are largely similar to the non-downsampled counterpart above.

library(iSEE)
app <- iSEE(sce.pbmc, initial=
    list(RowDataTable(Selected="ENSG00000251562", Search="MALAT1"), 
        FeatureAssayPlot(XAxis="Column data", XAxisColumnData="Library",
            YAxisFeatureSource="RowDataTable1", 
            VisualChoices="Point", Downsample=TRUE,
            VisualBoxOpen=TRUE
        )
    )
)

Downsampling is possible in all iSEE() plotting panels that represent features or samples as points. We can turn on downsampling for all such panels using the relevant field in panelDefaults(), which spares us the hassle of setting Downsample= individually in each panel constructor.

panelDefaults(Downsample=TRUE)

The downsampling only affects the visualization and the speed of the plot rendering. Any interactions with other panels occur as if all of the points were still there. For example, if one makes a brush, all of the points therein will be selected regardless of whether they were downsampled.

The downsampling resolution determines the degree to which points are considered to be overlapping. Decreasing the resolution will downsample more aggressively, improving plotting speed but potentially affecting the fidelity of the visualization. This may compromise the aesthetics of the plot when the size of the points is small, in which case an increase in resolution may be required at the cost of speed.

Obviously, downsampling will not preserve overlays for partially transparent points, but any reliance on partial transparency is probably not a good idea in the first place when there are many points.

4 Changing the interface

One can generally improve the speed of the iSEE() interface by only initializing the app with the desired panels. For example, it makes little sense to spend time rendering a RowDataPlot when only the ReducedDimensionPlot is of interest. Specification of the initial state is straightforward with the initial= argument, as described in a previous vignette.

On occasion, there may be alternative panels with more efficient visualizations for the same data. The prime example is the ReducedDimensionHexPlot class from the iSEEu package; this will create a hexplot rather than a scatter plot, thus avoiding the need to render each point in the latter.

5 Comments on deployment

It is straightforward to host iSEE applications on hosting platforms like Shiny Server or Rstudio Connect. All one needs to do is to create an app.R file that calls iSEE() with the desired parameters, and then follow the instructions for the target platform. For a better user experience, we suggest setting a minimum number of processes to avoid the initial delay from R start-up.

It is also possible to deploy and host Shiny app on shinyapps.io, a platform as a service (PaaS) provided by RStudio. In many cases, users will need to configure the settings of their deployed apps, in particular selecting larger instances to provide sufficient memory for the app. The maximum amount of 1GB available to free accounts may not be sufficient to deploy large datasets; in which case you may consider using out-of-memory matrices, filtering your dataset (e.g., removing lowly detected features), or going for a paid account. Detailed instructions to get started are available at https://shiny.rstudio.com/articles/shinyapps.html. For example, see the isee-shiny-contest app, winner of the 1st Shiny Contest.

Session Info

sessionInfo()
#> R version 4.3.1 (2023-06-16)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.18-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] TENxPBMCData_1.19.0         HDF5Array_1.30.0           
#>  [3] rhdf5_2.46.0                DelayedArray_0.28.0        
#>  [5] SparseArray_1.2.0           S4Arrays_1.2.0             
#>  [7] abind_1.4-5                 Matrix_1.6-1.1             
#>  [9] scater_1.30.0               ggplot2_3.4.4              
#> [11] scuttle_1.12.0              scRNAseq_2.15.0            
#> [13] iSEE_2.14.0                 SingleCellExperiment_1.24.0
#> [15] SummarizedExperiment_1.32.0 Biobase_2.62.0             
#> [17] GenomicRanges_1.54.0        GenomeInfoDb_1.38.0        
#> [19] IRanges_2.36.0              S4Vectors_0.40.0           
#> [21] BiocGenerics_0.48.0         MatrixGenerics_1.14.0      
#> [23] matrixStats_1.0.0           BiocStyle_2.30.0           
#> 
#> loaded via a namespace (and not attached):
#>   [1] splines_4.3.1                 later_1.3.1                  
#>   [3] BiocIO_1.12.0                 bitops_1.0-7                 
#>   [5] filelock_1.0.2                tibble_3.2.1                 
#>   [7] XML_3.99-0.14                 lifecycle_1.0.3              
#>   [9] doParallel_1.0.17             lattice_0.22-5               
#>  [11] ensembldb_2.26.0              magrittr_2.0.3               
#>  [13] sass_0.4.7                    rmarkdown_2.25               
#>  [15] jquerylib_0.1.4               yaml_2.3.7                   
#>  [17] httpuv_1.6.12                 DBI_1.1.3                    
#>  [19] RColorBrewer_1.1-3            zlibbioc_1.48.0              
#>  [21] Rtsne_0.16                    purrr_1.0.2                  
#>  [23] AnnotationFilter_1.26.0       RCurl_1.98-1.12              
#>  [25] rappdirs_0.3.3                circlize_0.4.15              
#>  [27] GenomeInfoDbData_1.2.11       ggrepel_0.9.4                
#>  [29] irlba_2.3.5.1                 DelayedMatrixStats_1.24.0    
#>  [31] codetools_0.2-19              DT_0.30                      
#>  [33] xml2_1.3.5                    tidyselect_1.2.0             
#>  [35] shape_1.4.6                   ScaledMatrix_1.10.0          
#>  [37] viridis_0.6.4                 shinyWidgets_0.8.0           
#>  [39] BiocFileCache_2.10.0          GenomicAlignments_1.38.0     
#>  [41] jsonlite_1.8.7                GetoptLong_1.0.5             
#>  [43] BiocNeighbors_1.20.0          ellipsis_0.3.2               
#>  [45] iterators_1.0.14              foreach_1.5.2                
#>  [47] tools_4.3.1                   progress_1.2.2               
#>  [49] Rcpp_1.0.11                   glue_1.6.2                   
#>  [51] gridExtra_2.3                 xfun_0.40                    
#>  [53] mgcv_1.9-0                    dplyr_1.1.3                  
#>  [55] shinydashboard_0.7.2          withr_2.5.1                  
#>  [57] BiocManager_1.30.22           fastmap_1.1.1                
#>  [59] rhdf5filters_1.14.0           fansi_1.0.5                  
#>  [61] shinyjs_2.1.0                 digest_0.6.33                
#>  [63] rsvd_1.0.5                    R6_2.5.1                     
#>  [65] mime_0.12                     colorspace_2.1-0             
#>  [67] biomaRt_2.58.0                RSQLite_2.3.1                
#>  [69] utf8_1.2.4                    generics_0.1.3               
#>  [71] rtracklayer_1.62.0            prettyunits_1.2.0            
#>  [73] httr_1.4.7                    htmlwidgets_1.6.2            
#>  [75] pkgconfig_2.0.3               gtable_0.3.4                 
#>  [77] blob_1.2.4                    ComplexHeatmap_2.18.0        
#>  [79] XVector_0.42.0                htmltools_0.5.6.1            
#>  [81] bookdown_0.36                 ProtGenerics_1.34.0          
#>  [83] rintrojs_0.3.2                clue_0.3-65                  
#>  [85] scales_1.2.1                  png_0.1-8                    
#>  [87] knitr_1.44                    rjson_0.2.21                 
#>  [89] nlme_3.1-163                  curl_5.1.0                   
#>  [91] shinyAce_0.4.2                cachem_1.0.8                 
#>  [93] GlobalOptions_0.1.2           stringr_1.5.0                
#>  [95] BiocVersion_3.18.0            parallel_4.3.1               
#>  [97] miniUI_0.1.1.1                vipor_0.4.5                  
#>  [99] AnnotationDbi_1.64.0          restfulr_0.0.15              
#> [101] pillar_1.9.0                  grid_4.3.1                   
#> [103] vctrs_0.6.4                   promises_1.2.1               
#> [105] BiocSingular_1.18.0           dbplyr_2.3.4                 
#> [107] beachmat_2.18.0               xtable_1.8-4                 
#> [109] cluster_2.1.4                 beeswarm_0.4.0               
#> [111] evaluate_0.22                 GenomicFeatures_1.54.0       
#> [113] cli_3.6.1                     compiler_4.3.1               
#> [115] Rsamtools_2.18.0              rlang_1.1.1                  
#> [117] crayon_1.5.2                  ggbeeswarm_0.7.2             
#> [119] stringi_1.7.12                viridisLite_0.4.2            
#> [121] BiocParallel_1.36.0           munsell_0.5.0                
#> [123] Biostrings_2.70.0             lazyeval_0.2.2               
#> [125] colourpicker_1.3.0            ExperimentHub_2.10.0         
#> [127] hms_1.1.3                     sparseMatrixStats_1.14.0     
#> [129] bit64_4.0.5                   Rhdf5lib_1.24.0              
#> [131] KEGGREST_1.42.0               shiny_1.7.5.1                
#> [133] interactiveDisplayBase_1.40.0 AnnotationHub_3.10.0         
#> [135] fontawesome_0.5.2             igraph_1.5.1                 
#> [137] memoise_2.0.1                 bslib_0.5.1                  
#> [139] bit_4.0.5
# devtools::session_info()