Plot Heatmap of Gene Expression or Factor Loading
Usage
plotGeneHeatmap(
object,
features,
cellIdx = NULL,
slot = c("normData", "rawData", "scaleData", "scaleUnsharedData"),
useCellMeta = NULL,
cellAnnotation = NULL,
featureAnnotation = NULL,
cellSplitBy = NULL,
featureSplitBy = NULL,
viridisOption = "C",
...
)
plotFactorHeatmap(
object,
factors = NULL,
cellIdx = NULL,
slot = c("H.norm", "H"),
useCellMeta = NULL,
cellAnnotation = NULL,
factorAnnotation = NULL,
cellSplitBy = NULL,
factorSplitBy = NULL,
trim = c(0, 0.03),
viridisOption = "D",
...
)
Arguments
- object
A liger object, with data to be plot available.
- features, factors
Character vector of genes of interests or numeric index of factor to be involved.
features
is required, whilefactors
is by default all the factors (reads object recorded k value inuns
slot).- cellIdx
Valid index to subscribe cells to be included. See
subsetLiger
. DefaultNULL
use all cells.- slot
Use the chosen matrix for heatmap. For
plotGeneHeatmap
, default"normData"
, alternatively"rawData"
,"scaleData"
or"scaleUnsharedData"
. ForplotFactorHeatmap
, default"H.norm"
, alternatively"H"
.- useCellMeta
Character vector of available variable names in
cellMeta
, variables will be added as annotation to the heatmap. DefaultNULL
.- cellAnnotation
data.frame object for using external annotation, with each column a variable and each row is a cell. Row names of this data.frame will be used for matching cells involved in heatmap. For cells not found in this data.frame,
NA
s will be added with warning. DefaultNULL
.- featureAnnotation, factorAnnotation
Similar as
cellAnnotation
, while each row would be a gene or factor, respectively. DefaultNULL
.- cellSplitBy
Character vector of variable names available in annotation given by
useCellMeta
andcellAnnotation
. This slices the heatmap by specified variables. DefaultNULL
.- featureSplitBy, factorSplitBy
Similar as
cellSplitBy
. DefaultNULL
- viridisOption
See
option
argument ofviridis
. Default"C"
(plasma) forplotGeneHeatmap
and"D"
(viridis) forplotFactorHeatmap
.- ...
Additional arguments passed to general function
.plotHeatmap
andHeatmap
.- trim
Numeric vector of two numbers. Higher value limits the maximum value and lower value limits the minimum value. Default
c(0, 0.03)
.
Value
HeatmapList-class
object
Examples
# \donttest{
plotGeneHeatmap(pbmcPlot, varFeatures(pbmcPlot))
#> ℹ Subsetting dataset: "ctrl"
#> ℹ Subsetting dataset: "stim"
#> ✔ Subsetting dataset: "stim" ... done
#>
#> ℹ Subsetting dataset: "ctrl"
#> ✔ Subsetting dataset: "ctrl" ... done
#>
plotGeneHeatmap(pbmcPlot, varFeatures(pbmcPlot),
useCellMeta = c("leiden_cluster", "dataset"),
cellSplitBy = "leiden_cluster")
#> ℹ Subsetting dataset: "ctrl"
#> ℹ Subsetting dataset: "stim"
#> ✔ Subsetting dataset: "stim" ... done
#>
#> ℹ Subsetting dataset: "ctrl"
#> ✔ Subsetting dataset: "ctrl" ... done
#>
plotFactorHeatmap(pbmcPlot)
plotFactorHeatmap(pbmcPlot, cellIdx = pbmcPlot$leiden_cluster %in% 1:3,
useCellMeta = c("leiden_cluster", "dataset"),
cellSplitBy = "leiden_cluster")
# }