Skip to content
Snippets Groups Projects
Commit 53246012 authored by Luke Zappia's avatar Luke Zappia
Browse files

Add single mode to splat

parent 5f8679c0
No related branches found
No related tags found
No related merge requests found
Package: splatter
Type: Package
Title: Simple Simulation of Single-cell RNA Sequencing Data
Version: 0.3.21
Date: 2016-10-11
Version: 0.3.22
Date: 2016-10-12
Author: Luke Zappia
Authors@R: as.person(c(
"Luke Zappia <luke.zappia@mcri.edu.au> [aut, cre]",
......
......@@ -13,6 +13,7 @@ export(splat)
export(splatGroups)
export(splatParams)
export(splatPaths)
export(splatSingle)
importFrom(Biobase,"assayData<-")
importFrom(Biobase,"fData<-")
importFrom(Biobase,"pData<-")
......
......@@ -15,9 +15,10 @@
#'
#' @param params splatParams object containing parameters for the simulation.
#' See \code{\link{splatParams}} for details.
#' @param method which simulation method to use. Options are "groups" which
#' produces distinct groups (eg. cell types) or "paths" which selects
#' cells from a continuous trajectory (eg. differentiation process).
#' @param method which simulation method to use. Options are "single" which
#' produces a single population, "groups" which produces distinct groups
#' (eg. cell types) or "paths" which selects cells from continuous
#' trajectories (eg. differentiation process).
#' @param verbose logical. Whether to print progress messages.
#' @param ... any additional parameter settings to override what is provided in
#' \code{params}.
......@@ -106,7 +107,8 @@
#' @importFrom Biobase fData pData pData<- assayData
#' @importFrom scater newSCESet counts
#' @export
splat <- function(params = defaultParams(), method = c("groups", "paths"),
splat <- function(params = defaultParams(),
method = c("single", "groups", "paths"),
verbose = TRUE, ...) {
method <- match.arg(method)
......@@ -126,6 +128,11 @@ splat <- function(params = defaultParams(), method = c("groups", "paths"),
nGroups <- getParams(params, "nGroups")
group.cells <- getParams(params, "groupCells")
if (nGroups == 1 && method == "groups") {
warning("nGroups is 1, switching to single mode")
method <- "single"
}
if (verbose) {message("Creating simulation object...")}
# Set up name vectors
cell.names <- paste0("Cell", 1:nCells)
......@@ -149,16 +156,20 @@ splat <- function(params = defaultParams(), method = c("groups", "paths"),
# Make groups vector which is the index of param$groupCells repeated
# params$groupCells[index] times
groups <- lapply(1:nGroups, function(i, g) {rep(i, g[i])},
g = group.cells)
groups <- unlist(groups)
pData(sim)$Group <- group.names[groups]
if (method != "single") {
groups <- lapply(1:nGroups, function(i, g) {rep(i, g[i])},
g = group.cells)
groups <- unlist(groups)
pData(sim)$Group <- group.names[groups]
}
if (verbose) {message("Simulating library sizes...")}
sim <- simLibSizes(sim, params)
if (verbose) {message("Simulating gene means...")}
sim <- simGeneMeans(sim, params)
if (method == "groups") {
if (method == "single") {
sim <- simSingleCellMeans(sim, params)
} else if (method == "groups") {
if (verbose) {message("Simulating group DE...")}
sim <- simGroupDE(sim, params)
if (verbose) {message("Simulating cell means...")}
......@@ -193,11 +204,17 @@ splat <- function(params = defaultParams(), method = c("groups", "paths"),
return(sce)
}
#' @rdname splat
#' @export
splatSingle <- function(params = defaultParams(), verbose = TRUE, ...) {
sim <- splat(params = params, method = "single", verbose = verbose, ...)
return(sim)
}
#' @rdname splat
#' @export
splatGroups <- function(params = defaultParams(), verbose = TRUE, ...) {
sim <- splat(params = params, method = "groups", verbose = verbose, ...)
return(sim)
}
......@@ -205,7 +222,6 @@ splatGroups <- function(params = defaultParams(), verbose = TRUE, ...) {
#' @export
splatPaths <- function(params = defaultParams(), verbose = TRUE, ...) {
sim <- splat(params = params, method = "paths", verbose = verbose, ...)
return(sim)
}
......@@ -346,6 +362,35 @@ simPathDE <- function(sim, params) {
return(sim)
}
#' Simulate single population cell means
#'
#' Simulate a gene by cell matrix giving the mean expression for each gene in
#' each cell.
#'
#' @param sim SCESet to add cell means to.
#' @param params splatParams object with simulation parameters.
#'
#' @return SCESet with added cell means.
#'
#' @importFrom Biobase fData pData assayData assayData<-
simSingleCellMeans <- function(sim, params) {
nCells <- getParams(params, "nCells")
cell.names <- pData(sim)$Cell
gene.names <- fData(sim)$Gene
exp.lib.sizes <- pData(sim)$ExpLibSize
cell.means.gene <- as.matrix(fData(sim)[, rep("GeneMean", nCells)])
cell.props.gene <- t(t(cell.means.gene) / colSums(cell.means.gene))
base.means.cell <- t(t(cell.props.gene) * exp.lib.sizes)
colnames(base.means.cell) <- cell.names
rownames(base.means.cell) <- gene.names
assayData(sim)$BaseCellMeans <- base.means.cell
return(sim)
}
#' Simulate group cell means
#'
#' Simulate a gene by cell matrix giving the mean expression for each gene in
......@@ -368,10 +413,6 @@ simGroupCellMeans <- function(sim, params) {
exp.lib.sizes <- pData(sim)$ExpLibSize
group.means.gene <- fData(sim)[, paste0("GeneMean", group.names)]
if (nGroups == 1) {
group.means.gene <- matrix(group.means.gene)
colnames(group.means.gene) <- "GeneMeanGroup1"
}
cell.means.gene <- as.matrix(group.means.gene[, factor(groups)])
cell.props.gene <- t(t(cell.means.gene) / colSums(cell.means.gene))
base.means.cell <- t(t(cell.props.gene) * exp.lib.sizes)
......
% Generated by roxygen2: do not edit by hand
% Please edit documentation in R/simulate.R
\name{simSingleCellMeans}
\alias{simSingleCellMeans}
\title{Simulate single population cell means}
\usage{
simSingleCellMeans(sim, params)
}
\arguments{
\item{sim}{SCESet to add cell means to.}
\item{params}{splatParams object with simulation parameters.}
}
\value{
SCESet with added cell means.
}
\description{
Simulate a gene by cell matrix giving the mean expression for each gene in
each cell.
}
......@@ -4,11 +4,14 @@
\alias{splat}
\alias{splatGroups}
\alias{splatPaths}
\alias{splatSingle}
\title{Simulate scRNA-seq data}
\usage{
splat(params = defaultParams(), method = c("groups", "paths"),
splat(params = defaultParams(), method = c("single", "groups", "paths"),
verbose = TRUE, ...)
splatSingle(params = defaultParams(), verbose = TRUE, ...)
splatGroups(params = defaultParams(), verbose = TRUE, ...)
splatPaths(params = defaultParams(), verbose = TRUE, ...)
......@@ -17,9 +20,10 @@ splatPaths(params = defaultParams(), verbose = TRUE, ...)
\item{params}{splatParams object containing parameters for the simulation.
See \code{\link{splatParams}} for details.}
\item{method}{which simulation method to use. Options are "groups" which
produces distinct groups (eg. cell types) or "paths" which selects
cells from a continuous trajectory (eg. differentiation process).}
\item{method}{which simulation method to use. Options are "single" which
produces a single population, "groups" which produces distinct groups
(eg. cell types) or "paths" which selects cells from continuous
trajectories (eg. differentiation process).}
\item{verbose}{logical. Whether to print progress messages.}
......
context("splatter simulations")
test_that("one group switches to single mode", {
expect_warning(splat(method = "groups"),
"nGroups is 1, switching to single mode")
expect_silent(splat(method = "paths", verbose = FALSE))
})
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment