Title: | Streamlining the creation of initial states for starting an iSEE instance |
---|---|
Description: | iSEEfier provides a set of functionality to quickly and intuitively create, inspect, and combine initial configuration objects. These can be conveniently passed in a straightforward manner to the function call to launch iSEE() with the specified configuration. This package currently works seamlessly with the sets of panels provided by the iSEE and iSEEu packages, but can be extended to accommodate the usage of any custom panel (e.g. from iSEEde, iSEEpathways, or any panel developed independently by the user). |
Authors: | Najla Abassi [aut, cre] , Federico Marini [aut] |
Maintainer: | Najla Abassi <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.3.0 |
Built: | 2024-10-31 06:17:50 UTC |
Source: | https://github.com/bioc/iSEEfier |
Constant values used throughout iSEEfier
iSEE_panel_colors
iSEE_panel_colors
An object of class character
of length 17.
color values (as string character or hex value) for the panels included by
default in iSEE
and iSEEu
Glue a set of initial
configuration objects, combining them into a single
valid initial
set.
glue_initials( ..., remove_duplicate_panels = TRUE, verbose = TRUE, custom_panels_allowed = NULL )
glue_initials( ..., remove_duplicate_panels = TRUE, verbose = TRUE, custom_panels_allowed = NULL )
... |
A set of |
remove_duplicate_panels |
Logical, defaults to |
verbose |
Logical, defaults to |
custom_panels_allowed |
Character vector, defaults to |
The usage of custom_panels_allowed
can be especially relevant when one creates
one or more custom panels, with a specific name that needs to be indicated in
this parameter.
For example, if using a panel of class FancyPlotPanel
and one called
FancyTablePanel
, the value for custom_panels_allowed
should be set to
c("FancyPlotPanel", "FancyTablePanel")
.
It is worth mentioning that iSEE::iSEE()
is actually able to handle the
automatic renaming of panels that could be detected as duplicated. This can
basically relax the requirement on the "uniqueness" of the configured panels, with
the only caveat of having to think of how the transmissions between panels
will be handled; nevertheless, most users might not even need to face this
situation.
A single initial
list object, in the format that is required to
be passed as a parameter in the call to iSEE::iSEE()
, concatenating the
values provided as input.
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list_1 <- c("ENSMUSG00000026581") gene_list_2 <- c("ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial1 <- iSEEinit(sce = sce, features = gene_list_1, clusters = cluster, groups = group) initial2 <- iSEEinit(sce = sce, features = gene_list_2, clusters = cluster, groups = group) initials_merged <- glue_initials(initial1, initial2) view_initial_tiles(initial1) view_initial_tiles(initial2) view_initial_tiles(initials_merged) ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial_merged)
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list_1 <- c("ENSMUSG00000026581") gene_list_2 <- c("ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial1 <- iSEEinit(sce = sce, features = gene_list_1, clusters = cluster, groups = group) initial2 <- iSEEinit(sce = sce, features = gene_list_2, clusters = cluster, groups = group) initials_merged <- glue_initials(initial1, initial2) view_initial_tiles(initial1) view_initial_tiles(initial2) view_initial_tiles(initials_merged) ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial_merged)
iSEEinit()
defines the initial setup of an iSEE instance, recommending
tailored visual elements to effortlessly illustrate the expression of a gene
list in a single view.
iSEEinit( sce, features, reddim_type = "TSNE", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], gene_id = "id", add_markdown_panel = FALSE )
iSEEinit( sce, features, reddim_type = "TSNE", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], gene_id = "id", add_markdown_panel = FALSE )
sce |
SingleCellExperiment object |
features |
A character vector or a data.frame containing a list of genes.
If |
reddim_type |
A string vector containing the dimensionality reduction type |
clusters |
A character string containing the name of the clusters/cell-type/state...(as listed in the colData of the sce) |
groups |
A character string of the groups/conditions...(as it appears in the colData of the sce) |
gene_id |
A character string containing the name of the column name containing gene names/ids, when 'features' is a data.frame |
add_markdown_panel |
A logical indicating whether or not to include the MarkdownBoard panel in the initial configuration |
A list of "Panel" objects specifying the initial state of iSEE instance
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group)
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group)
iSEEmarker()
creates an initial state of an iSEE instance for interactive
exploration of marker genes through the DynamicMarkerTable
panel, synchronizing
selections with a ReducedDimensionPlot
and a FeatureAssayPlot
to visualize
the expression of selected marker genes
iSEEmarker( sce, reddim_type = "TSNE", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], selection_plot_format = c("ColumnDataPlot", "ReducedDimensionPlot") )
iSEEmarker( sce, reddim_type = "TSNE", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], selection_plot_format = c("ColumnDataPlot", "ReducedDimensionPlot") )
sce |
SingleCellExperiment object |
reddim_type |
A string vector containing the dimensionality reduction |
clusters |
A character string containing the name of the clusters/cell-type/state...(as listed in the colData of the sce) |
groups |
A character string of the groups/conditions...(as it appears in the colData of the sce) |
selection_plot_format |
A string character containing the class of the panel.
It can be either |
A list of "Panel" objects specifying the initial state of iSEE instance
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) cluster <- "stimulus" group <- "single cell quality" initial <- iSEEmarker(sce = sce, clusters = cluster, groups = group, selection_plot_format = "ColumnDataPlot")
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) cluster <- "stimulus" group <- "single cell quality" initial <- iSEEmarker(sce = sce, clusters = cluster, groups = group, selection_plot_format = "ColumnDataPlot")
iSEEnrich()
creates an initial state of an iSEE instance for interactive
exploration of feature sets extracted from GO and KEGG database, displaying
all associated genes in a RowDataTable
panel.
iSEEnrich( sce, collection = c("GO", "KEGG"), organism = "org.Hs.eg.db", gene_identifier = "ENTREZID", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], reddim_type = "PCA" )
iSEEnrich( sce, collection = c("GO", "KEGG"), organism = "org.Hs.eg.db", gene_identifier = "ENTREZID", clusters = colnames(colData(sce))[1], groups = colnames(colData(sce))[1], reddim_type = "PCA" )
sce |
SingleCellExperiment object |
collection |
A character vector specifying the gene set collections of interest (GO,KEGG) |
organism |
A character string of the org.*.eg.db package to use to extract mappings of gene sets to gene IDs. |
gene_identifier |
A character string specifying the identifier to use to extract gene IDs for the organism package |
clusters |
A character string containing the name of the clusters/cell-type/state...(as listed in the colData of the sce) |
groups |
A character string of the groups/conditions...(as it appears in the colData of the sce) |
reddim_type |
A string vector containing the dimensionality reduction type |
A list of "Panel" objects specifying the initial state of iSEE instance
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) GO_collection <- "GO" Mm_organism <- "org.Mm.eg.db" gene_id <- "SYMBOL" clusters <- "stimulus" groups <- "single cell quality" reddim_type <- "PCA" results <- iSEEnrich(sce = sce, collection = GO_collection, organism = Mm_organism, gene_identifier = gene_id, clusters = clusters, groups = groups, reddim_type = reddim_type)
sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) GO_collection <- "GO" Mm_organism <- "org.Mm.eg.db" gene_id <- "SYMBOL" clusters <- "stimulus" groups <- "single cell quality" reddim_type <- "PCA" results <- iSEEnrich(sce = sce, collection = GO_collection, organism = Mm_organism, gene_identifier = gene_id, clusters = clusters, groups = groups, reddim_type = reddim_type)
Translates the layout of the initial
configuration object as a networks,
representing panels as nodes and links between them as edges.
view_initial_network(initial, plot_format = c("igraph", "visNetwork", "none"))
view_initial_network(initial, plot_format = c("igraph", "visNetwork", "none"))
initial |
An |
plot_format |
Character string, one of |
Panels are the nodes, with color and names to identify them easily. The connections among panels are represented through directed edges. This can be a compact visualization to obtain an overview for the configuration, without the need of fully launching the app and loading the content of all panels
This function is particularly useful with mid-to-large initial
objects, as
they can be quickly generated in a programmatic manner via the iSEEinit()
provided in this package.
An igraph
object, underlying the visual representation provided.
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group) g_init <- view_initial_network(initial) g_init view_initial_network(initial, plot_format = "visNetwork") ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial)
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group) g_init <- view_initial_network(initial) g_init view_initial_network(initial, plot_format = "visNetwork") ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial)
Previews the layout of the initial
configuration object in a graphical form.
view_initial_tiles(initial)
view_initial_tiles(initial)
initial |
An |
Tiles are used to represent the panel types, and reflect the values of their width. This can be a compact visualization to obtain an overview for the configuration, without the need of fully launching the app and loading the content of all panels
This function is particularly useful with mid-to-large initial
objects, as
they can be quickly generated in a programmatic manner via the iSEEinit()
provided in this package.
A ggplot
object, representing a schematic view for the initial
object.
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group) view_initial_tiles (initial) ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial)
## Load a dataset and preprocess this quickly sce <- scRNAseq::RichardTCellData() sce <- scuttle::logNormCounts(sce) sce <- scater::runPCA(sce) sce <- scater::runTSNE(sce) ## Select some features and aspects to focus on gene_list <- c("ENSMUSG00000026581", "ENSMUSG00000005087", "ENSMUSG00000015437") cluster <- "stimulus" group <- "single cell quality" initial <- iSEEinit(sce = sce, features = gene_list, clusters = cluster, groups = group) view_initial_tiles (initial) ## Continue your exploration directly within iSEE! if (interactive()) iSEE(sce, initial = initial)