Title: | Single-Cell Data Visualisation Using Kernel Gene-Weighted Density Estimation |
---|---|
Description: | This package provides a enhanced visualization of single-cell data based on gene-weighted density estimation. Nebulosa recovers the signal from dropped-out features and allows the inspection of the joint expression from multiple features (e.g. genes). Seurat and SingleCellExperiment objects can be used within Nebulosa. |
Authors: | Jose Alquicira-Hernandez [aut, cre] |
Maintainer: | Jose Alquicira-Hernandez <[email protected]> |
License: | GPL-3 |
Version: | 1.17.0 |
Built: | 2024-10-30 08:32:19 UTC |
Source: | https://github.com/bioc/Nebulosa |
Estimate weighted kernel density
calculate_density(w, x, method, adjust = 1, map = TRUE)
calculate_density(w, x, method, adjust = 1, map = TRUE)
w |
Vector with weights for each observation |
x |
Matrix with dimensions where to calculate the density from. Only the first two dimensions will be used |
method |
Kernel density estimation method:
|
adjust |
Numeric value to adjust to bandwidth. Default: 1. Not available
for |
map |
Whether to map densities to individual observations |
If map
is TRUE
, a vector with corresponding densities
for each observation is returned. Otherwise,
a list with the density estimates from the selected method is returned.
Jose Alquicira-Hernandez
dens <- Nebulosa:::calculate_density(iris[, 3], iris[, 1:2], method = "wkde")
dens <- Nebulosa:::calculate_density(iris[, 3], iris[, 1:2], method = "wkde")
Plot gene-weighted 2D kernel density
plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... ) ## S4 method for signature 'Seurat' plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... ) ## S4 method for signature 'SingleCellExperiment' plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... )
plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... ) ## S4 method for signature 'Seurat' plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... ) ## S4 method for signature 'SingleCellExperiment' plot_density( object, features, slot = NULL, joint = FALSE, reduction = NULL, dims = c(1, 2), method = c("ks", "wkde"), adjust = 1, size = 1, shape = 16, combine = TRUE, pal = "viridis", raster = TRUE, ... )
object |
Seurat or SingleCellExperiment object |
features |
Features (e.g. genes) to visualize |
slot |
Type of data: |
joint |
Return joint density plot? By default |
reduction |
Name of the reduction to visualize. If not provided, last computed reduction is visualized |
dims |
Vector of length 2 specifying the dimensions to be plotted. By default, the first two dimensions are considered. |
method |
Kernel density estimation method:
|
adjust |
Numeric value to adjust to bandwidth. Default: 1. Not available
for |
size |
Size of the geom to be plotted (e.g. point size) |
shape |
Shape of the geom to be plotted |
combine |
Create a single plot? If |
pal |
String specifying the viridis color palette to use. |
raster |
Rasterise plot |
... |
Further scale arguments passed to scale_color_viridis_c Options:
|
A scatterplot from a given reduction showing the gene-weighted density
plot_density(Seurat)
: Plot gene-weighted 2D kernel density
plot_density(SingleCellExperiment)
: Plot gene-weighted 2D kernel density
Jose Alquicira-Hernandez
data <- SeuratObject::pbmc_small plot_density(data, "CD3E")
data <- SeuratObject::pbmc_small plot_density(data, "CD3E")
Plot density estimates
plot_density_( z, feature, cell_embeddings, dim_names, shape, size, legend_title, pal = c("viridis", "magma", "cividis", "inferno", "plasma"), raster, ... )
plot_density_( z, feature, cell_embeddings, dim_names, shape, size, legend_title, pal = c("viridis", "magma", "cividis", "inferno", "plasma"), raster, ... )
z |
Vector with density values for each cells |
feature |
Name of the feature being plotted |
cell_embeddings |
Matrix with cell embeddings |
dim_names |
Names of the dimensions from the cell embeddings |
shape |
Geom shape |
size |
Geom size |
legend_title |
String used as legend title |
pal |
String specifying the viridis color palette to use |
raster |
Rasterise plot |
... |
Further scale arguments passed to scale_color_viridis_c |
A ggplot object
Jose Alquicira-Hernandez
Weighted 2D kernel density estimation
wkde2d(x, y, w, h, adjust = 1, n = 100, lims = c(range(x), range(y)))
wkde2d(x, y, w, h, adjust = 1, n = 100, lims = c(range(x), range(y)))
x |
Dimension 1 |
y |
Dimension 2 |
w |
Weight variable |
h |
vector of bandwidths for x and y directions. Defaults to normal reference bandwidth (ks::hpi). A scalar value will be taken to apply to both directions. |
adjust |
Bandwidth adjustment |
n |
Number of grid points in each direction. Can be scalar or a length-2 integer vector. |
lims |
The limits of the rectangle covered by the grid as c(xl, xu, yl, yu). |
A list of three components.
x, y
The x and y coordinates of the grid points, vectors of
length n.
z
An n[1] by n[2] matrix of the weighted estimated density:
rows correspond to the value of x, columns to the value of y.
Jose Alquicira-Hernandez
set.seed(1) x <- rnorm(100) set.seed(2) y <- rnorm(100) set.seed(3) w <- sample(c(0, 1), 100, replace = TRUE) dens <- Nebulosa:::wkde2d(x, y, w)
set.seed(1) x <- rnorm(100) set.seed(2) y <- rnorm(100) set.seed(3) w <- sample(c(0, 1), 100, replace = TRUE) dens <- Nebulosa:::wkde2d(x, y, w)