diff --git a/DESCRIPTION b/DESCRIPTION index e38608d..9332d76 100644 --- a/DESCRIPTION +++ b/DESCRIPTION @@ -1,7 +1,7 @@ Package: SingleR Title: Reference-Based Single-Cell RNA-Seq Annotation -Version: 2.7.8 -Date: 2024-10-26 +Version: 2.7.9 +Date: 2024-10-28 Authors@R: c(person("Dvir", "Aran", email="dvir.aran@ucsf.edu", role=c("aut", "cph")), person("Aaron", "Lun", email="infinite.monkeys.with.keyboards@gmail.com", role=c("ctb", "cre")), person("Daniel", "Bunis", role="ctb"), diff --git a/R/plotMarkerHeatmap.R b/R/plotMarkerHeatmap.R index 4a37c90..a9a6eaa 100644 --- a/R/plotMarkerHeatmap.R +++ b/R/plotMarkerHeatmap.R @@ -5,7 +5,7 @@ #' @param results A \linkS4class{DataFrame} containing the output from \code{\link{SingleR}}, #' \code{\link{classifySingleR}}, \code{\link{combineCommonResults}}, or \code{\link{combineRecomputedResults}}. #' @param test A numeric matrix of log-normalized expression values where rows are genes and columns are cells. -#' Each row should be named with the gene name. +#' Each row should be named with the same gene name that was used to compute \code{results}. #' #' Alternatively, a \linkS4class{SummarizedExperiment} object containing such a matrix. #' @param label String specifying the label of interest. @@ -14,6 +14,9 @@ #' @param assay.type Integer scalar or string specifying the matrix of expression values to use if \code{test} is a \linkS4class{SummarizedExperiment}. #' @param use.pruned Logical scalar indicating whether the pruned labels should be used instead. #' @param order.by String specifying the column of the output of \code{\link[scran]{scoreMarkers}} with which to sort for interesting markers. +#' @param display.row.names Character vector of length equal to the number of rows of \code{test}, +#' containing the names of the features to show on the heatmap (e.g., to replace IDs with symbols). +#' If \code{NULL}, the existing row names of \code{test} are used. #' @param top Integer scalar indicating the top most interesting markers to show in the heatmap. #' @param BPPARAM A \linkS4class{BiocParallelParam} object specifying the parallelization scheme to use for marker scoring. #' @@ -48,6 +51,7 @@ plotMarkerHeatmap <- function( label, other.labels=NULL, assay.type="logcounts", + display.row.names=NULL, use.pruned=FALSE, order.by="rank.logFC.cohen", top=20, @@ -92,13 +96,19 @@ plotMarkerHeatmap <- function( to.show <- names(abundance)[order(abundance, decreasing=TRUE)] } - to.show <- head(to.show, top) - limits <- range(test, na.rm=TRUE) + to.show <- match(head(to.show, top), rownames(test)) colnames(test) <- seq_len(ncol(test)) col.order <- order(predictions) + test <- test[to.show,col.order,drop=FALSE] + predictions <- predictions[col.order] + + if (!is.null(display.row.names)) { + rownames(test) <- display.row.names[rkeep][to.show] + } + limits <- range(test, na.rm=TRUE) pheatmap::pheatmap( - test[to.show,col.order,drop=FALSE], + test, breaks=seq(limits[1], limits[2], length.out=26), color=viridis::viridis(25), annotation_col=data.frame(labels=predictions, row.names=colnames(test)), diff --git a/inst/NEWS.Rd b/inst/NEWS.Rd index a99e3b9..4c7c6bc 100644 --- a/inst/NEWS.Rd +++ b/inst/NEWS.Rd @@ -21,6 +21,8 @@ Rather, filtering should be done prior to \code{trainSingleR()}, as is done in t \item \code{combineRecomputedResults()} now supports fine-tuning to resolve closely-related labels from different references. This is similar to the fine-tuning in \code{classifySingleR()} where the feature space is iterately redefined as the union of markers of labels with near-highest scores. + +\item Added the \code{plotMarkerHeatmap()} function to plot a diagnostic heatmap of the most interesting markers for each label. }} \section{Version 2.0.0}{\itemize{ diff --git a/man/plotMarkerHeatmap.Rd b/man/plotMarkerHeatmap.Rd index d020a7a..4603b22 100644 --- a/man/plotMarkerHeatmap.Rd +++ b/man/plotMarkerHeatmap.Rd @@ -10,6 +10,7 @@ plotMarkerHeatmap( label, other.labels = NULL, assay.type = "logcounts", + display.row.names = NULL, use.pruned = FALSE, order.by = "rank.logFC.cohen", top = 20, @@ -21,7 +22,7 @@ plotMarkerHeatmap( \code{\link{classifySingleR}}, \code{\link{combineCommonResults}}, or \code{\link{combineRecomputedResults}}.} \item{test}{A numeric matrix of log-normalized expression values where rows are genes and columns are cells. -Each row should be named with the gene name. +Each row should be named with the same gene name that was used to compute \code{results}. Alternatively, a \linkS4class{SummarizedExperiment} object containing such a matrix.} @@ -32,6 +33,10 @@ Defaults to all available labels.} \item{assay.type}{Integer scalar or string specifying the matrix of expression values to use if \code{test} is a \linkS4class{SummarizedExperiment}.} +\item{display.row.names}{Character vector of length equal to the number of rows of \code{test}, +containing the names of the features to show on the heatmap (e.g., to replace IDs with symbols). +If \code{NULL}, the existing row names of \code{test} are used.} + \item{use.pruned}{Logical scalar indicating whether the pruned labels should be used instead.} \item{order.by}{String specifying the column of the output of \code{\link[scran]{scoreMarkers}} with which to sort for interesting markers.} diff --git a/tests/testthat/test-plotMarkerHeatmap.R b/tests/testthat/test-plotMarkerHeatmap.R index 64407a9..1e559ab 100644 --- a/tests/testthat/test-plotMarkerHeatmap.R +++ b/tests/testthat/test-plotMarkerHeatmap.R @@ -18,6 +18,11 @@ test_that("marker heatmaps work with a subset of other labels", { expect_s3_class(plotMarkerHeatmap(pred, test, pred$labels[1], other.labels=unique(pred$labels)[1:2]), "pheatmap") }) +test_that("marker heatmaps work with other names", { + alt.names <- paste0("X_", rownames(test)) + expect_s3_class(plotMarkerHeatmap(pred, test, pred$labels[1], display.row.names=alt.names), "pheatmap") +}) + test_that("marker heatmap falls back to average abundances", { lab <- pred$labels[1] keep <- pred$labels == lab