plot_hexbin_fc: Plot of fold change of selected gene in single cell data...

View source: R/plot_hexbin_fc.R

plot_hexbin_fcR Documentation

Plot of fold change of selected gene in single cell data using bivariate hexagon cells.

Description

Plot of fold change of selected gene in single cell data using bivariate hexagon cells.

Usage

plot_hexbin_fc(
  sce,
  col,
  mod = "RNA",
  type,
  feature,
  title = NULL,
  xlab = NULL,
  ylab = NULL,
  colors
)

Arguments

sce

A SingleCellExperiment object.

col

A string referring to the name of one column in the meta data of sce by which to compare. Note this factor can only contain two levels.

mod

A string referring to the name of one column in the meta data of sce by which to compare. Note this factor can only contain two levels.

type

A string referring to the name of one column in the meta data of sce by which to compare. Note this factor can only contain two levels.

feature

A string referring to the name of one feature.

title

A string containing the title of the plot.

xlab

A string containing the title of the x axis.

ylab

A string containing the title of the y axis.

colors

A vector of strings specifying which colors to use for plotting the different levels in the selected column of the meta data.

Details

This function plots fold change within each hexagon, which are calculated with make_hexbin. Note that the fold change is only accurate if the condition investigated is within the same individual. For conditions across different individuals different methods that account for individual-specific effects are required.

Value

A ggplot2{ggplot} object.

Examples

# For SingleCellExperiment
library(TENxPBMCData)
library(scater)
tenx_pbmc3k <- TENxPBMCData(dataset = "pbmc3k")
rm_ind <- calculateAverage(tenx_pbmc3k) < 0.1
tenx_pbmc3k <- tenx_pbmc3k[!rm_ind, ]
colData(tenx_pbmc3k) <- cbind(colData(tenx_pbmc3k), perCellQCMetrics(tenx_pbmc3k))
tenx_pbmc3k <- logNormCounts(tenx_pbmc3k)
tenx_pbmc3k <- runPCA(tenx_pbmc3k)
tenx_pbmc3k <- make_hexbin(tenx_pbmc3k, 20, dimension_reduction = "PCA")
tenx_pbmc3k$random <- factor(sample(1:2, ncol(tenx_pbmc3k), replace = TRUE))
plot_hexbin_fc(tenx_pbmc3k, col = "random", feature = "ENSG00000187608", type = "counts")

SaskiaFreytag/schex documentation built on Feb. 4, 2024, 7:49 p.m.