vignettes/iSEE-lab.Rmd
iSEE-lab.Rmd
This workshop demonstrates how to use the iSEE package to create and configure interactive applications for the exploration of various types of genomics data sets (e.g., bulk and single-cell RNA-seq, CyTOF, gene expression microarray).
This workshop will be presented as a lab session that combines an instructor-led live demo, followed by hands-on experimentation guided by completely worked examples and stand-alone notes that participants may continue to use after the workshop.
The instructor-led live demo comprises three parts:
The hands-on lab comprises three parts:
Participants are encouraged to ask questions at any time during the workshop.
Additional background reading about the programming environment, relevant packages, and example use cases:
Students will participate by following along an R markdown (https://rmarkdown.rstudio.com/) document, and asking questions throughout the workshop. There is also scope for participants to apply iSEE to their own data sets, and fuel the discussion with more questions about specific use cases.
In this workshop, we use example data from the TENxPBMCData package. This package provides an R / Bioconductor resource for representing and manipulating different single-cell RNA-seq data sets profiling peripheral blood mononuclear cells (PBMC) generated by 10x Genomics (https://support.10xgenomics.com/single-cell-gene-expression/datasets).
library(TENxPBMCData)
The man page for the TENxPBMCData()
function gives an idea of the datasets that are available from this package. It can be opened with the following command.
help(TENxPBMCData)
Here, we use the pbmc3k
dataset, which contains gene expression profiles for 2,700 single peripheral blood mononuclear cells. The first time this dataset is loaded, this command downloads the dataset to a local cache, which takes some time, depending on the speed of your internet connection. Subsequent times, the same command loads the dataset directly from the local cache.
sce <- TENxPBMCData(dataset = "pbmc3k")
At this point we can inspect the dataset in the console.
sce #> class: SingleCellExperiment #> dim: 32738 2700 #> 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): #> altExpNames(0):
The dataset is provided as an object of the SingleCellExperiment
class. In particular, this summary view indicates that the following pieces of information are available:
counts
NULL
Barcode
) and the donor identifier (Individual
).Note that a SingleCellExperiment
object (or, more generally, any SummarizedExperiment
object) like this one already contains sufficient information to launch an interactive application instance to visualize the available data and metadata, using the iSEE()
function.
For the purpose of this workshop, we first apply some preprocessing to the SingleCellExperiment
object, in order to populate it with more information that can be visualized with iSEE
.
We start by adding column names to the object, and use gene symbols instead of Ensembl IDs as row names. In the case where multiple Ensembl identifiers correspond to the same gene symbol, the scater::uniquifyFeatureNames
function concatenates the Ensembl ID and the gene symbol in order to generate unique feature names.
library(scater) colnames(sce) <- paste0("Cell", seq_len(ncol(sce))) rownames(sce) <- scater::uniquifyFeatureNames( ID = rowData(sce)$ENSEMBL_ID, names = rowData(sce)$Symbol_TENx ) head(rownames(sce)) #> [1] "MIR1302-10" "FAM138A" "OR4F5" "RP11-34P13.7" "RP11-34P13.8" #> [6] "AL627309.1"
Next, we use the scater package to calculate gene- and cell-level quality metrics. These metrics are added as columns to the rowData
and colData
slots of the SingleCellExperiment
object, respectively. We also add some additional metrics that are not automatically computed by the scater package.
MT <- rownames(sce)[grep("^MT-", rownames(sce))] sce <- scater::addPerCellQC(sce, subsets = list(MT = MT)) #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") sce <- scater::addPerFeatureQC(sce) #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") sce$log10_total <- log10(sce$total) rowData(sce)$n_cells <- as.integer(rowData(sce)$detected * ncol(sce)) rowData(sce)$log10_total <- log10(rowSums(assay(sce, "counts")) + 1) sce #> class: SingleCellExperiment #> dim: 32738 2700 #> metadata(0): #> assays(1): counts #> rownames(32738): MIR1302-10 FAM138A ... AC002321.2 AC002321.1 #> rowData names(7): ENSEMBL_ID Symbol_TENx ... n_cells log10_total #> colnames(2700): Cell1 Cell2 ... Cell2699 Cell2700 #> colData names(18): Sample Barcode ... total log10_total #> reducedDimNames(0): #> altExpNames(0):
We filter out a few cells with a large fraction of the counts coming from mitochondrial genes, since these may be damaged cells. Notice the reduced number of columns in the dataset below.
(sce <- sce[, sce$subsets_MT_percent < 5]) #> class: SingleCellExperiment #> dim: 32738 2643 #> metadata(0): #> assays(1): counts #> rownames(32738): MIR1302-10 FAM138A ... AC002321.2 AC002321.1 #> rowData names(7): ENSEMBL_ID Symbol_TENx ... n_cells log10_total #> colnames(2643): Cell1 Cell2 ... Cell2699 Cell2700 #> colData names(18): Sample Barcode ... total log10_total #> reducedDimNames(0): #> altExpNames(0):
Next, we calculate size factors and normalized and log-transformed expression values, using the scran and scater packages. Note that it is typically recommended to pre-cluster the cells before computing the size factors, as follows:
# set.seed(1000) # clusters <- scran::quickCluster(sce, BSPARAM = IrlbaParam()) # sce <- scran::computeSumFactors(sce, cluster = clusters, min.mean = 0.1)
However, for time reasons, we will skip the pre-clustering step in this workshop.
library(scran) sce <- scran::computeSumFactors(sce, min.mean = 0.1) #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") summary(sizeFactors(sce)) #> Min. 1st Qu. Median Mean 3rd Qu. Max. #> 0.1164 0.7309 0.9004 1.0000 1.1635 9.8580 sce <- scater::logNormCounts(sce)
In order to extract the most informative genes, we first model the mean-variance trend and decompose the variance into biological and technical components.
dec <- scran::modelGeneVar(sce) #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") top.dec <- dec[order(dec$bio, decreasing = TRUE), ] head(top.dec) #> DataFrame with 6 rows and 6 columns #> mean total tech bio p.value FDR #> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric> #> LYZ 1.63950 3.92209 0.884955 3.03713 4.59557e-117 2.54472e-113 #> S100A9 1.07323 3.36899 0.778205 2.59079 2.72222e-110 9.04431e-107 #> HLA-DRA 1.56937 3.16577 0.880146 2.28562 5.73028e-68 7.93262e-65 #> FTL 3.67846 2.88499 0.708226 2.17676 2.46520e-94 5.85028e-91 #> CD74 2.23653 2.75839 0.854474 1.90392 1.34728e-50 1.01732e-47 #> CST3 1.10995 2.59172 0.790231 1.80149 7.24671e-53 6.01912e-50
Next, we apply Principal Components Analysis (PCA) and t-distributed Stochastic Neighbor Embedding (t-SNE) to generate low-dimensional representations of the cells in our data set. These low-dimensional representations are added to the reducedDim
slot of the SingleCellExperiment
object.
library(BiocSingular) set.seed(1000) sce <- scran::denoisePCA(sce, technical = dec) #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") ncol(reducedDim(sce, "PCA")) #> [1] 5 set.seed(1000) sce <- scater::runTSNE(sce, dimred = "PCA", perplexity = 30) sce <- scater::runUMAP(sce, dimred = "PCA") sce #> class: SingleCellExperiment #> dim: 32738 2643 #> metadata(0): #> assays(2): counts logcounts #> rownames(32738): MIR1302-10 FAM138A ... AC002321.2 AC002321.1 #> rowData names(7): ENSEMBL_ID Symbol_TENx ... n_cells log10_total #> colnames(2643): Cell1 Cell2 ... Cell2699 Cell2700 #> colData names(19): Sample Barcode ... log10_total sizeFactor #> reducedDimNames(3): PCA TSNE UMAP #> altExpNames(0):
Finally, we cluster the cells using a graph-based algorithm, and find ‘marker genes’ for each cluster as the genes that are significantly upregulated in the cluster compared to each of the other inferred clusters. The adjusted p-values from this test, for each cluster, are added to the rowData
slot of the object.
snn.gr <- scran::buildSNNGraph(sce, use.dimred = "PCA") clusters <- igraph::cluster_walktrap(snn.gr) sce$Cluster <- factor(clusters$membership) table(sce$Cluster) #> #> 1 2 3 4 5 6 7 8 9 10 11 12 13 14 #> 420 354 143 212 161 318 324 189 150 144 117 75 25 11 markers <- scran::findMarkers(sce, groups = sce$Cluster, test.type = "t", direction = "up", pval.type = "all") #> Warning: 'rowGrid' is deprecated. #> Use 'rowAutoGrid' instead. #> See help("Deprecated") #> Warning: 'colGrid' is deprecated. #> Use 'colAutoGrid' instead. #> See help("Deprecated") for (i in names(markers)) { rowData(sce)[, paste0("FDR_cluster", i)] <- markers[[i]]$FDR[match(rownames(sce), rownames(markers[[i]]))] } sce #> class: SingleCellExperiment #> dim: 32738 2643 #> metadata(0): #> assays(2): counts logcounts #> rownames(32738): MIR1302-10 FAM138A ... AC002321.2 AC002321.1 #> rowData names(21): ENSEMBL_ID Symbol_TENx ... FDR_cluster13 #> FDR_cluster14 #> colnames(2643): Cell1 Cell2 ... Cell2699 Cell2700 #> colData names(20): Sample Barcode ... sizeFactor Cluster #> reducedDimNames(3): PCA TSNE UMAP #> altExpNames(0):
This concludes the preparation of the data. We have now a SingleCellExperiment
object that contains different types of abundance values, representations in reduced dimensions, as well as a range of row (feature) and column (cell) metadata.
We can launch an iSEE instance for exploring this data set using the iSEE()
function:
shiny::runApp(app)
The main argument to the iSEE()
function is a SummarizedExperiment
object, or an object of any class extending SummarizedExperiment
(such as SingleCellExperiment
, in this case). No other restrictions are made on the type of data stored in the object, and iSEE can be used for interactive visualization of many different types of data. It is also worth noting that for various types of data, Bioconductor packages provides functionality for directly importing quantifications generated by external software packages into a SummarizedExperiment
object. For example, the DropletUtils package can read quantifications from the 10x Genomics CellRanger pipeline for single-cell RNA-seq data, and the tximeta package can be used to read data from transcript quantification pipelines into a SummarizedExperiment
object.
While we will not make explicit use of it in this workshop, we note that it is also possible to call iSEE()
without providing a SummarizedExperiment
object. In that case, the user will be prompted to upload such an object, serialized into a .rds
file. It is also possible to import a specification of the initial panel setup.
app <- iSEE()
shiny::runApp(app)
This section provides an overview of the graphical interface of iSEE applications. To follow along, make sure that you have launched an iSEE instance as described in the code block at the end of the previous section. Note that in the default configuration, the panels do not look exactly like the ones shown in the screenshots below. For example, data points are not immediately colored, and the default annotation variables displayed by each panel may differ. Later sections of this workshop demonstrate how to modify the content of the panels and how they are displayed.
Note that for simplicity, we typically refer to a SummarizedExperiment
in this workshop; however, iSEE works seamlessly for objects of any class extending SummarizedExperiment
as well (e.g., SingleCellExperiment
, DESeqDataSet
). That said, some types of panels – such as the Reduced dimension plot
– are only available for objects that contain a reducedDim
slot (in particular, SingleCellExperiment
objects); the basic SummarizedExperiment
class does not contain this slot. In this workshop, we refer to the rows of the SummarizedExperiment
object as ‘features’ (these can be genes, transcripts, genomic regions, etc) and to the columns as ‘samples’ (which, in our example data set, are single cells).
The iSEE user interface consists of a number of panels, each displaying the data provided in the SummarizedExperiment
from a specific perspective. There are 8 standard panel types; 6 plot panels and 2 table panels, all showcased in the figure shown at the end of section Preparation of example scRNAseq data. In the default configuration, one panel of each type is included when launching the iSEE user interface. However, users are free to rearrange, resize, remove or add panels freely, as described in the later section Configuration of the app interface. We provide a brief overview of each panel type in the following subsections.
In addition to the 8 standard panel types included in iSEE, users can create custom panels (both plots and tables). Moreover, the iSEEu (iSEE universe) package contains additional panel types. The creation and configuration of custom panels is discussed later, in the dedicated section Custom panels.
The reduced dimension plot can display any reduced dimension representation that is present in the reducedDim
slot of the SingleCellExperiment
object.
Note that this slot is not defined for the base SummarizedExperiment
class, in which case the user interface does not allow the inclusion of panels of this type.
The column data plot can display one or two of the provided column annotations (from the colData
slot). Depending on the class of the selected annotations, the panel shows either a scatter plot, a violin plot, or a Hinton diagram (Hinton and Shallice 1991; Bremner, Gotts, and Denham 1994).
Analogous to the column data plot above, the row data plot displays one or two of the provided row annotations (from the rowData
slot). Depending on the class of the selected annotations, the panel displays either a scatter plot, a violin plot, or a Hinton plot.
The complex heatmap panel displays, for any assay, the observed values for a subset of the features across the samples.
The feature assay plot displays the observed values for one feature across the samples. It is also possible to plot the observed values for two features, in a scatter plot.
Analogous to the feature assay plot above, the sample assay plot shows the observed values for all features, for one of the samples. It is also possible to plot the observed values for two samples, in a scatter plot.
The row data table displays all information provided in the rowData
slot of the SummarizedExperiment
object, leveraging the interactivity provided by the DT package.
Analogous to the row data table above, the column data table displays all information provided in the colData
slot of the SummarizedExperiment
object.
Each plot panel type has a Data parameters
collapsible box. This box has different content for each panel type, but in all cases it lets the user control the data that is displayed in the plot.
In contrast to the Data parameters
collapsible box that lets users control what is displayed in the plot, the Visual parameters
box lets users control how the information is displayed.
This collapsible box contains the controls to change the size, shape, opacity, and color of the points, to facet the plot by any available categorical annotation, to subsample points for increased speed of plot rendering, and to control how legends are displayed.
The Selection parameters
collapsible box provides controls to transfer selections of points (features or samples) between panels.
We discuss point transmission in more detail later in this workshop.
At the top-right corner of the iSEE application, users can find additional controls for reproducibility, configuration, and help. Each of these is discussed in more detail later in this workshop. Links are provided in the last column of the table below.
Organization |
Organize panels Examine panel chart |
|
Export |
Download panel output Extract R code Display panel settings |
|
Documentation |
Quick tour Open the vignette |
|
Additional info |
About this session About iSEE |
As mentioned above, the default behaviour of the iSEE()
function is to launch an instance of the user interface that displays one panel of each of the standard types (provided the underlying data is available, e.g. for reduced dimension plots the reducedDim
slot is required). However, in some cases it is desirable to have multiple panels of the same type, and/or exclude some panel types.
In order to accommodate such situations, users can add, remove, change the order of and resize all panels via the Organization
menu in the top-right corner.
Action: Using your open iSEE instance, open the
Organization
menu () and click on theOrganize panels
() button. Try to add and remove panels, and resize the existing ones. Remember to click onApply settings
() to apply the changes you have made to the interface.
Clicking in the selectize box listing all current panels will present you with a drop-down menu from which you can choose additional panels to add. Similarly, panels can be removed by clicking on the icon associated with the panel name.
Each panel can be individually resized by changing the width and height. Note that the total width of a row in the interface is 12 units. When the width of a panel is greater than the space available, the panel is moved to a new row.
It can take a great amount of time to achieve a satisfactory panel configuration. To avoid the need to manually organize the panels each time the app is opened, iSEE offers the possibility to specify how the app is initialized, as well as to inspect and export the current panel settings for future use.
As a first example, let’s assume that we would like to start iSEE
with the default set of eight panels, but we want each panel to have width equal to 3 units rather than the default 4, so that we can fit all the panels in two rows. To achieve this, we define a list that we name initial
, which contains the specifications for all the panels that should be included at startup.
initial <- list( ReducedDimensionPlot(PanelWidth = 3L), ColumnDataPlot(PanelWidth = 3L), RowDataPlot(PanelWidth = 3L), ComplexHeatmapPlot(PanelWidth = 3L), FeatureAssayPlot(PanelWidth = 3L), SampleAssayPlot(PanelWidth = 3L), RowDataTable(PanelWidth = 3L), ColumnDataTable(PanelWidth = 3L) )
Action: Close the application, run the code above in your R session, and relaunch the application, now with the additional argument
initial
set:
app <- iSEE(sce, initial = initial)
shiny::runApp(app)
Notice how the panel configuration at startup is now different from the default one. The setup is fully customizable - it is possible to include multiple panels of each types, and the size specifications of a given panel are independent of those for the other panels. It is even possible to start an app without any panels at all, by simply providing an empty list to initial
(see also iSEEu::modeEmpty()
). Regardless of the initialization choice, once the app is launched, the user by default retains full flexibility to reorganize the interface interactively (it is possible to limit the types of panels that can be added via the extra
argument to iSEE
).
Now that we have an application consisting of a suitable collection of panels for our data set, we’ll explore the collapsible boxes below each plot, which lets us control their appearance.
By default, each Reduced dimension plot displays the first reduced dimension representation provided (contained in the reducedDim
slot of the SingleCellExperiment
object). If additional reduced dimension representations are included, the one chosen for the display can be changed in the Data parameters
collapsible box in the Reduced dimension plot panel.
Action: Use this control to show the t-SNE representation instead of the PCA.
Notably, the instance can be preconfigured to immediately display the t-SNE representation when the application is launched.
Action: Close the application and run the following block of code to show two Reduced dimension plot panels; one showing the PCA representation, one showing the t-SNE representation.
initial <- list( ReducedDimensionPlot(Type = "PCA", PanelWidth = 6L), ReducedDimensionPlot(Type = "TSNE", PanelWidth = 6L) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
With default settings, all points are colored black. However, any provided row (for features) or column (for samples) annotation can be used to color the points.
Action: To demonstrate this, open the
Visual parameters
collapsible box in one of the Reduced dimension plot panels, and setColor by:
toColumn data
. Scroll through the dropdown menu that appears in order to see all the sample annotations that can be used to color the points. All these values are taken from thecolData
slot of the provided object.
Action: First, color the cells by
Cluster
, which contains the cluster labels that were assigned to the cells in the preprocessing step in a previous section. Note how a discrete color scheme is selected, since the variable is categorical.
Action: Next, color the cells by the log10-transformed number of detected genes (
log10_total
). Notice how the color scale is now continuous.
The instance can be preconfigured to control the coloring of data points, to display the same result on startup.
Action: Close the application and run the following block of code to edit the configuration of the two Reduced dimension plot panels prepared in the previous section.
initial <- list( ReducedDimensionPlot(Type = "PCA", ColorBy = "Column data", ColorByColumnData = "log10_total", PanelWidth = 6L), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", PanelWidth = 6L) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
iSEE allows precise control over the color schemes used. More information about this can be found in the ExperimentColorMap vignette.
In addition to coloring cells by the provided annotations as above, we can also color them by the observed values for a given feature (as provided in one of the assays
).
Action: Set
Color by:
toFeature name
in one of the Reduced dimension plot panels, and select one of the genes from the dropdown menu. You can search for a gene of interest by typing in the dropdown box. You can also choose which assay should be used to extract the values to color according to.
Again, the same result can be achieved by preconfiguring the application as follows.
Action: Close the application and run the following block of code to set up 3 panels, all displaying the t-SNE representation. In the first panel, we color cells by their
Cluster
label. In the second and third panel, we display CD79A and CD74, a combination of gene markers typically observed in B cells.
## The displayed features IDs must be row names of the object c("CD79A", "CD74") %in% rownames(sce) #> [1] TRUE TRUE initial <- list( ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster"), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Feature name", ColorByFeatureName = "CD79A", ColorByFeatureNameAssay = "logcounts"), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Feature name", ColorByFeatureName = "CD74", ColorByFeatureNameAssay = "logcounts") ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
The Color by: Sample name
option in the Reduced dimension plot can also be used to highlight a single sample.
Action: Set
Color by:
toSample name
in one of the Reduced dimension plots, and use the dropdown menu underneath to change the highlighted cell.
Again, the same result can be achieved by preconfiguring the application.
Action: Close the app and run the following block of code to continue editing the configuration of the Reduced dimension plot panels prepared in the previous section. We also open the Visual parameters collapsible box to highlight the preconfigured options.
initial <- list( ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", VisualBoxOpen = TRUE), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Sample name", ColorBySampleName = "Cell2439", VisualBoxOpen = TRUE), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Sample name", ColorBySampleName = "Cell673", VisualBoxOpen = TRUE) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
By default, only the Color
options are displayed in the Visual parameters
box. Additional controls can be shown by clicking the check boxes for Shape
, Size
, Point
, Facet
and Text
. The Size
, Shape
and Point
panels let you set the shape and size of points according to provided annotations (or to resize all points), and to set the point opacity.
Here, you can also downsample points to increase the speed of rendering, which can be useful for data sets with large numbers of points. Crucially, downsampling is purely a visual effect; data points hidden by the downsampling process are still included in selections captured by Shiny brushes and lasso selections. We illustrate this point by applying the same Shiny brush to select a visually well defined cluster visually downsampled to a resolution of 100 and 50 bins, respectively; note how the information message for both selection declares the same number of data points selected, despite the clearly different number of data points displayed in the corresponding area of each plot.
Action: Close the app and run the following block of code to continue editing the configuration of the Reduced dimension plot panels prepared in the previous section. Here, we reduce the point size and increase the transparency (i.e. reduce the value for the alpha attribute). We also apply a downsampling grid of 100 horizontal and vertical bins to the plot, where only a single data point is plotted (if any) for each bin. Finally, we display the
Shape
andPoint
options.
brushData <- list(xmin = -16.44, xmax = 5.6976316150235, ymin = -41.147401279341, ymax = -15.99494636763, coords_css = list(xmin = 231.34375, xmax = 383.34375, ymin = 306L, ymax = 417L), coords_img = list(xmin = 231.34375, xmax = 383.34375, ymin = 306L, ymax = 417L), img_css_ratio = list(x = 1L, y = 1L), mapping = list(x = "X", y = "Y", colour = "ColorBy"), domain = list(left = -44.5050628659896, right = 38.791275779239, bottom = -43.4753804355979, top = 47.8385432957266), range = list(left = 38.7190443065069, right = 610.520547945205, bottom = 427.273577161815, top = 24.2971850062453), log = list(x = NULL, y = NULL), direction = "xy", brushId = "ReducedDimensionPlot2_Brush", outputId = "ReducedDimensionPlot2") initial <- list( ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", VisualBoxOpen = TRUE, Downsample = TRUE, DownsampleResolution = 200), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Sample name", ColorBySampleName = "Cell2439", VisualBoxOpen = TRUE, PointSize = 0.4, PointAlpha = 0.4, Downsample = TRUE, DownsampleResolution = 100, VisualChoices = c("Shape", "Point"), BrushData = brushData), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Sample name", ColorBySampleName = "Cell673", VisualBoxOpen = TRUE, PointSize = 0.4, PointAlpha = 0.4, Downsample = TRUE, DownsampleResolution = 50, VisualChoices = c("Shape", "Point"), BrushData = brushData) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
The Column data plot and Row data plot panels are used to display the information in the colData
and rowData
slots of the SummarizedExperiment
object. These panels can display either one or two annotations. By default, the first annotation in the respective metadata table is selected.
The plot layout adapts automatically to the type of annotation(s) displayed.
As a side note, categorical variables with “too many” unique values (which is typically the case for sample IDs, for instance) are represented as numeric variables, by replacing each value by an index. The limit for what is considered “too many” unique values can be set using options(iSEE.maxlevels = 30)
.
Action: Start an instance with the default set of panels and explore the different types of plots provided in the Column data or Row data plot panels by selecting different annotations as the
Column of interest (Y-axis)
and/or settingX-axis:
toColumn data
and selecting an additional annotation.
In the following code chunk we demonstrate the preconfiguration of an app instance that immediately displays two panels.
Action: Close the app and run the following block of code to set up one Column data plot and one Row data plot panel. We set up the Column data plot panel to immediately examine the proportion of counts assigned to mitochondrial features in each cluster. We set up the Row data plot panel to examine the proportion of cells with detectable expression against the log-transformed total counts for each gene.
initial <- list( ColumnDataPlot(YAxis = "subsets_MT_percent", XAxis = "Column data", XAxisColumnData = "Cluster", PanelWidth = 6L), RowDataPlot(YAxis = "n_cells", XAxis = "Row data", XAxisRowData = "log10_total", PanelWidth = 6L) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
The Feature assay plot is used to display the observed values for one feature across the samples on the Y-axis. By default, the first feature in the rownames
of the dataset is selected, displaying its distribution across the entire dataset.
While continuous measurements are more common (e.g., gene expression), discrete assays are also supported (e.g., single nucleotide polymorphisms), and displayed as a Hinton plot, consisting of rectangles with an area proportional to the number of samples they represent.
It is also possible to plot the observed values of a second feature on the X-axis. This will result in a scatter plot or Hinton plot, depending on the continuous or discrete nature of the data, respectively.
Alternatively, the X-axis can also be used to stratify gene expression according to a discrete sample metadata in the form of a violin plot, while a continuous metadata will generate a scatter plot.
Action: Start an instance with the default set of panels and explore the different types of plots generated by the Feature assay plot panel using the choices available for the two axes.
Similarly to all other panels, the Feature assay plot can be preconfigured to immediately display specific information when the app is launched.
Action: Close the app and run the following block of code to set up three Feature assay plot panels demonstrating different configurations. All three panels display the same gene on the Y-axis and color data points according to their cluster label. On the X-axis, the first panel displays a discrete sample metadata variable, the cluster label; the second panel display a continuous sample metadata variable, the log-transformed library size; the third panel displays the expression of a second gene.
initial <- list( FeatureAssayPlot(ColorBy = "Column data", ColorByColumnData = "Cluster", YAxisFeatureName = "CD3D", XAxis = "Column data", XAxisColumnData = "Cluster", DataBoxOpen = TRUE, Assay = "logcounts"), FeatureAssayPlot(ColorBy = "Column data", ColorByColumnData = "Cluster", YAxisFeatureName = "CD3D", XAxis = "Column data", XAxisColumnData = "log10_total", DataBoxOpen = TRUE, Assay = "logcounts"), FeatureAssayPlot(ColorBy = "Column data", ColorByColumnData = "Cluster", YAxisFeatureName = "CD3D", XAxis = "Feature name", XAxisFeatureName = "CD40LG", DataBoxOpen = TRUE, Assay = "logcounts") ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
When exploring data, it is often useful to be able to select a subset of points and investigate different aspects of their characteristic features. In iSEE, this can be achieved by selecting points in one panel, and transmitting this selection to one or more other panels.
Action: Select a set of points in one of the Feature assay plot panels, by drawing a rectangle around them. Then open the Selection parameters collapsible box of one of the other Feature assay plot panels, and select the Feature assay plot panel where you made the point selection from the dropdown menu under Receive column selection from. Note how the points corresponding to the cells that you selected in the first panel are highlighted in the receiving panel. You can highlight the points in different ways by changing the Selection effect.
Also the brushing and point selection can be preconfigured. In the example below, we configure three types of panels: a Feature assay plot, a Reduced dimension plot, and a Column data plot. For visual convenience, we color data points in all three panels by the cluster label.
To demonstrate the feature, we define a Shiny brush in the Feature assay plot panel to select cells with detectable levels of the CD3D gene. For the other two panels, pay particular attention to the "ColumnSelectionSource"
instruction which declares that they should both receive the data points selected in the Feature assay plot panel. In addition, we specify that the Reduced dimension plot should highlight the selection by applying a transparency effect to unselected data points. Alternatively, we use the Column data plot to demonstrate how the selection can also be highlighted by applying a given color (here, red) to the selected data points.
Action: Close the app and run the following block of code to launch the example application described above.
initial <- list( FeatureAssayPlot(YAxisFeatureName = "CD3D", BrushData = list( xmin = 0.5, xmax = 14.5, ymin = 0.5, ymax = 5.7, mapping = list(x = "X", y = "Y", group = "GroupBy"), direction = "xy", brushId = "featAssayPlot1_Brush", outputId = "featAssayPlot1"), XAxis = "Column data", XAxisColumnData = "Cluster", ColorBy = "Column data", ColorByColumnData = "Cluster"), ReducedDimensionPlot(ColumnSelectionSource = "FeatureAssayPlot1", Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", SelectionBoxOpen = TRUE), ColumnDataPlot(ColumnSelectionSource = "FeatureAssayPlot1", SelectionEffect = "Color", YAxis = "log10_total", XAxis = "Column data", XAxisColumnData = "subsets_MT_percent", ColorBy = "Column data", ColorByColumnData = "Cluster", SelectionBoxOpen = TRUE) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
Notably, iSEE can display active links transmitting point selections between panel using the “Examine panel chart” menu.
iSEE allows you to zoom into a subregion of a plot, by simply drawing a rectangular selection around the area of interest and double-clicking inside the rectangle. To zoom back out to the full plot, double-click anywhere in the plot.
To demonstrate the feature, we set up two Reduced dimension plot panels. In the first panel, we draw a Shiny brush to indicate an area of interest in the reduced dimension plot. We preconfigure the second panel to display only that area. In both panels, we color data points by cluster label using the same color map, to help visualize the zoom relationship between the two panels.
Action: Close the app and run the following block of code to launch the example application described above. Double-click on the area selected in the first panel to zoom and display the same view as the second panel. Double click again anywhere in the panel to zoom out.
initial <- list( ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", PointSize = 1, BrushData = list( xmin = 3.4, xmax = 27.3, ymin = 0.8, ymax = 30.0, log = list(x = NULL, y = NULL), mapping = list(x = "X", y = "Y", colour = "ColorBy"), direction = "xy", brushId = "redDimPlot1_Brush", outputId = "redDimPlot1")), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster", PointSize = 2, ZoomData = c(xmin = 3.4, xmax = 27.3, ymin = 0.8, ymax = 30.0)) ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
Sometimes, rectangular selections may not offer enough flexibility to include or exclude certain data points. While RStudio Shiny does not offer freehand/lasso selections yet, iSEE implements this feature through a server-side click-based functionality to draw lasso selections in plot panels.
To start a lasso selection, click once on the plot at the location where you wish to create the first lasso waypoint, and let the plot refresh to show the point. To draw subsequent waypoints, click once at the desired locations, each time allowing the plot to refresh. To close the lasso and activate the selection, click on the starting waypoint of the lasso (a tolerance of 1% of the axis range is applied).
To demonstrate the feature, we define one open (i.e. incomplete) and one closed (i.e., active) lasso selection in two separate Row data plot panels, to select features (i.e., genes) that show a higher total UMI count across all cells than the general trend for other features detected in the same number of cells. In addition, we specify that the the Row data table should display the features captured by the closed lasso selection, along with their associated metadata.
Action: Close the app and run the following block of code to launch the example application described above. Close the incomplete lasso selection in the first panel by clicking on the initial waypoint.
initial <- list( RowDataPlot(XAxis = "Row data", YAxis = "log10_total", XAxisRowData = "detected", BrushData = list(lasso = NULL, closed = FALSE, panelvar1 = NULL, panelvar2 = NULL, mapping = list(x = "X", y = "Y"), coord = structure(c( -0.12, 11.03, 20.69, 39.08, 60.5, 85.53, 93.68, 95.79, 60.8, 34.26, -0.12, 2.29, 2.77, 3.03, 3.32, 3.54, 3.98, 4.24, 4.72, 4.69, 4.43, 3.41), .Dim = c(11L, 2L)))), RowDataPlot(XAxis = "Row data", YAxis = "log10_total", XAxisRowData = "detected", BrushData = list(lasso = NULL, closed = TRUE, panelvar1 = NULL, panelvar2 = NULL, mapping = list(x = "X", y = "Y"), coord = structure(c( -0.12, 11.03, 20.69, 39.08, 60.5, 85.53, 93.68, 95.79, 60.8, 34.26, -0.12, -0.12, 2.29, 2.77, 3.03, 3.32, 3.54, 3.98, 4.24, 4.72, 4.69, 4.43, 3.41, 2.29), .Dim = c(12L, 2L)))), RowDataTable(RowSelectionSource = "RowDataPlot2") ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
Tables may also be preconfigured to immediately display a subset of rows according to their filters or select a row different from the first one when the app is launched.
For simplicity, we demonstrate this feature using a subset of 3 column metadata–namely “Cluster”, “total”, and “subsets_MT_percent”–out of the set of 20 column metadata computed by the workflow above.
In the example below, we preconfigure a Column statistics table to display only cells with cluster labels “3” and “5”, and a library size between 5,000 and 15,844 (the maximum library size in the dataset). In addition, we also preselect the cell named “Cell943”.
Action: Close the app and run the following block of code to launch the example application described above.
sce1 <- sce colData(sce1) <- colData(sce1)[, c("Cluster", "total", "subsets_MT_percent")] initial <- list( ColumnDataTable(SearchColumns = c('["3", "5"]", "5000 ... 15844"', ""), Selected = "Cell943", PanelWidth = 12L) ) app <- iSEE(sce1, initial = initial)
shiny::runApp(app)
Notably the preselected cell can be used to immediately establish a relationship with another panel and highlighted the selection in a different panel as demonstrated earlier.
Action: Close the app and run the following block of code to launch an application that includes a Reduced dimension plot panel dynamically higlighting the cell selected in the Column statistics table. Select different rows of the cell to see the Reduced dimension plot panel update accordingly.
initial <- list( ColumnDataTable(SearchColumns = c('["3", "5"]', "5000 ... 15844", ""), Selected = "Cell943", PanelId = 1L, PanelWidth = 8L), ReducedDimensionPlot(Type = "TSNE", VisualBoxOpen = TRUE, ColorBy = "Sample name", ColorBySampleSource = "ColumnDataTable1") ) app <- iSEE(sce1, initial = initial)
shiny::runApp(app)
In the previous sections we have seen how we can specify the initial state of each panel. Often, the initial exploration of a data set, including reorganization and configuration of panels, is done interactively. Once a satisfactory panel setup has been achieved, the corresponding settings can be exported and used to restore the state of the app at startup.
Action: Click on the export icon () in the top-right corner, and select ‘Display panel settings’. Copy all the code shown in the pop-up window. Close the app and paste this code in your R session. This defines the ‘initial’ list that we have composed manually in the examples above. Then launch a new instance using the following code. Note how the app starts in the same configuration as it was before closing. Of course, it is still possible to continue exploring the data interactively - we have only changed the starting configuration.
iSEE(sce1, initial = initial)
In addition to the standard panel types, iSEE allows the user to create custom panels, displaying any type of plot or table. The custom panels can receive row and column selections from the standard panels. In order to create a custom panel, we need to define a function that generates the desired output. This function must receive a SummarizedExperiment
object as it first input. In addition, it should take a list of character vectors containing row names as the second argument (can be NULL
if no selection is made), and similarly a list of character vectors containing column names as the third argument.
We will first show how to create a custom panel that performs dimension reduction on a subset of the rows and columns, selected in and transmitted from other panels in the interface.
Action: Copy the code below in your R session in order to define the custom function.
library(scater) CUSTOM_DIMRED <- function(se, rows, columns, ntop=500, scale=TRUE, mode=c("PCA", "TSNE", "UMAP")) { if (is.null(columns)) { return( ggplot() + theme_void() + geom_text( aes(x, y, label=label), data.frame(x=0, y=0, label="No column data selected."), size=5) ) } mode <- match.arg(mode) if (mode=="PCA") { calcFUN <- runPCA } else if (mode=="TSNE") { calcFUN <- runTSNE } else if (mode=="UMAP") { calcFUN <- runUMAP } set.seed(1000) kept <- se[, unique(unlist(columns))] kept <- calcFUN(kept, ncomponents=2, ntop=ntop, scale=scale, subset_row=unique(unlist(rows))) plotReducedDim(kept, mode) }
Action: Next, launch an iSEE session with the custom reduced dimension plot panel, as well as a Row data plot panel and a Column data plot panel, from which the row and column selections for the custom reduced dimension representation are obtained. Select subsets of the points in each of these two panels and note how the plot in the custom panel is updated.
GENERATOR <- createCustomPlot(CUSTOM_DIMRED) initial=list( ReducedDimensionPlot( ColorBy = "Column data", ColorByColumnData = "Cluster", PanelId=1L ), RowDataPlot( XAxis="Row data", YAxis="detected", XAxisRowData="log10_total", PanelId=1L), GENERATOR(mode="TSNE", ntop=1000, ColumnSelectionSource="ReducedDimensionPlot1", RowSelectionSource="RowDataPlot1") ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
In the example below, we demonstrate how to transmit multiple selections (active and saved). This information is used to compute the log-fold change for each gene between the active selection and each of the other saved selections.
Action: Copy the code below in your R session in order to define the custom functions.
CUSTOM_DIFFEXP <- function(se, ri, ci, assay="logcounts") { ri <- ri$active if (is.null(ri)) { ri <- rownames(se) } assayMatrix <- assay(se, assay)[ri, , drop=FALSE] if (is.null(ci$active) || length(ci)<2L) { return(data.frame(row.names=character(0), LogFC=integer(0))) # dummy value. } active <- rowMeans(assayMatrix[,ci$active,drop=FALSE]) all_saved <- ci[grep("saved", names(ci))] lfcs <- vector("list", length(all_saved)) for (i in seq_along(lfcs)) { saved <- rowMeans(assayMatrix[,all_saved[[i]],drop=FALSE]) lfcs[[i]] <- active - saved } names(lfcs) <- sprintf("LogFC/%i", seq_along(lfcs)) do.call(data.frame, lfcs) } CUSTOM_HEAT <- function(se, ri, ci, assay="logcounts") { everything <- CUSTOM_DIFFEXP(se, ri, ci, assay=assay) if (nrow(everything) == 0L) { return(ggplot()) # empty ggplot if no genes reported. } everything <- as.matrix(everything) top <- head(order(rowMeans(abs(everything)), decreasing=TRUE), 50) topFC <- everything[top, , drop=FALSE] dfFC <- data.frame( gene=rep(rownames(topFC), ncol(topFC)), contrast=rep(colnames(topFC), each=nrow(topFC)), value=as.vector(topFC) ) ggplot(dfFC, aes(contrast, gene)) + geom_raster(aes(fill = value)) }
Action: Next, launch an iSEE session with the custom table and plot panels, as well as a Reduced dimension plot from which the column selections for the custom panels are obtained. Select and save multiple subsets of the points in the reduced dimension plot (e.g., clusters) and note how the DE table in the custom panel is updated.
TAB_GEN <- createCustomTable(CUSTOM_DIFFEXP) HEAT_GEN <- createCustomPlot(CUSTOM_HEAT) rdp <- ReducedDimensionPlot( PanelId=1L) initial <- list(rdp, TAB_GEN(ColumnSelectionSource="ReducedDimensionPlot1"), HEAT_GEN(ColumnSelectionSource="ReducedDimensionPlot1")) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
More examples of custom panels can be found at https://github.com/iSEE/iSEE_custom and in the Custom panels vignette.
The fact that data exploration is done interactively is no reason to forego reproducibility! To this end, iSEE lets you export the exact R code used to create each of the currently visible plots.
Action: Close the application and run the following block of code to set up 3 panels, all displaying the t-SNE representation.
initial <- list( ReducedDimensionPlot(Type = "TSNE", ColorBy = "Column data", ColorByColumnData = "Cluster"), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Feature name", ColorByFeatureName = "CD79A", ColorByFeatureNameAssay = "logcounts"), ReducedDimensionPlot(Type = "TSNE", ColorBy = "Feature name", ColorByFeatureName = "CD74", ColorByFeatureNameAssay = "logcounts") ) app <- iSEE(sce, initial = initial)
shiny::runApp(app)
Then click on the export icon () in the top-right corner, and select ‘Extract the R code’. Copy the preamble and the code for the first Reduced dimension plot from the pop-up window. Close the app and paste the code in your R session. Note how this recreates the first reduced dimension plot from the app.
## The following list of commands will generate the plots created in iSEE ## Copy them into a script or an R session containing your SingleCellExperiment. ## All commands below refer to your SingleCellExperiment object as `se`. se <- sce colormap <- ExperimentColorMap() colormap <- synchronizeAssays(colormap, se) all_contents <- list() ################################################################################ # Defining brushes ################################################################################ all_active <- list() all_active[['ReducedDimensionPlot1']] <- list() all_active[['ReducedDimensionPlot2']] <- list() all_active[['ReducedDimensionPlot3']] <- list() ################################################################################ ## Reduced dimension plot 1 ################################################################################ red.dim <- reducedDim(se, "TSNE"); plot.data <- data.frame(X=red.dim[, 1], Y=red.dim[, 2], row.names=colnames(se)); plot.data$ColorBy <- colData(se)[, "Cluster"]; # Avoid visual biases from default ordering by shuffling the points set.seed(2643); plot.data <- plot.data[sample(nrow(plot.data)),,drop=FALSE]; dot.plot <- ggplot() + geom_point(aes(x=X, y=Y, color=ColorBy), alpha=1, plot.data, size=1) + labs(x="Dimension 1", y="Dimension 2", color="Cluster", title="TSNE") + coord_cartesian(xlim=range(plot.data$X, na.rm=TRUE), ylim=range(plot.data$Y, na.rm=TRUE), expand=TRUE) + scale_color_manual(values=colDataColorMap(colormap, "Cluster", discrete=TRUE)(14), na.value='grey50', drop=FALSE) + scale_fill_manual(values=colDataColorMap(colormap, "Cluster", discrete=TRUE)(14), na.value='grey50', drop=FALSE) + theme_bw() + theme(legend.position='bottom', legend.box='vertical', legend.text=element_text(size=9), legend.title=element_text(size=11), axis.text=element_text(size=10), axis.title=element_text(size=12), title=element_text(size=12)) # Saving data for transmission all_contents[['ReducedDimensionPlot1']] <- plot.data
dot.plot
The iSEEu (“iSEE universe”) Bioconductor package defines additional custom panels and predefined ‘modes’ (startup configurations) that may be useful for specific applications. Here we illustrate the use of the reduced dimension mode, which will start an application with one reduced dimension panel for each reduced dimension representation in the input object.
library(iSEEu) app <- modeReducedDim(sce)
shiny::runApp(app)
Furthermore, since iSEE version 2.0.0, users can leverage the implementation of panels as a hierarchy of S4 classes to rapidly extend the framework and develop new types of panels with virtually unlimited freedom of functionality. In this framework, new panel types can be readily integrated in apps alongside built-in panel types and immediately benefit of most panel functionality with minimal effort, including the capacity to transmit selections to and from other panels. New panels may be implemented by extending the core virtual class Panel
directly for full control over the user interface and reactive observers; alternatively, it is often desirable for developers to inherit from one of the concrete panel classes to get most of the essential functionality from the parent class “for free”.
Use cases demonstrating the implementation of new panels with various levels of complexity are available in the Extending iSEE bookdown. For example, the Hexagonal reduced dimension plot - implemented in the iSEEu package - demonstrates an alternative to the downsampling strategy, by summarizing data points into hexagonal bins.
app <- iSEE(sce, initial = list( ReducedDimensionPlot(PanelWidth = 6L), ReducedDimensionHexPlot(PanelWidth = 6L) ))
shiny::runApp(app)
New panels and modes may be defined in independent R packages that declare Imports: iSEE
in their DESCRIPTION file, and distributed online to create an “iSEE-verse” of interoperable packages for interactive visualization, akin to the tidyverse collection of R packages for data science. If you have ideas about additional panels and/or modes that might be useful, we welcome issues and pull requests at https://github.com/iSEE/iSEEu.
iSEE
bookdown book: https://isee.github.io/iSEE-book/
sessionInfo() #> R version 4.0.0 (2020-04-24) #> Platform: x86_64-pc-linux-gnu (64-bit) #> Running under: Ubuntu 18.04.4 LTS #> #> Matrix products: default #> BLAS/LAPACK: /usr/lib/x86_64-linux-gnu/libopenblasp-r0.2.20.so #> #> locale: #> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C #> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8 #> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=C #> [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] parallel stats4 stats graphics grDevices utils datasets #> [8] methods base #> #> other attached packages: #> [1] iSEEu_1.1.3 iSEE_2.1.10 #> [3] BiocSingular_1.5.0 scran_1.17.10 #> [5] scater_1.17.4 ggplot2_3.3.2 #> [7] TENxPBMCData_1.7.0 HDF5Array_1.17.3 #> [9] rhdf5_2.33.7 SingleCellExperiment_1.11.6 #> [11] SummarizedExperiment_1.19.6 DelayedArray_0.15.7 #> [13] matrixStats_0.56.0 Matrix_1.2-18 #> [15] Biobase_2.49.0 GenomicRanges_1.41.5 #> [17] GenomeInfoDb_1.25.8 IRanges_2.23.10 #> [19] S4Vectors_0.27.12 BiocGenerics_0.35.4 #> [21] BiocStyle_2.17.0 #> #> loaded via a namespace (and not attached): #> [1] backports_1.1.8 circlize_0.4.10 #> [3] AnnotationHub_2.21.1 BiocFileCache_1.13.0 #> [5] igraph_1.2.5 shinydashboard_0.7.1 #> [7] splines_4.0.0 BiocParallel_1.23.2 #> [9] digest_0.6.25 htmltools_0.5.0 #> [11] viridis_0.5.1 magrittr_1.5 #> [13] memoise_1.1.0 cluster_2.1.0 #> [15] limma_3.45.9 ComplexHeatmap_2.5.3 #> [17] pkgdown_1.5.1 colorspace_1.4-1 #> [19] blob_1.2.1 rappdirs_0.3.1 #> [21] xfun_0.15 dplyr_1.0.0 #> [23] crayon_1.3.4 RCurl_1.98-1.2 #> [25] jsonlite_1.7.0 glue_1.4.1 #> [27] gtable_0.3.0 zlibbioc_1.35.0 #> [29] XVector_0.29.3 GetoptLong_1.0.2 #> [31] Rhdf5lib_1.11.3 shape_1.4.4 #> [33] scales_1.1.1 DBI_1.1.0 #> [35] edgeR_3.31.4 miniUI_0.1.1.1 #> [37] Rcpp_1.0.5 viridisLite_0.3.0 #> [39] xtable_1.8-4 clue_0.3-57 #> [41] dqrng_0.2.1 bit_1.1-15.2 #> [43] rsvd_1.0.3 DT_0.14 #> [45] htmlwidgets_1.5.1 httr_1.4.2 #> [47] FNN_1.1.3 RColorBrewer_1.1-2 #> [49] shinyAce_0.4.1 ellipsis_0.3.1 #> [51] farver_2.0.3 pkgconfig_2.0.3 #> [53] scuttle_0.99.11 uwot_0.1.8 #> [55] dbplyr_1.4.4 locfit_1.5-9.4 #> [57] tidyselect_1.1.0 labeling_0.3 #> [59] rlang_0.4.7 later_1.1.0.1 #> [61] AnnotationDbi_1.51.1 munsell_0.5.0 #> [63] BiocVersion_3.12.0 tools_4.0.0 #> [65] generics_0.0.2 RSQLite_2.2.0 #> [67] ExperimentHub_1.15.0 rintrojs_0.2.2 #> [69] evaluate_0.14 stringr_1.4.0 #> [71] fastmap_1.0.1 yaml_2.2.1 #> [73] knitr_1.29 bit64_0.9-7.1 #> [75] fs_1.4.2 purrr_0.3.4 #> [77] nlme_3.1-148 mime_0.9 #> [79] compiler_4.0.0 beeswarm_0.2.3 #> [81] curl_4.3 png_0.1-7 #> [83] interactiveDisplayBase_1.27.5 tibble_3.0.3 #> [85] statmod_1.4.34 stringi_1.4.6 #> [87] desc_1.2.0 RSpectra_0.16-0 #> [89] lattice_0.20-41 shinyjs_1.1 #> [91] vctrs_0.3.2 pillar_1.4.6 #> [93] lifecycle_0.2.0 rhdf5filters_1.1.2 #> [95] BiocManager_1.30.10 GlobalOptions_0.1.2 #> [97] BiocNeighbors_1.7.0 bitops_1.0-6 #> [99] irlba_2.3.3 httpuv_1.5.4 #> [101] R6_2.4.1 promises_1.1.1 #> [103] gridExtra_2.3 vipor_0.4.5 #> [105] colourpicker_1.0 MASS_7.3-51.6 #> [107] assertthat_0.2.1 rprojroot_1.3-2 #> [109] rjson_0.2.20 shinyWidgets_0.5.3 #> [111] withr_2.2.0 GenomeInfoDbData_1.2.3 #> [113] mgcv_1.8-31 grid_4.0.0 #> [115] beachmat_2.5.0 rmarkdown_2.3 #> [117] DelayedMatrixStats_1.11.1 Rtsne_0.15 #> [119] shiny_1.5.0 ggbeeswarm_0.6.0
Bremner, Frederick J, Stephen J Gotts, and Dina L Denham. 1994. “Hinton Diagrams: Viewing Connection Strengths in Neural Networks.” Behav. Res. Methods Instrum. Comput. 26 (2): 215–18.
Hinton, Geoffrey E, and Tim Shallice. 1991. “Lesioning an Attractor Network: Investigations of Acquired Dyslexia.” Psychological Review 98 (1): 74–95.
Huber, Wolfgang, Vincent J Carey, Robert Gentleman, Simon Anders, Marc Carlson, Benilton S Carvalho, Hector Corrada Bravo, et al. 2015. “Orchestrating High-Throughput Genomic Analysis with Bioconductor.” Nat. Methods 12 (2): 115–21.
Rue-Albrecht, Kevin, Federico Marini, Charlotte Soneson, and Aaron T L Lun. 2018. “iSEE: Interactive SummarizedExperiment Explorer.” F1000Res. 7 (June).
University of Oxford, Oxford, UK↩︎
Friedrich Miescher Institute for Biomedical Research and SIB Swiss Institute of Bioinformatics↩︎
Institute of Medical Biostatistics, Epidemiology and Informatics (IMBEI), Mainz and Center for Thrombosis and Hemostasis (CTH), Mainz↩︎
Genentech, Bioinformatics and Computational Biology↩︎