Skip to content

Commit de25043

Browse files
committed
Pass along additional arguments to scoreMarkers in plotMarkerHeatmap.
1 parent a5acd72 commit de25043

3 files changed

Lines changed: 33 additions & 25 deletions

File tree

R/plotMarkerHeatmap.R

Lines changed: 26 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,7 @@
1414
#' @param use.pruned Logical scalar indicating whether the pruned labels should be used instead.
1515
#' @param order.by.effect String specifying the effect size from \code{\link[scrapper]{scoreMarkers}} with which to sort for interesting markers.
1616
#' @param order.by.summary String specifying the summary statistic from \code{\link[scrapper]{scoreMarkers}} with which to sort for interesting markers.
17+
#' @param score.args Named list of additional arguments to pass to \code{\link[scrapper]{scoreMarkers}}, e.g., \code{block}, \code{threshold}.
1718
#' @param display.row.names Character vector of length equal to the number of rows of \code{test},
1819
#' containing the names of the features to show on the heatmap (e.g., to replace IDs with symbols).
1920
#' If \code{NULL}, the existing row names of \code{test} are used.
@@ -77,15 +78,16 @@ plotMarkerHeatmap <- function(
7778
results,
7879
test,
7980
label,
80-
other.labels=NULL,
81-
assay.type="logcounts",
82-
display.row.names=NULL,
83-
use.pruned=FALSE,
84-
order.by.effect="cohens.d",
85-
order.by.summary="min.rank",
86-
average=FALSE,
87-
center=FALSE,
88-
top=20,
81+
other.labels = NULL,
82+
assay.type = "logcounts",
83+
display.row.names = NULL,
84+
use.pruned = FALSE,
85+
order.by.effect = "cohens.d",
86+
order.by.summary = "min.rank",
87+
score.args = list(),
88+
average = FALSE,
89+
center = FALSE,
90+
top = 20,
8991
num.threads = 1,
9092
BPPARAM = NULL,
9193
...
@@ -101,6 +103,7 @@ plotMarkerHeatmap <- function(
101103
use.pruned=use.pruned,
102104
order.by.effect=order.by.effect,
103105
order.by.summary=order.by.summary,
106+
score.args=score.args,
104107
num.threads=num.threads
105108
)
106109

@@ -149,13 +152,14 @@ configureMarkerHeatmap <- function(
149152
results,
150153
test,
151154
label,
152-
other.labels=NULL,
153-
assay.type="logcounts",
154-
use.pruned=FALSE,
155-
order.by.effect="cohens.d",
156-
order.by.summary="min.rank",
157-
num.threads=1)
158-
{
155+
other.labels = NULL,
156+
assay.type = "logcounts",
157+
use.pruned = FALSE,
158+
order.by.effect = "cohens.d",
159+
order.by.summary = "min.rank",
160+
score.args = list(),
161+
num.threads=1
162+
) {
159163
test <- .to_clean_matrix(test, assay.type, check.missing=FALSE, num.threads=num.threads)
160164
all.markers <- metadata(results)$de.genes[[label]]
161165

@@ -186,16 +190,13 @@ configureMarkerHeatmap <- function(
186190
# Prioritize the markers with interesting variation in the test data for
187191
# visualization. If we only have one label, we use the most abundant markers.
188192
if (length(unique(predictions)) > 1L) {
189-
interesting <- scrapper::scoreMarkers(
190-
test,
191-
predictions,
192-
num.threads=num.threads,
193-
compute.auc=(order.by.effect=="auc"),
194-
compute.cohens.d=(order.by.effect=="cohens.d"),
195-
compute.delta.mean=(order.by.effect=="delta.mean"),
196-
compute.delta.detected=(order.by.effect=="delta.detected")
197-
)
193+
score.args$num.threads <- num.threads
194+
score.args$compute.auc <- (order.by.effect == "auc")
195+
score.args$compute.cohens.d <- (order.by.effect == "cohens.d")
196+
score.args$compute.delta.mea <- (order.by.effect == "delta.mean")
197+
score.args$compute.delta.detected <- (order.by.effect == "delta.detected")
198198

199+
interesting <- do.call(scrapper::scoreMarkers, c(list(test, predictions), score.args))
199200
stats <- interesting[[order.by.effect]][[label]][[order.by.summary]]
200201
decreasing <- (order.by.summary!="min.rank")
201202

inst/NEWS.Rd

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -11,7 +11,10 @@ Similarly, suggest the use of \code{aggr.ref=TRUE} for large single-cell referen
1111
\item Updated to the latest version of \pkg{singlepp}, which provides direct support for sparse test/reference matrices.
1212
This also eliminates the need for \pkg{BiocNeighbors} so all \code{BNPARAM=} arguments are now soft-deprecated.
1313

14+
\item Added \code{score.args=} to \code{configureMarkerHeatmap()} to pass along to \code{scoreMarkers()}, mainly for \code{block=} and \code{threshold=} options.
15+
1416
\item Added \code{center=} and \code{average=} options to \code{plotMarkerHeatmap()}, for centering of genes and averaging of labels, respectively.
17+
Also pass along \code{score.args=}.
1518

1619
\item Support \code{block=} in \code{getClassicMarkers()} to define a blocking factor within a single reference matrix.
1720
This serves as an alternative to supplying multiple reference matrices.

man/plotMarkerHeatmap.Rd

Lines changed: 4 additions & 0 deletions
Some generated files are not rendered by default. Learn more about customizing how changed files appear on GitHub.

0 commit comments

Comments
 (0)