R/runFindMarker.R
getFindMarkerTopTable.RdFetch the table of top markers that pass the filtering
getFindMarkerTopTable(
inSCE,
log2fcThreshold = 0,
fdrThreshold = 0.05,
minClustExprPerc = 0.5,
maxCtrlExprPerc = 0.5,
minMeanExpr = 0,
topN = 1
)
findMarkerTopTable(
inSCE,
log2fcThreshold = 1,
fdrThreshold = 0.05,
minClustExprPerc = 0.7,
maxCtrlExprPerc = 0.4,
minMeanExpr = 1,
topN = 10
)SingleCellExperiment inherited object.
Only use DEGs with the absolute values of log2FC
larger than this value. Default 1
Only use DEGs with FDR value smaller than this value.
Default 0.05
A numeric scalar. The minimum cutoff of the
percentage of cells in the cluster of interests that expressed the marker
gene. Default 0.7.
A numeric scalar. The maximum cutoff of the
percentage of cells out of the cluster (control group) that expressed the
marker gene. Default 0.4.
A numeric scalar. The minimum cutoff of the mean
expression value of the marker in the cluster of interests. Default 1.
An integer. Only to fetch this number of top markers for each
cluster in maximum, in terms of log2FC value. Use NULL to cancel the
top N subscription. Default 10.
An organized data.frame object, with the top marker gene
information.
Users have to run runFindMarker prior to using this
function to extract a top marker table.
data("mouseBrainSubsetSCE", package = "singleCellTK")
mouseBrainSubsetSCE <- runFindMarker(mouseBrainSubsetSCE,
useAssay = "logcounts",
cluster = "level1class")
#> Wed Apr 23 11:37:30 2025 ... Identifying markers for cluster 'microglia', using DE method 'wilcox'
#> Wed Apr 23 11:37:31 2025 ... Identifying markers for cluster 'oligodendrocytes', using DE method 'wilcox'
#> Wed Apr 23 11:37:31 2025 ... Organizing findMarker result
getFindMarkerTopTable(mouseBrainSubsetSCE)
#> Gene Log2_FC Pvalue FDR level1class clusterExprPerc
#> 1228 Apoe 4.439115 3.642815e-06 0.0001796628 microglia 1
#> 1002 Mal 5.361403 7.495067e-06 0.0002501837 oligodendrocytes 1
#> ControlExprPerc clusterAveExpr
#> 1228 0.3333333 5.076598
#> 1002 0.4000000 6.370450