plotHeatmap: Plot heatmap

View source: R/plotHeatmap.R

plotHeatmapR Documentation

Plot heatmap

Description

Plot heatmap showing top clusters or cluster-marker combinations

Usage

plotHeatmap(
  out = NULL,
  analysis_type = c("DA", "DS"),
  top_n = 20,
  threshold = 0.1,
  res = NULL,
  d_se = NULL,
  d_counts = NULL,
  d_medians = NULL,
  d_medians_by_cluster_marker = NULL,
  sample_order = NULL
)

Arguments

out

Output object from diffcyt wrapper function, containing results object res and data objects d_se, d_counts, d_medians, and d_medians_by_cluster_marker. Alternatively, the results and data objects can be provided individually.

analysis_type

Whether to plot heatmap for differential abundance (DA) or differential state (DS) test results.

top_n

Number of top clusters (DA tests) or cluster-marker combinations (DS tests) to display. Default = 20.

threshold

Threshold for significant adjusted p-values. Default = 0.1.

res

Object containing differential test results. Alternatively, the combined output object from the wrapper function diffcyt can be provided.

d_se

Data object. Alternatively, the combined output object from the wrapper function diffcyt can be provided.

d_counts

Data object. Alternatively, the combined output object from the wrapper function diffcyt can be provided.

d_medians

Data object. (Required for DS tests only.) Alternatively, the combined output object from the wrapper function diffcyt can be provided.

d_medians_by_cluster_marker

Data object. Alternatively, the combined output object from the wrapper function diffcyt can be provided.

sample_order

(Optional) Custom ordering for samples (columns) in right-hand panel of heatmap. (This is useful when the default ordering does not group samples by condition; e.g. samples are ordered alphabetically by sample IDs instead.)

Details

Display heatmap to visualize results for the top (most highly significant) detected clusters or cluster-marker combinations.

For DA tests, the heatmap consists of the following panels:

  • median (arcsinh-transformed) expression (across all samples) for 'cell type' markers

  • cluster abundances by sample

  • row annotation indicating significant detected clusters

For DS tests, the heatmap consists of:

  • median (arcsinh-transformed) expression (across all samples) for 'cell type' markers

  • median (arcsinh-transformed) expression (across all samples) for 'cell state' markers

  • median (arcsinh-transformed) expression (by sample) for 'cell state' markers for the top cluster-marker combinations

  • row annotation indicating significant detected cluster-marker combinations

Heatmaps are generated using the ComplexHeatmap package (Gu et al., 2016), and color scales are generated using the circlize package (Gu et al., 2014). Both packages are available from Bioconductor.

Value

Displays a heatmap.

Examples

# For a complete workflow example demonstrating each step in the 'diffcyt' pipeline, 
# see the package vignette.

# Function to create random data (one sample)
d_random <- function(n = 20000, mean = 0, sd = 1, ncol = 20, cofactor = 5) {
  d <- sinh(matrix(rnorm(n, mean, sd), ncol = ncol)) * cofactor
  colnames(d) <- paste0("marker", sprintf("%02d", 1:ncol))
  d
}

# Create random data (without differential signal)
set.seed(123)
d_input <- list(
  sample1 = d_random(), 
  sample2 = d_random(), 
  sample3 = d_random(), 
  sample4 = d_random()
)

# Add differential abundance (DA) signal
ix_DA <- 801:900
ix_cols_type <- 1:10
d_input[[3]][ix_DA, ix_cols_type] <- d_random(n = 1000, mean = 2, ncol = 10)
d_input[[4]][ix_DA, ix_cols_type] <- d_random(n = 1000, mean = 2, ncol = 10)

# Add differential states (DS) signal
ix_DS <- 901:1000
ix_cols_DS <- 19:20
d_input[[1]][ix_DS, ix_cols_type] <- d_random(n = 1000, mean = 3, ncol = 10)
d_input[[2]][ix_DS, ix_cols_type] <- d_random(n = 1000, mean = 3, ncol = 10)
d_input[[3]][ix_DS, c(ix_cols_type, ix_cols_DS)] <- d_random(n = 1200, mean = 3, ncol = 12)
d_input[[4]][ix_DS, c(ix_cols_type, ix_cols_DS)] <- d_random(n = 1200, mean = 3, ncol = 12)

experiment_info <- data.frame(
  sample_id = factor(paste0("sample", 1:4)), 
  group_id = factor(c("group1", "group1", "group2", "group2")), 
  stringsAsFactors = FALSE
)

marker_info <- data.frame(
  channel_name = paste0("channel", sprintf("%03d", 1:20)), 
  marker_name = paste0("marker", sprintf("%02d", 1:20)), 
  marker_class = factor(c(rep("type", 10), rep("state", 10)), 
                        levels = c("type", "state", "none")), 
  stringsAsFactors = FALSE
)

# Create design matrix
design <- createDesignMatrix(experiment_info, cols_design = "group_id")

# Create contrast matrix
contrast <- createContrast(c(0, 1))

# Test for differential abundance (DA) of clusters (using default method 'diffcyt-DA-edgeR')
out_DA <- diffcyt(d_input, experiment_info, marker_info, 
                  design = design, contrast = contrast, 
                  analysis_type = "DA", method_DA = "diffcyt-DA-edgeR", 
                  seed_clustering = 123, verbose = FALSE)

# Test for differential states (DS) within clusters (using default method 'diffcyt-DS-limma')
out_DS <- diffcyt(d_input, experiment_info, marker_info, 
                  design = design, contrast = contrast, 
                  analysis_type = "DS", method_DS = "diffcyt-DS-limma", 
                  seed_clustering = 123, verbose = FALSE)

# Display results for top DA clusters
topTable(out_DA, format_vals = TRUE)

# Display results for top DS cluster-marker combinations
topTable(out_DS, format_vals = TRUE)

# Plot heatmap for DA tests
plotHeatmap(out_DA, analysis_type = "DA")

# Plot heatmap for DS tests
plotHeatmap(out_DS, analysis_type = "DS")


lmweber/diffcyt documentation built on March 19, 2024, 5:24 a.m.