View source: R/medianByClusterMarker.R
medianByClusterMarker | R Documentation |
Calculate median value of each marker in each cluster.
medianByClusterMarker(
SE,
assay = 1,
marker_in_column = TRUE,
column_cluster = "cluster_id",
use_marker = NULL
)
SE |
A |
assay |
A numeric index or assay name indicating with assay of |
marker_in_column |
A logical scalar, indicating whether markers (genes,
features) are in the columns of |
column_cluster |
The name of the column of |
use_marker |
A logical or numeric vector such that
|
A SummarizedExperiment
object
containing the median value of each marker in each cluster.
Ruizhu Huang, Charlotte Soneson
suppressPackageStartupMessages({
library(SummarizedExperiment)
})
## Simulate data with 100 cells and 10 markers (5 type, 5 state markers)
set.seed(1)
count <- matrix(rpois(n = 1000, lambda = 10), nrow = 100)
colnames(count) <- paste0("mk", 1:10)
rowD <- data.frame("cluster" = sample(seq_len(6), 100, replace = TRUE))
colD <- data.frame(type_marker = rep(c(FALSE, TRUE), each = 5))
## SE with markers in columns
d_se <- SummarizedExperiment(assays = list(counts = count),
rowData = rowD,
colData = colD)
medianByClusterMarker(SE = d_se, marker_in_column = TRUE,
column_cluster = "cluster",
use_marker = colData(d_se)$type_marker)
## SE with markers in rows
d_se <- SummarizedExperiment(assays = list(counts = t(count)),
rowData = colD,
colData = rowD)
medianByClusterMarker(SE = d_se, marker_in_column = FALSE,
column_cluster = "cluster",
use_marker = rowData(d_se)$type_marker)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.