Title: | iSEE extension for summarising data points in hexagonal bins |
---|---|
Description: | This package provides panels summarising data points in hexagonal bins for `iSEE`. It is part of `iSEEu`, the iSEE universe of panels that extend the `iSEE` package. |
Authors: | Kevin Rue-Albrecht [aut, cre] , Charlotte Soneson [aut] , Federico Marini [aut] , Aaron Lun [aut] |
Maintainer: | Kevin Rue-Albrecht <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.9.0 |
Built: | 2024-11-29 08:16:06 UTC |
Source: | https://github.com/bioc/iSEEhex |
The ReducedDimensionHexPlot is a ReducedDimensionPlot subclass that is dedicated to creating a reduced dimension plot summarising data points in hexagonal bins.
The following slots control the parameters used in the visualization:
BinResolution
, a numeric positive scalar specifying the number of hexagonal bins in both vertical and horizontal directions.
Defaults to 100.
In addition, this class inherits all slots from its parent ReducedDimensionPlot, ColumnDotPlot, DotPlot and Panel classes.
ReducedDimensionHexPlot(...)
creates an instance of a ReducedDimensionHexPlot class,
where any slot and its value can be passed to ...
as a named argument.
In the following code snippets, x
is an instance of a ReducedDimensionHexPlot class.
Refer to the documentation for each method for more details on the remaining arguments.
For defining the interface:
.panelColor(x)
will return the specified default color for this panel class.
.fullName(x)
will return "Hexagonal reduced dimension plot"
.
.hideInterface(x, field)
will return TRUE
for field="Downsample"
as downsampling is not applicable to this panel that summarizes all data points in each hexagonal bin;
otherwise this function will call the ReducedDimensionPlot method.
.defineVisualShapeInterface(x)
will return NULL
for this panel, as the shape aesthetic is not applicable to this panel that does not display individual data points.
.defineVisualSizeInterface(x)
overrides the equivalent method inherited from all parents classes and will return instead an HTML tag definition that contains a user input controlling the number of hexagonal bins in both vertical and horizontal directions.
.defineVisualOtherInterface(x)
will return NULL
, as there are no additional visual parameters for this panel.
.allowableColorByDataChoices(x, se)
will return a character vector with the names of all continuous fields in colData(se)
,
where se
is the input SummarizedExperiment object.
For monitoring reactive expressions:
.createObservers(x, se, input, session, pObjects, rObjects)
sets up observers for all new slots described above, as well as in the parent classes via the ReducedDimensionPlot method.
For creating the plot:
.generateDotPlot(x, envir)
will return a list with plot
, a ggplot2::ggplot()
object; and commands
, a character vector of commands to produce that object when evaluated inside envir
.
For documentation:
.definePanelTour(x)
returns an data.frame containing the steps of a panel-specific tour.
.getDotPlotColorHelp(x, color_choices)
returns a function that generates an rintrojs tour for the color choice UI.
Kevin Rue-Albrecht
ReducedDimensionPlot, for the base class.
library(scRNAseq) # Example data ---- sce <- ReprocessedAllenData(assays="tophat_counts") class(sce) library(scater) sce <- logNormCounts(sce, exprs_values="tophat_counts") sce <- runPCA(sce, ncomponents=4) sce <- runTSNE(sce) rowData(sce)$ave_count <- rowMeans(assay(sce, "tophat_counts")) rowData(sce)$n_cells <- rowSums(assay(sce, "tophat_counts") > 0) # launch the app itself ---- if (interactive()) { iSEE(sce, initial=list( ReducedDimensionHexPlot(BinResolution=50), ReducedDimensionPlot() )) }
library(scRNAseq) # Example data ---- sce <- ReprocessedAllenData(assays="tophat_counts") class(sce) library(scater) sce <- logNormCounts(sce, exprs_values="tophat_counts") sce <- runPCA(sce, ncomponents=4) sce <- runTSNE(sce) rowData(sce)$ave_count <- rowMeans(assay(sce, "tophat_counts")) rowData(sce)$n_cells <- rowSums(assay(sce, "tophat_counts") > 0) # launch the app itself ---- if (interactive()) { iSEE(sce, initial=list( ReducedDimensionHexPlot(BinResolution=50), ReducedDimensionPlot() )) }