Title: | Visualizing Single Cell and Spatial Transcriptomics |
---|---|
Description: | Useful functions to visualize single cell and spatial data. It supports visualizing 'Seurat', 'SingleCellExperiment' and 'SpatialExperiment' objects through grammar of graphics syntax implemented in 'ggplot2'. |
Authors: | Guangchuang Yu [aut, cre, cph] , Shuangbin Xu [aut] , Noriaki Sato [ctb] |
Maintainer: | Guangchuang Yu <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.5.0 |
Built: | 2024-11-18 03:59:43 UTC |
Source: | https://github.com/bioc/ggsc |
Two-Dimensional Weighted Kernel Density Estimation And Mapping the Result To Original Dimension
CalWkdeCpp(x, w, l, h, adjust = 1, n = 400L)
CalWkdeCpp(x, w, l, h, adjust = 1, n = 400L)
x |
The 2-D coordinate matrix |
w |
The weighted sparse matrix, the number columns the same than the number rows than x. |
l |
The limits of the rectangle covered by the grid as c(xl, xu, yl, yu) |
h |
The vector of bandwidths for x and y directions, defaults to normal reference bandwidth (see bandwidth.nrd), A scalar value will be taken to apply to both directions (see ks::hpi). |
adjust |
numeric value to adjust to bandwidth, default is 1. |
n |
number of grid points in the two directions, default is 400. |
Each Geom has an associated function that draws the key when the geom needs to be displayed in a legend. These are the options built into ggplot2.
draw_key_bgpoint(data, params, size)
draw_key_bgpoint(data, params, size)
data |
A single row data frame containing the scaled aesthetics to display in this key |
params |
A list of additional parameters supplied to the geom. |
size |
Width and height of key in mm. |
A grid grob.
this add the background color for geom_point
geom_bgpoint( mapping = NULL, data = NULL, stat = "identity", position = "identity", ..., na.rm = FALSE, show.legend = NA, inherit.aes = TRUE, gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1, pointsize = NULL )
geom_bgpoint( mapping = NULL, data = NULL, stat = "identity", position = "identity", ..., na.rm = FALSE, show.legend = NA, inherit.aes = TRUE, gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1, pointsize = NULL )
mapping |
Set of aesthetic mappings created by |
data |
The data to be displayed in this layer. There are three options: If A A |
stat |
The statistical transformation to use on the data for this layer.
When using a
|
position |
A position adjustment to use on the data for this layer. This
can be used in various ways, including to prevent overplotting and
improving the display. The
|
... |
Other arguments passed on to |
na.rm |
If |
show.legend |
logical. Should this layer be included in the legends?
|
inherit.aes |
If |
gap_colour |
colour of gap background between the bottom background
and top point point layer, default is |
gap_alpha |
numeric the transparency of gap background colour, default is 1. |
bg_line_width |
numeric the line width of background point layer,
default is |
gap_line_width |
numeric the line width of gap between the background and
top point point layer, default is |
pointsize |
numeric the size of point, default is NULL, will use the
internal size aesthetics of |
colour
the colour of point, default is black
.
bg_colour
the colour of background point, default is NA
.
alpha
the transparency of colour, default is 1.
subset
subset the data frame which meet conditions to display.
geom_bgpoint()
understands the following aesthetics (required aesthetics are in bold):
Learn more about setting these aesthetics in vignette("ggplot2-specs")
.
Shuangbin Xu
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species, bg_colour=Species) ) + geom_bgpoint(pointsize=4, gap_line_width = .1, bg_line_width = .3)
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species, bg_colour=Species) ) + geom_bgpoint(pointsize=4, gap_line_width = .1, bg_line_width = .3)
this add the background colour for the geom_scattermore
geom_scattermore2( mapping = NULL, data = NULL, stat = "identity", position = "identity", ..., na.rm = FALSE, show.legend = NA, inherit.aes = TRUE, interpolate = FALSE, pointsize = 0, pixels = c(512, 512), gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1 )
geom_scattermore2( mapping = NULL, data = NULL, stat = "identity", position = "identity", ..., na.rm = FALSE, show.legend = NA, inherit.aes = TRUE, interpolate = FALSE, pointsize = 0, pixels = c(512, 512), gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1 )
mapping |
Set of aesthetic mappings created by |
data |
The data to be displayed in this layer. There are three options: If A A |
stat |
The statistical transformation to use on the data for this layer.
When using a
|
position |
A position adjustment to use on the data for this layer. This
can be used in various ways, including to prevent overplotting and
improving the display. The
|
... |
Other arguments passed on to |
na.rm |
If |
show.legend |
logical. Should this layer be included in the legends?
|
inherit.aes |
If |
interpolate |
A logical value indicating whether to linearly interpolate
the image (the alternative is to use nearest-neighbour interpolation,
which gives a more blocky result). Default |
pointsize |
Radius of rasterized point. Use ‘0’ for single pixels (fastest). |
pixels |
Vector with X and Y resolution of the raster, default |
gap_colour |
colour of gap background between the bottom background
and top point point layer, default is |
gap_alpha |
numeric the transparency of gap background colour, default is 1. |
bg_line_width |
numeric the line width of background point layer,
default is |
gap_line_width |
numeric the line width of gap between the background and
top point point layer, default is |
colour
the colour of point, default is black
.
bg_colour
the colour of background point, default is NA
.
alpha
the transparency of colour, default is 1.
subset
subset the data frame which meet conditions to display.
polygonal point layer
geom_scattermore2()
understands the following aesthetics (required aesthetics are in bold):
Learn more about setting these aesthetics in vignette("ggplot2-specs")
.
Shuangbin Xu
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species, bg_colour=Species) ) + geom_scattermore2(pointsize=4, gap_line_width = .1, bg_line_width = .3)
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species, bg_colour=Species) ) + geom_scattermore2(pointsize=4, gap_line_width = .1, bg_line_width = .3)
plot_lisa_feature
plot_lisa_feature( spe, lisa.res, features = NULL, assay.type = "logcounts", geom = geom_bgpoint, pointsize = 2, hlpointsize = 1.8, clustertype = "High", hlcolor = c("black"), gap_line_width = 0.1, bg_line_width = 0.3, facet_name = NULL, reduction = NULL, image.plot = FALSE, label_wrap_width = 30, ... )
plot_lisa_feature( spe, lisa.res, features = NULL, assay.type = "logcounts", geom = geom_bgpoint, pointsize = 2, hlpointsize = 1.8, clustertype = "High", hlcolor = c("black"), gap_line_width = 0.1, bg_line_width = 0.3, facet_name = NULL, reduction = NULL, image.plot = FALSE, label_wrap_width = 30, ... )
spe |
SpatialExperiment or SingleCellExperiment object. |
lisa.res |
the result returned by |
features |
selected features to be visualized, default is NULL. |
assay.type |
the assay name where data will be used from
(e.g., 'data', 'counts'), default is |
geom |
the function of geometric layer, default is |
pointsize |
numeric the size of point, default is |
hlpointsize |
numeric the size of point which contains corresbonding
spatially variable gene(i.e., SVG), default is |
clustertype |
cell type which is from the result of |
hlcolor |
the color of circular line which enfolds the point
that contains SVG, default is |
gap_line_width |
numeric the line width of gap between the background and
top point point layer, default is |
bg_line_width |
numeric the line width of background point layer,
default is |
facet_name |
the name of facet used in |
reduction |
reduction method, default is |
image.plot |
logical whether display the image of spatial experiment, default is FALSE. |
label_wrap_width |
numeric maximum number of characters before wrapping the strip.
default is |
... |
additional parameters pass to
|
ggplot object
## Not run: library(ggplot2) library(SingleCellExperiment) |> suppressPackageStartupMessages() library(SpatialExperiment) |> suppressPackageStartupMessages() library(STexampleData) # create ExperimentHub instance eh <- ExperimentHub() # query STexampleData datasets myfiles <- query(eh, "STexampleData") ah_id <- myfiles$ah_id[myfiles$title == 'Visium_humanDLPFC'] spe <- myfiles[[ah_id]] spe <- spe[, colData(spe)$in_tissue == 1] spe <-scater::logNormCounts(spe) genes <- c('MOBP', 'PCP4', 'SNAP25', 'HBB', 'IGKC', 'NPY') target.features <- rownames(spe)[match(genes, rowData(spe)$gene_name)] library(SVP) lisa.res1 <- runLISA(spe, assay.type='logcounts', features=target.features[seq(2)], weight.method='knn', k=50) plot_lisa_feature(spe, lisa.res=lisa.res1, features=target.features[seq(2)], pointsize=2, hlpointsize=2, gap_line_width=.1) ## End(Not run)
## Not run: library(ggplot2) library(SingleCellExperiment) |> suppressPackageStartupMessages() library(SpatialExperiment) |> suppressPackageStartupMessages() library(STexampleData) # create ExperimentHub instance eh <- ExperimentHub() # query STexampleData datasets myfiles <- query(eh, "STexampleData") ah_id <- myfiles$ah_id[myfiles$title == 'Visium_humanDLPFC'] spe <- myfiles[[ah_id]] spe <- spe[, colData(spe)$in_tissue == 1] spe <-scater::logNormCounts(spe) genes <- c('MOBP', 'PCP4', 'SNAP25', 'HBB', 'IGKC', 'NPY') target.features <- rownames(spe)[match(genes, rowData(spe)$gene_name)] library(SVP) lisa.res1 <- runLISA(spe, assay.type='logcounts', features=target.features[seq(2)], weight.method='knn', k=50) plot_lisa_feature(spe, lisa.res=lisa.res1, features=target.features[seq(2)], pointsize=2, hlpointsize=2, gap_line_width=.1) ## End(Not run)
sc_dim
sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... )
sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_dim( object, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, geom = sc_geom_point, ... )
object |
Seurat object or SingleCellExperiment object |
dims |
selected dimensions (must be a two-length vector) that are used in visualization |
reduction |
reduction method, default is NULL and will use the default setting store in the object |
cells |
selected cells to plot (default is all cells) |
slot |
slot to pull expression data from (e.g., 'count' or 'data') |
mapping |
aesthetic mapping, the |
geom |
the function of geometric layer, default is sc_geom_point,
other geometric layer, such as |
... |
additional parameters pass to
|
dimension reduction plot
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_label() f2 <- p2 + sc_dim_geom_label( geom = shadowtext::geom_shadowtext, color='black', bg.color='white' )
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_label() f2 <- p2 + sc_dim_geom_label( geom = shadowtext::geom_shadowtext, color='black', bg.color='white' )
sc_dim_count
sc_dim_count(sc_dim_plot)
sc_dim_count(sc_dim_plot)
sc_dim_plot |
dimension reduction plot of single cell data |
a bar plot to present the cell numbers of different clusters
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p <- sc_dim(sce, reduction = 'UMAP') p1 <- sc_dim_count(p)
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p <- sc_dim(sce, reduction = 'UMAP') p1 <- sc_dim_count(p)
sc_dim_geom_ellipse
sc_dim_geom_ellipse(geom = stat_ellipse, mapping = NULL, level = 0.95, ...)
sc_dim_geom_ellipse(geom = stat_ellipse, mapping = NULL, level = 0.95, ...)
geom |
the layer function, default is |
mapping |
aesthetic mapping |
level |
the level at which to draw an ellipse |
... |
additional parameters pass to the stat_ellipse |
layer of ellipse
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_ellipse()
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_ellipse()
sc_dim_geom_feature
sc_dim_geom_feature( object, features, dims = c(1, 2), ncol = 3, ..., .fun = function(.data) dplyr::filter(.data, .data$value > 0) )
sc_dim_geom_feature( object, features, dims = c(1, 2), ncol = 3, ..., .fun = function(.data) dplyr::filter(.data, .data$value > 0) )
object |
Seurat or SingleCellExperiment object |
features |
selected features (i.e., genes) |
dims |
selected dimensions (must be a two-length vector) that are used in visualization |
ncol |
number of facet columns if 'length(features) > 1' |
... |
additional parameters pass to 'scattermore::geom_scattermore()' |
.fun |
user defined function that will be applied to selected features (default is to filter out genes with no expression values) |
layer of points for selected features
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') set.seed(123) genes <- rownames(sce) |> sample(6) f1 <- p1 + sc_dim_geom_feature( object = sce, features = genes )
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') set.seed(123) genes <- rownames(sce) |> sample(6) f1 <- p1 + sc_dim_geom_feature( object = sce, features = genes )
sc_dim_geom_label
sc_dim_geom_label(geom = ggplot2::geom_text, mapping = NULL, ...)
sc_dim_geom_label(geom = ggplot2::geom_text, mapping = NULL, ...)
geom |
geometric layer (default: geom_text) to display the lables |
mapping |
aesthetic mapping |
... |
additional parameters pass to the geom |
layer of labels
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_label()
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP', mapping = aes(colour = Cell_Cycle)) p2 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_label()
sc_dim_geom_subset
sc_dim_geom_sub(mapping = NULL, subset, .column = "ident", ...)
sc_dim_geom_sub(mapping = NULL, subset, .column = "ident", ...)
mapping |
aesthetic mapping |
subset |
subset of clusters to be displayed |
.column |
which column represents cluster (e.g., 'ident') |
... |
additional parameters pass to sc_geom_point |
plot with a layer of specified clusters
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_sub(subset = c(1, 2), .column = 'label', bg_colour='black')
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_geom_sub(subset = c(1, 2), .column = 'label', bg_colour='black')
sc_dim_sub
sc_dim_sub(subset, .column = "ident")
sc_dim_sub(subset, .column = "ident")
subset |
subset of clusters to be displayed |
.column |
which column represents cluster (e.g., 'ident') |
update plot with only subset displayed
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_sub(subset = c(1, 2), .column = 'label')
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') p1 <- sc_dim(sce, reduction = 'UMAP') f1 <- p1 + sc_dim_sub(subset = c(1, 2), .column = 'label')
sc_dot
sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... ) ## S4 method for signature 'Seurat' sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... ) ## S4 method for signature 'SingleCellExperiment' sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... )
sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... ) ## S4 method for signature 'Seurat' sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... ) ## S4 method for signature 'SingleCellExperiment' sc_dot( object, features, group.by = NULL, split.by = NULL, cols = c("lightgrey", "blue"), col.min = -2.5, col.max = 2.5, dot.min = 0, dot.scale = 6, slot = "data", .fun = NULL, mapping = NULL, scale = TRUE, scale.by = "radius", scale.min = NA, scale.max = NA, cluster.idents = FALSE, ... )
object |
Seurat or SingleCellExperiment object |
features |
selected features |
group.by |
grouping factor |
split.by |
additional split factor |
cols |
colors of the points |
col.min |
minimum scaled averaged expression threshold |
col.max |
maximum scaled averaged expression threshold |
dot.min |
the threshold of percentage of cells for the the smallest dot |
dot.scale |
Scaling factor for size of points |
slot |
slot to pull expression data from (e.g., 'count' or 'data') |
.fun |
user defined function that will be applied to selected features (default is NULL and there is no data operation) |
mapping |
aesthetic mapping |
scale |
whether to scale the expression value (default to TRUE) |
scale.by |
scale the size of the points by |
scale.min |
lower limit of scaling |
scale.max |
upper limit of scaling |
cluster.idents |
Order identities by hierarchical clusters based on average expression and perventage of expression (default is FALSE) |
... |
additional parameters pass to 'ggplot2::geom_point()' |
dot plot to visualize feature expression distribution
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) set.seed(123) genes <- rownames(sce) |> sample(6) sc_dot(sce, genes[1:5], 'Treatment', slot = 'logcounts')
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) set.seed(123) genes <- rownames(sce) |> sample(6) sc_dot(sce, genes[1:5], 'Treatment', slot = 'logcounts')
sc_feature
sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... )
sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_feature( object, features, dims = c(1, 2), reduction = NULL, cells = NULL, slot = "data", mapping = NULL, ncol = 3, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, geom = sc_geom_point, ... )
object |
Seurat object |
features |
selected features (i.e., genes) |
dims |
selected dimensions (must be a two-length vector) that are used in visualization |
reduction |
reduction method, default is NULL and will use the default setting store in the object |
cells |
selected cells to plot (default is all cells) |
slot |
slot to pull expression data from (e.g., 'count' or 'data') |
mapping |
aesthetic mapping |
ncol |
number of facet columns if 'length(features) > 1' |
density |
whether plot the 2D weighted kernel density, default is FALSE. |
grid.n |
number of grid points in the two directions to estimate 2D weighted kernel density, default is 100. |
joint |
whether joint the multiple features with |
joint.fun |
how to joint the multiple features if |
common.legend |
whether to use |
geom |
the function of geometric layer, default is sc_geom_point,
other geometric layer, such as |
... |
additional parameters pass to 'scattermore::geom_scattermore()'
|
dimension reduction plot colored by selected features
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runTSNE(sce, assay.type = 'logcounts') set.seed(123) genes <- rownames(sce) |> sample(6) p1 <- sc_feature(sce, genes[1], slot='logcounts', reduction = 'TSNE') p2 <- sc_feature(sce, genes, slot='logcounts', reduction = 'TSNE') f1 <- sc_dim(sce, slot='logcounts', reduction = 'TSNE') + sc_dim_geom_feature(sce, genes[1], color='black') f2 <- sc_dim(sce, alpha=.3, slot='logcounts', reduction = 'TSNE') + ggnewscale::new_scale_color() + sc_dim_geom_feature(sce, genes, mapping=aes(color=features)) + scale_color_viridis_d() p1 + p2 + f1 + f2 # The features can also be specified the variables from # colData or reducedDims pp <- sc_feature(sce, features = 'sizeFactor', reduction='TSNE', geom=geom_bgpoint) pp
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runTSNE(sce, assay.type = 'logcounts') set.seed(123) genes <- rownames(sce) |> sample(6) p1 <- sc_feature(sce, genes[1], slot='logcounts', reduction = 'TSNE') p2 <- sc_feature(sce, genes, slot='logcounts', reduction = 'TSNE') f1 <- sc_dim(sce, slot='logcounts', reduction = 'TSNE') + sc_dim_geom_feature(sce, genes[1], color='black') f2 <- sc_dim(sce, alpha=.3, slot='logcounts', reduction = 'TSNE') + ggnewscale::new_scale_color() + sc_dim_geom_feature(sce, genes, mapping=aes(color=features)) + scale_color_viridis_d() p1 + p2 + f1 + f2 # The features can also be specified the variables from # colData or reducedDims pp <- sc_feature(sce, features = 'sizeFactor', reduction='TSNE', geom=geom_bgpoint) pp
add the annotation layer for ggsc object
sc_geom_annot( data = NULL, mapping = NULL, pointsize = 2, pixels = c(512, 512), gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1, show.legend = NA, na.rm = FALSE, ... )
sc_geom_annot( data = NULL, mapping = NULL, pointsize = 2, pixels = c(512, 512), gap_colour = "white", gap_alpha = 1, bg_line_width = 0.3, gap_line_width = 0.1, show.legend = NA, na.rm = FALSE, ... )
data |
The data to be displayed in this layer. There are three
options:
If |
mapping |
Set of aesthetic mappings created by |
pointsize |
Radius of rasterized point. Use ‘0’ for single pixels (fastest). |
pixels |
Vector with X and Y resolution of the raster, default |
gap_colour |
colour of gap background between the bottom background
and top point point layer, default is |
gap_alpha |
numeric the transparency of gap background colour, default is 1. |
bg_line_width |
numeric the line width of background point layer,
default is |
gap_line_width |
numeric the line width of gap between the background and
top point point layer, default is |
show.legend |
logical. Should this layer be included in the legends?
|
na.rm |
If |
... |
Other arguments passed on to |
layer object
sc_geom_point
sc_geom_point(mapping = NULL, ...)
sc_geom_point(mapping = NULL, ...)
mapping |
aesthetic mapping |
... |
additional parameters pass to 'scattermore::geom_scattermore()' |
layer of points
sc_dim()
and sc_feature()
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species) ) + sc_geom_point()
library(ggplot2) ggplot(iris, aes(x= Sepal.Length, y = Petal.Width, color=Species) ) + sc_geom_point()
sc_spatial
sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = "data", plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = NULL, remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = "data", plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = NULL, remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = 1, plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = "v", remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... )
sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = "data", plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = NULL, remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... ) ## S4 method for signature 'Seurat' sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = "data", plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = NULL, remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... ) ## S4 method for signature 'SingleCellExperiment' sc_spatial( object, features = NULL, sample.id = NULL, image.id = NULL, slot = 1, plot.pie = FALSE, pie.radius.scale = 0.3, image.plot = TRUE, image.first.operation = "rotate", image.rotate.degree = NULL, image.mirror.axis = "v", remove.point = FALSE, mapping = NULL, ncol = 6, density = FALSE, grid.n = 100, joint = FALSE, joint.fun = prod, common.legend = TRUE, pointsize = 5, geom = sc_geom_point, ... )
object |
Seurat object |
features |
selected features to be visualized |
sample.id |
the index name of sample id, which only work with SingleCellExperiment or SpatialExperiment. |
image.id |
the index name of image id, which only work with SingleCellExperiment or SpatialExperiment. |
slot |
if plotting a feature, which data will be used (e.g., 'data', 'counts'), the assay name if object is SingleCellExperiment or SpatialExperiment. |
plot.pie |
logical whether plot the features with pie, default is |
pie.radius.scale |
numeric scale to the radius of pie only work with |
image.plot |
whether to display the issue image as background. |
image.first.operation |
character which the first operation to image, 'rotate' or 'mirror', default is 'rotate'. |
image.rotate.degree |
integer the degree to ratate image, default is NULL. |
image.mirror.axis |
character the direction to mirror the image, default is 'h'. |
remove.point |
whether to remove the spot points, it is nice if your just view the issue image, default is FALSE. |
mapping |
aesthetic mapping, default is NULL. |
ncol |
integer number of facet columns if 'length(features) > 1', default is 6. |
density |
whether plot the 2D weighted kernel density, default is FALSE. |
grid.n |
number of grid points in the two directions to estimate 2D weighted kernel density, default is 100. |
joint |
whether joint the multiple features with |
joint.fun |
how to joint the multiple features if |
common.legend |
whether to use |
pointsize |
the size of point, default is 5. |
geom |
the layer of point, default is |
... |
additional parameters, see also
|
ggplot object
## Not run: library(STexampleData) # create ExperimentHub instance eh <- ExperimentHub() # query STexampleData datasets myfiles <- query(eh, "STexampleData") ah_id <- myfiles$ah_id[myfiles$title == 'Visium_humanDLPFC'] spe <- myfiles[[ah_id]] spe <- spe[, colData(spe)$in_tissue == 1] set.seed(123) genes <- rownames(spe) |> sample(6) p <- sc_spatial(spe, features = genes, image.rotate.degree = -90, image.mirror.axis = NULL, ncol = 3) # The features also can be specified # the variables from colData or reducedDims. p1 <- sc_spatial(spe, features = 'cell_count', image.rotate.degree = -90, image.mirror.axis = NULL) ## End(Not run)
## Not run: library(STexampleData) # create ExperimentHub instance eh <- ExperimentHub() # query STexampleData datasets myfiles <- query(eh, "STexampleData") ah_id <- myfiles$ah_id[myfiles$title == 'Visium_humanDLPFC'] spe <- myfiles[[ah_id]] spe <- spe[, colData(spe)$in_tissue == 1] set.seed(123) genes <- rownames(spe) |> sample(6) p <- sc_spatial(spe, features = genes, image.rotate.degree = -90, image.mirror.axis = NULL, ncol = 3) # The features also can be specified # the variables from colData or reducedDims. p1 <- sc_spatial(spe, features = 'cell_count', image.rotate.degree = -90, image.mirror.axis = NULL) ## End(Not run)
sc_violin
sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... ) ## S4 method for signature 'Seurat' sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... ) ## S4 method for signature 'SingleCellExperiment' sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... )
sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... ) ## S4 method for signature 'Seurat' sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... ) ## S4 method for signature 'SingleCellExperiment' sc_violin( object, features, cells = NULL, slot = "data", .fun = NULL, mapping = NULL, ncol = 3, geom = geom_violin, ... )
object |
Seurat object |
features |
selected features |
cells |
selected cells to plot (default is all cells) |
slot |
slot to pull expression data from (e.g., 'count' or 'data') |
.fun |
user defined function that will be applied to selected features (default is NULL and there is no data operation) |
mapping |
aesthetic mapping |
ncol |
number of facet columns if 'length(features) > 1' |
geom |
the geom function, default is geom_violin, other option is geom_boxplot |
... |
additional parameters pass to 'ggplot2::geom_geom_violin()' |
violin plot to visualize feature expression distribution
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') set.seed(123) genes <- rownames(sce) |> sample(6) sc_violin(sce, genes[1], slot = 'logcounts') sc_violin(sce, genes[1], slot = 'logcounts', .fun=function(d) dplyr::filter(d, value > 0) ) + ggforce::geom_sina(size=.1) sc_violin(sce, genes, slot = 'logcounts') + theme(axis.text.x = element_text(angle=45, hjust=1))
library(scuttle) library(scater) library(scran) library(ggplot2) sce <- mockSCE() sce <- logNormCounts(sce) clusters <- clusterCells(sce, assay.type = 'logcounts') colLabels(sce) <- clusters sce <- runUMAP(sce, assay.type = 'logcounts') set.seed(123) genes <- rownames(sce) |> sample(6) sc_violin(sce, genes[1], slot = 'logcounts') sc_violin(sce, genes[1], slot = 'logcounts', .fun=function(d) dplyr::filter(d, value > 0) ) + ggforce::geom_sina(size=.1) sc_violin(sce, genes, slot = 'logcounts') + theme(axis.text.x = element_text(angle=45, hjust=1))
Create your own discrete scale
scale_bg_colour_identity( name = waiver(), ..., guide = "none", aesthetics = "bg_colour" ) scale_bg_colour_manual( ..., values, aesthetics = "bg_colour", breaks = waiver(), na.value = "grey50" )
scale_bg_colour_identity( name = waiver(), ..., guide = "none", aesthetics = "bg_colour" ) scale_bg_colour_manual( ..., values, aesthetics = "bg_colour", breaks = waiver(), na.value = "grey50" )
... |
Arguments passed on to
|
name |
The name of the scale. Used as the axis or legend title. If
|
guide |
A function used to create a guide or its name. See
|
aesthetics |
The names of the aesthetics that this scale works with. |
values |
a set of aesthetic values to map data values to. If this
is a named vector, then the values will be matched based on the names.
If unnamed, values will be matched in order (usually alphabetical) with
the limits of the scale. Any data values that don't match will be
given |
breaks |
One of:
|
na.value |
If |
bg_colour scale constructor