With an input SingleCellExperiment object and specifying the
clustering labels, this function iteratively call the differential expression
analysis on each cluster against all the others. runFindMarker
will be deprecated in the future.
runFindMarker(
inSCE,
useAssay = "logcounts",
useReducedDim = NULL,
method = c("wilcox", "MAST", "DESeq2", "Limma", "ANOVA"),
cluster = "cluster",
covariates = NULL,
log2fcThreshold = NULL,
fdrThreshold = 0.05,
minClustExprPerc = NULL,
maxCtrlExprPerc = NULL,
minMeanExpr = NULL,
detectThresh = 0
)
findMarkerDiffExp(
inSCE,
useAssay = "logcounts",
useReducedDim = NULL,
method = c("wilcox", "MAST", "DESeq2", "Limma", "ANOVA"),
cluster = "cluster",
covariates = NULL,
log2fcThreshold = NULL,
fdrThreshold = 0.05,
minClustExprPerc = NULL,
maxCtrlExprPerc = NULL,
minMeanExpr = NULL,
detectThresh = 0
)
SingleCellExperiment inherited object.
character. A string specifying which assay to use for the
MAST calculations. Default "logcounts"
.
character. A string specifying which reducedDim to use
for MAST calculations. Set useAssay
to NULL
when using.
Required.
A single character for specific differential expression
analysis method. Choose from 'wilcox'
, 'MAST'
, 'DESeq2'
,
'Limma'
, and 'ANOVA'
. Default "wilcox"
.
One single character to specify a column in
colData(inSCE)
for the clustering label. Alternatively, a vector or
a factor is also acceptable. Default "cluster"
.
A character vector of additional covariates to use when
building the model. All covariates must exist in
names(colData(inSCE))
. Not applicable when method
is
"MAST"
method. Default NULL
.
Only out put DEGs with the absolute values of log2FC
larger than this value. Default NULL
Only out put DEGs with FDR value smaller than this
value. Default NULL
A numeric scalar. The minimum cutoff of the
percentage of cells in the cluster of interests that expressed the marker
gene. From 0 to 1. Default NULL
.
A numeric scalar. The maximum cutoff of the
percentage of cells out of the cluster (control group) that expressed the
marker gene. From 0 to 1. Default NULL
.
A numeric scalar. The minimum cutoff of the mean
expression value of the marker in the cluster of interests. Default
NULL
.
A numeric scalar, above which a matrix value will be
treated as expressed when calculating cluster/control expression percentage.
Default 0
.
The input SingleCellExperiment object with
metadata(inSCE)$findMarker
updated with a data.table of the up-
regulated DEGs for each cluster.
The returned marker table, in the metadata
slot, consists of 8
columns: "Gene"
, "Log2_FC"
, "Pvalue"
, "FDR"
,
cluster
, "clusterExprPerc"
, "ControlExprPerc"
and
"clusterAveExpr"
.
"clusterExprPerc"
is the fraction of cells,
that has marker value (e.g. gene expression counts) larger than
detectThresh
, in the cell population of the cluster. As for each
cluster, we set all cells out of this cluster as control. Similarly,
"ControlExprPerc"
is the fraction of cells with marker value larger
than detectThresh
in the control cell group.
data("mouseBrainSubsetSCE", package = "singleCellTK")
mouseBrainSubsetSCE <- runFindMarker(mouseBrainSubsetSCE,
useAssay = "logcounts",
cluster = "level1class")
#> Sat Mar 18 10:31:00 2023 ... Identifying markers for cluster 'microglia', using DE method 'wilcox'
#> Sat Mar 18 10:31:01 2023 ... Identifying markers for cluster 'oligodendrocytes', using DE method 'wilcox'
#> Sat Mar 18 10:31:02 2023 ... Organizing findMarker result