diff --git a/DESCRIPTION b/DESCRIPTION index 989c525..13f2934 100644 --- a/DESCRIPTION +++ b/DESCRIPTION @@ -1,6 +1,7 @@ Package: velociraptor Title: Toolkit for Single-Cell Velocity -Version: 0.99.8 +Version: 0.99.9 +Date: 2020-10-10 Authors@R: c(person("Kevin", "Rue-Albrecht", role = c("aut", "cre"), email = "kevinrue67@gmail.com", comment = c(ORCID = "0000-0003-3899-3872")), person("Aaron", "Lun", role="aut", email="infinite.monkeys.with.keyboards@gmail.com", comment = c(ORCID = '0000-0002-3564-4813')), person("Charlotte", "Soneson", role="aut", email="charlottesoneson@gmail.com", comment = c(ORCID = '0000-0003-3833-2169'))) diff --git a/NAMESPACE b/NAMESPACE index ffdaa2a..d6bad73 100644 --- a/NAMESPACE +++ b/NAMESPACE @@ -1,8 +1,11 @@ # Generated by roxygen2: do not edit by hand +export(.grid_vectors) export(embedVelocity) export(gridVectors) export(scvelo) +exportMethods(embedVelocity) +exportMethods(gridVectors) exportMethods(scvelo) import(SummarizedExperiment) import(basilisk) diff --git a/NEWS.md b/NEWS.md index c211082..d2488ac 100644 --- a/NEWS.md +++ b/NEWS.md @@ -1,3 +1,7 @@ +# velociraptor 0.99.9 + +* Converted various functions to S4 generics for easier use with `SingleCellExperiment` objects. + # velociraptor 0.99.8 * Trigger new build to repeat `ExperimentHub` download. diff --git a/R/embedVelocity.R b/R/embedVelocity.R index c8dec36..74f0cbf 100644 --- a/R/embedVelocity.R +++ b/R/embedVelocity.R @@ -3,9 +3,15 @@ #' Project the velocity vector for each cell onto an existing low-dimensional embedding. #' #' @param x A numeric matrix of low-dimensional coordinates, e.g., after t-SNE. -#' @param v A \linkS4class{SingleCellExperiment} containing the output of the velocity calculations, +#' Alternatively, a \linkS4class{SingleCellExperiment} containing such coordinates in its \code{\link{reducedDims}}. +#' @param vobj A \linkS4class{SingleCellExperiment} containing the output of the velocity calculations, #' typically after running \code{\link{scvelo}}. -#' @param ... Further arguments to pass to the \code{velocity_embedding} Python function from \pkg{scVelo}. +#' @param ... For the generic, further arguments to pass to specific methods. +#' +#' For the ANY method, further arguments to pass to the \code{velocity_embedding} Python function from \pkg{scVelo}. +#' +#' For the SingleCellExperiment method, further arguments to pass to the ANY method. +#' @param use.dimred String or integer scalar specifying the reduced dimensions to retrieve from \code{x}. #' #' @details #' This is a simple wrapper around the \code{scvelo.tools.velocity_embedding} function. @@ -28,21 +34,38 @@ #' projected <- embedVelocity(tsne.results, out) #' #' @export +#' @name embedVelocity +NULL + #' @importFrom SingleCellExperiment reducedDim<- -embedVelocity <- function(x, v, ...) { - reducedDim(v, "X_target") <- as.matrix(x) - basiliskRun(env=velo.env, fun=.run_embedder, v=v, ...) +.embed_velocity <- function(x, vobj, ...) { + reducedDim(vobj, "X_target") <- as.matrix(x) + basiliskRun(env=velo.env, fun=.run_embedder, vobj=vobj, ...) } #' @importFrom reticulate import #' @importFrom zellkonverter SCE2AnnData -.run_embedder <- function(v, ...) { +.run_embedder <- function(vobj, ...) { scv <- import("scvelo") args <- list(..., basis="target", autoscale=FALSE) - adata <- SCE2AnnData(v) + adata <- SCE2AnnData(vobj) do.call(scv$tl$velocity_embedding, c(list(data=adata), args)) adata$obsm["velocity_target"] } + +#' @export +#' @rdname embedVelocity +setGeneric("embedVelocity", function(x, vobj, ...) standardGeneric("embedVelocity")) + +#' @export +#' @rdname embedVelocity +setMethod("embedVelocity", "ANY", .embed_velocity) + +#' @export +#' @rdname embedVelocity +setMethod("embedVelocity", "SingleCellExperiment", function(x, vobj, ..., use.dimred=1) { + .embed_velocity(reducedDim(x, use.dimred), vobj, ...) +}) diff --git a/R/gridVectors.R b/R/gridVectors.R index 68f033b..fccc711 100644 --- a/R/gridVectors.R +++ b/R/gridVectors.R @@ -2,13 +2,19 @@ #' #' Summarize the velocity vectors into a grid, usually for easy plotting. #' -#' @param x A low-dimensional embedding of the dataset where each cell is a row. -#' @param v A low-dimensional projection of the velocity vectors, of the same dimensions as \code{x}. +#' @param x A numeric matrix of low-dimensional coordinates, e.g., after t-SNE. +#' Alternatively, a \linkS4class{SingleCellExperiment} containing such coordinates in its \code{\link{reducedDims}}. +#' @param embedded A low-dimensional projection of the velocity vectors into the embedding of \code{x}. +#' This should be of the same dimensions as \code{x} and is typically produced by \code{\link{embedVelocity}}. #' @param resolution Integer scalar specifying the resolution of the grid, #' in terms of the number of grid intervals along each axis. #' @param scale Logical scalar indicating whether the averaged vectors should be scaled by the grid resolution. #' @param as.data.frame Logical scalar indicating whether the output should be a data.frame. #' If \code{FALSE}, a list of two matrices is returned. +#' @param ... For the generic, further arguments to pass to specific methods. +#' +#' For the SingleCellExperiment method, further arguments to pass to the ANY method. +#' @param use.dimred String or integer scalar specifying the reduced dimensions to retrieve from \code{x}. #' #' @details #' This partitions the bounding box of \code{x} into a grid with \code{resolution} units in each dimension. @@ -16,7 +22,7 @@ #' This is most obviously useful for visualization to avoid overplotting of velocity vectors. #' #' If \code{scale=TRUE}, per-block vectors are scaled so that the median vector length is comparable to the spacing between blocks. -#' This improves visualization when the scales of \code{x} and \code{v} are not immediately comparable. +#' This improves visualization when the scales of \code{x} and \code{embedded} are not immediately comparable. #' #' @return #' If \code{as.data.frame=FALSE}, a list is returned containing \code{start} and \code{end}, @@ -38,10 +44,15 @@ #' plot(tsne.results[,1], tsne.results[,2], col='grey') #' arrows(out$start.1, out$start.2, out$end.1, out$end.2, length=0.05) #' +#' @seealso +#' \code{\link{embedVelocity}}, to generate \code{embedded}. +#' @name gridVectors +NULL + #' @export #' @importFrom S4Vectors selfmatch DataFrame #' @importFrom stats median -gridVectors <- function(x, v, resolution=40, scale=TRUE, as.data.frame=TRUE) { +.grid_vectors <- function(x, embedded, resolution=40, scale=TRUE, as.data.frame=TRUE) { limits <- apply(x, 2, range) intercept <- limits[1,] delta <- (limits[2,] - intercept)/resolution @@ -53,7 +64,7 @@ gridVectors <- function(x, v, resolution=40, scale=TRUE, as.data.frame=TRUE) { tab <- as.integer(table(grp)) pos <- rowsum(x, grp)/tab - vec <- rowsum(v, grp)/tab + vec <- rowsum(embedded, grp)/tab if (scale) { d <- sqrt(rowSums(vec^2)) @@ -64,3 +75,17 @@ gridVectors <- function(x, v, resolution=40, scale=TRUE, as.data.frame=TRUE) { FUN <- if (as.data.frame) data.frame else list FUN(start=pos, end=pos + vec) } + +#' @export +#' @rdname gridVectors +setGeneric("gridVectors", function(x, embedded, ...) standardGeneric("gridVectors")) + +#' @export +#' @rdname gridVectors +setMethod("gridVectors", "ANY", .grid_vectors) + +#' @export +#' @rdname gridVectors +setMethod("gridVectors", "SingleCellExperiment", function(x, embedded, ..., use.dimred=1) { + .grid_vectors(reducedDim(x, use.dimred), embedded, ...) +}) diff --git a/man/embedVelocity.Rd b/man/embedVelocity.Rd index fd5aa64..49de801 100644 --- a/man/embedVelocity.Rd +++ b/man/embedVelocity.Rd @@ -2,17 +2,30 @@ % Please edit documentation in R/embedVelocity.R \name{embedVelocity} \alias{embedVelocity} +\alias{embedVelocity,ANY-method} +\alias{embedVelocity,SingleCellExperiment-method} \title{Project velocities onto an embedding} \usage{ -embedVelocity(x, v, ...) +embedVelocity(x, vobj, ...) + +\S4method{embedVelocity}{ANY}(x, vobj, ...) + +\S4method{embedVelocity}{SingleCellExperiment}(x, vobj, ..., use.dimred = 1) } \arguments{ -\item{x}{A numeric matrix of low-dimensional coordinates, e.g., after t-SNE.} +\item{x}{A numeric matrix of low-dimensional coordinates, e.g., after t-SNE. +Alternatively, a \linkS4class{SingleCellExperiment} containing such coordinates in its \code{\link{reducedDims}}.} -\item{v}{A \linkS4class{SingleCellExperiment} containing the output of the velocity calculations, +\item{vobj}{A \linkS4class{SingleCellExperiment} containing the output of the velocity calculations, typically after running \code{\link{scvelo}}.} -\item{...}{Further arguments to pass to the \code{velocity_embedding} Python function from \pkg{scVelo}.} +\item{...}{For the generic, further arguments to pass to specific methods. + +For the ANY method, further arguments to pass to the \code{velocity_embedding} Python function from \pkg{scVelo}. + +For the SingleCellExperiment method, further arguments to pass to the ANY method.} + +\item{use.dimred}{String or integer scalar specifying the reduced dimensions to retrieve from \code{x}.} } \value{ A numeric matrix of the same dimensions as \code{x}, containing the projected velocity vectors in that embedding. diff --git a/man/gridVectors.Rd b/man/gridVectors.Rd index 569a8a9..0502d2d 100644 --- a/man/gridVectors.Rd +++ b/man/gridVectors.Rd @@ -2,14 +2,26 @@ % Please edit documentation in R/gridVectors.R \name{gridVectors} \alias{gridVectors} +\alias{gridVectors,ANY-method} +\alias{gridVectors,SingleCellExperiment-method} \title{Summarize vectors into a grid} \usage{ -gridVectors(x, v, resolution = 40, scale = TRUE, as.data.frame = TRUE) +gridVectors(x, embedded, ...) + +\S4method{gridVectors}{ANY}(x, embedded, resolution = 40, scale = TRUE, as.data.frame = TRUE) + +\S4method{gridVectors}{SingleCellExperiment}(x, embedded, ..., use.dimred = 1) } \arguments{ -\item{x}{A low-dimensional embedding of the dataset where each cell is a row.} +\item{x}{A numeric matrix of low-dimensional coordinates, e.g., after t-SNE. +Alternatively, a \linkS4class{SingleCellExperiment} containing such coordinates in its \code{\link{reducedDims}}.} + +\item{embedded}{A low-dimensional projection of the velocity vectors into the embedding of \code{x}. +This should be of the same dimensions as \code{x} and is typically produced by \code{\link{embedVelocity}}.} -\item{v}{A low-dimensional projection of the velocity vectors, of the same dimensions as \code{x}.} +\item{...}{For the generic, further arguments to pass to specific methods. + +For the SingleCellExperiment method, further arguments to pass to the ANY method.} \item{resolution}{Integer scalar specifying the resolution of the grid, in terms of the number of grid intervals along each axis.} @@ -18,6 +30,8 @@ in terms of the number of grid intervals along each axis.} \item{as.data.frame}{Logical scalar indicating whether the output should be a data.frame. If \code{FALSE}, a list of two matrices is returned.} + +\item{use.dimred}{String or integer scalar specifying the reduced dimensions to retrieve from \code{x}.} } \value{ If \code{as.data.frame=FALSE}, a list is returned containing \code{start} and \code{end}, @@ -37,7 +51,7 @@ The locations and vectors of all cells in each block are averaged to obtain a re This is most obviously useful for visualization to avoid overplotting of velocity vectors. If \code{scale=TRUE}, per-block vectors are scaled so that the median vector length is comparable to the spacing between blocks. -This improves visualization when the scales of \code{x} and \code{v} are not immediately comparable. +This improves visualization when the scales of \code{x} and \code{embedded} are not immediately comparable. } \examples{ tsne.results <- matrix(rnorm(10000), ncol=2) @@ -49,6 +63,9 @@ out <- gridVectors(tsne.results, tsne.vectors) plot(tsne.results[,1], tsne.results[,2], col='grey') arrows(out$start.1, out$start.2, out$end.1, out$end.2, length=0.05) +} +\seealso{ +\code{\link{embedVelocity}}, to generate \code{embedded}. } \author{ Aaron Lun diff --git a/tests/testthat/test-embed.R b/tests/testthat/test-embed.R index b9109bf..0c4e0c9 100644 --- a/tests/testthat/test-embed.R +++ b/tests/testthat/test-embed.R @@ -15,6 +15,11 @@ test_that("embedVelocity works correctly", { tsne.results <- matrix(rnorm(2*ncol(out)), ncol=2) projected <- embedVelocity(tsne.results, out) expect_identical(dim(tsne.results), dim(projected)) + + # Same results inside an SCE. + reducedDim(sce1, "TSNE") <- tsne.results + projected2 <- embedVelocity(sce1, out, use.dimred="TSNE") + expect_identical(projected, projected2) }) test_that("gridVectors works correctly", { @@ -26,4 +31,9 @@ test_that("gridVectors works correctly", { out <- gridVectors(tsne.results, projected, as.data.frame=FALSE) expect_type(out, "list") + + # Same results inside an SCE. + reducedDim(sce1, "TSNE") <- tsne.results + out2 <- gridVectors(sce1, projected, use.dimred="TSNE", as.data.frame=FALSE) + expect_identical(out, out2) })