Plot highly variable genes
plotTopHVG(
inSCE,
method = c("vst", "mean.var.plot", "dispersion", "modelGeneVar"),
hvgNumber = NULL,
useFeatureSubset = NULL,
labelsCount = 20,
featureDisplay = metadata(inSCE)$featureDisplay,
labelSize = 2,
dotSize = 2,
textSize = 12
)
Input SingleCellExperiment
object containing the
computations.
Select either "vst"
, "mean.var.plot"
,
"dispersion"
or "modelGeneVar"
.
Specify the number of top genes to highlight in red. Default
NULL
. See details.
A character string for the rowData
variable
name to store a logical index of selected features. Default NULL
. See
details.
Specify the number of data points/genes to label. Should
be less than hvgNumber
. Default 20
. See details.
A character string for the rowData
variable name
to indicate what type of feature ID should be displayed. If set by
setSCTKDisplayRow
, will by default use it. If NULL
, will
use rownames(inSCE)
.
Numeric, size of the text label on top HVGs. Default
2
.
Numeric, size of the dots of the features. Default 2
.
Numeric, size of the text of axis title, axis label, etc.
Default 12
.
ggplot of HVG metrics and top HVG labels
When hvgNumber = NULL
and useFeature = NULL
, only plot
the mean VS variance/dispersion scatter plot. When only hvgNumber
set,
label the top hvgNumber
HVGs ranked by the metrics calculated by
method
. When useFeatureSubset
set, label the features in
the subset on the scatter plot created with method
and ignore
hvgNumber
.
data("mouseBrainSubsetSCE", package = "singleCellTK")
mouseBrainSubsetSCE <- runModelGeneVar(mouseBrainSubsetSCE)
#> Warning: collapsing to unique 'x' values
plotTopHVG(mouseBrainSubsetSCE, method = "modelGeneVar")