| topClusters {diffcyt} | R Documentation |
Show results for top (most highly significant) clusters or cluster-marker combinations
topClusters(res, order = TRUE, all = FALSE, top_n = 20)
res |
Output object containing results from one of the differential testing
functions, in |
order |
Whether to order results by adjusted p-value. Default = TRUE. |
all |
Whether to display all clusters or cluster-marker combinations (instead of
top |
top_n |
Number of clusters or cluster-marker combinations to display (if |
Summary function to display results for top (most highly significant) detected clusters or cluster-marker combinations.
The differential testing functions return results in the form of p-values and adjusted
p-values for each cluster (DA tests) or cluster-marker combination (DS tests), which
can be used to rank the clusters or cluster-marker combinations by their evidence for
differential abundance or differential states. The p-values and adjusted p-values are
stored in the rowData of the output SummarizedExperiment object
generated by the testing functions.
This summary function displays the rowData from the SummarizedExperiment
output object as a DataFrame, with clusters or cluster-marker
combinations ordered by the adjusted p-values.
Returns a DataFrame of results for the top_n clusters or
cluster-marker combinations, ordered by adjusted p-values.
# 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 = 2)
# 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, plot = FALSE, verbose = FALSE)
# Display results for top DA clusters
topClusters(out_DA$res)
# Display results for top DS cluster-marker combinations
topClusters(out_DS$res)