diff --git a/DESCRIPTION b/DESCRIPTION
index 4be8b4b0849846355db27205e005ee484f1fc90a..af80af7822f107b89bc93841b30e293bfd0fd438 100644
--- a/DESCRIPTION
+++ b/DESCRIPTION
@@ -53,7 +53,8 @@ Suggests:
     scDD,
     scran,
     mfa,
-    phenopath
+    phenopath,
+    zinbwave
 biocViews: SingleCell, RNASeq, Transcriptomics, GeneExpression, Sequencing,
     Software
 URL: https://github.com/Oshlack/splatter
diff --git a/NAMESPACE b/NAMESPACE
index e388c4739b1229bb1dd6937937beb9313068f8d5..4cfad39e5ef56f527991f6f4713875addfdcd440 100644
--- a/NAMESPACE
+++ b/NAMESPACE
@@ -15,6 +15,8 @@ S3method(simpleEstimate,SingleCellExperiment)
 S3method(simpleEstimate,matrix)
 S3method(splatEstimate,SingleCellExperiment)
 S3method(splatEstimate,matrix)
+S3method(zinbEstimate,SingleCellExperiment)
+S3method(zinbEstimate,matrix)
 export(addGeneLengths)
 export(compareSCEs)
 export(diffSCEs)
@@ -37,6 +39,7 @@ export(newPhenoParams)
 export(newSCDDParams)
 export(newSimpleParams)
 export(newSplatParams)
+export(newZINBParams)
 export(phenoEstimate)
 export(phenoSimulate)
 export(scDDEstimate)
@@ -51,6 +54,8 @@ export(splatSimulateGroups)
 export(splatSimulatePaths)
 export(splatSimulateSingle)
 export(summariseDiff)
+export(zinbEstimate)
+export(zinbSimulate)
 exportClasses(Lun2Params)
 exportClasses(LunParams)
 exportClasses(MFAParams)
@@ -58,6 +63,7 @@ exportClasses(PhenoParams)
 exportClasses(SCDDParams)
 exportClasses(SimpleParams)
 exportClasses(SplatParams)
+exportClasses(ZINBParams)
 importFrom(BiocParallel,SerialParam)
 importFrom(BiocParallel,bplapply)
 importFrom(SingleCellExperiment,"cpm<-")
diff --git a/R/AllClasses.R b/R/AllClasses.R
index e9a131c926fd7e0cd2bedae14601b4b5a814c837..ac6cdc094f8027b1d2a36c0a87713a81cd383ac9 100644
--- a/R/AllClasses.R
+++ b/R/AllClasses.R
@@ -540,3 +540,38 @@ setClass("PhenoParams",
                    n.de.pst.beta = "numeric"),
          prototype = prototype(n.de = 2500, n.pst = 2500, n.pst.beta = 2500,
                                n.de.pst.beta = 2500))
+
+
+#' The ZINBParams class
+#'
+#' S4 class that holds parameters for the ZINB-WaVE simulation.
+#'
+#' @section Parameters:
+#'
+#' The ZINB-WaVE simulation uses the following parameters:
+#'
+#' \describe{
+#'     \item{\code{nGenes}}{The number of genes to simulate.}
+#'     \item{\code{nCells}}{The number of cells to simulate.}
+#'     \item{\code{[seed]}}{Seed to use for generating random numbers.}
+#'     \item{\code{model}}{Object describing a ZINB model.}
+#' }
+#'
+#' The majority of the parameters for this simulation are stored in a
+#' \code{\link[zinbwave]{ZinbModel}} object. Please refer to the documentation
+#' for this class and its constructor(\code{\link[zinbwave]{zinbModel}}) for
+#' details about all the parameters.
+#'
+#' The parameters not shown in brackets can be estimated from real data using
+#' \code{\link{zinbEstimate}}. For details of the ZINB-WaVE simulation
+#' see \code{\link{zinbSimulate}}.
+#'
+#' @name ZINBParams
+#' @rdname ZINBParams
+#' @aliases ZINBParams-class
+#' @exportClass ZINBParams
+setClass("ZINBParams",
+         contains = "Params",
+         slots = c(model = "ZinbModel"),
+         prototype = prototype(nGenes = 100, nCells = 50,
+                               model = zinbwave::zinbModel()))
diff --git a/R/ZINBParams-methods.R b/R/ZINBParams-methods.R
new file mode 100644
index 0000000000000000000000000000000000000000..2cddd4da2f4d507dff03f5894ea073989add3679
--- /dev/null
+++ b/R/ZINBParams-methods.R
@@ -0,0 +1,124 @@
+#' @rdname newParams
+#' @importFrom methods new
+#' @export
+newZINBParams <- function(...) {
+
+    if (!requireNamespace("zinbwave", quietly = TRUE)) {
+        stop("The ZINB-WaVE simulation requires the 'zinbwave' package.")
+    }
+
+    params <- new("ZINBParams")
+
+    params <- setParams(params, ...)
+
+    return(params)
+}
+
+setValidity("ZINBParams", function(object) {
+
+    v <- getParams(object, slotNames(object))
+
+    checks <- c(nGenes = checkmate::checkInt(v$nGenes, lower = 1),
+                nCells = checkmate::checkInt(v$nCells, lower = 1),
+                seed = checkmate::checkInt(v$seed, lower = 0),
+                model = checkmate::checkClass(v$model, "ZinbModel"),
+                model_valid = validObject(v$model, test = TRUE))
+
+    if (all(checks == TRUE)) {
+        valid <- TRUE
+    } else {
+        valid <- checks[checks != TRUE]
+        valid <- paste(names(valid), valid, sep = ": ")
+    }
+
+    return(valid)
+})
+
+#' @rdname setParam
+setMethod("setParam", "ZINBParams", function(object, name, value) {
+    checkmate::assertString(name)
+
+    if (name %in% c("nGenes", "nCells")) {
+        stop(name, " cannot be set directly, set model instead")
+    }
+
+    if (name == "model") {
+        object <- setParamUnchecked(object, "nGenes",
+                                    zinbwave::nFeatures(value))
+        object <- setParamUnchecked(object, "nCells",
+                                    zinbwave::nSamples(value))
+    }
+
+    object <- callNextMethod()
+
+    return(object)
+})
+
+setMethod("show", "ZINBParams", function(object) {
+
+    pp <- list("Design:"         = c("(Samples)"       = "X",
+                                     "(Genes)"         = "V"),
+               "Offsets:"        = c("(Mu)"            = "O_mu",
+                                     "(Pi)"            = "O_pi"),
+               "Indices:"        = c("(Sample Mu)"     = "which_X_mu",
+                                     "(Gene Mu)"       = "which_V_mu",
+                                     "(Sample Pi)"     = "which_X_pi",
+                                     "(Gene Pi)"       = "which_V_pi"),
+               "Intercepts:"     = c("(Sample Mu)"     = "X_mu_intercept",
+                                     "(Gene Mu)"       = "V_mu_intercept",
+                                     "(Sample Pi)"     = "X_pi_intercept",
+                                     "(Gene Pi)"       = "V_pi_intercept"),
+               "Latent factors:" = c("(W)"             = "W"),
+               "Coefficients:"   = c("(Sample Mu)"     = "beta_mu",
+                                     "(Gene Mu)"       = "gamma_mu",
+                                     "(Latent Mu)"     = "alpha_mu",
+                                     "(Sample Pi)"     = "beta_pi",
+                                     "(Gene Pi)"       = "gamma_pi",
+                                     "(Latent Pi)"     = "alpha_pi"),
+               "Regularisation:" = c("(Sample Mu)"     = "epsilon_beta_mu",
+                                     "(Gene Mu)"       = "epsilon_gamma_mu",
+                                     "(Sample Pi)"     = "epsilon_beta_pi",
+                                     "(Gene Pi)"       = "epsilon_gamma_pi",
+                                     "(Latent)"        = "epsilon_W",
+                                     "(Latent coeffs)" = "epsilon_alpha",
+                                     "(Zeta)"          = "epsilon_zeta",
+                                     "(Logit)"         = "epsilon_min_logit"))
+
+    callNextMethod()
+
+    model <- getParam(object, "model")
+    cat("Model:", "\n")
+    cat("ZinbModel with", zinbwave::nFeatures(model), "features,",
+        zinbwave::nSamples(model), "samples,", zinbwave::nFactors(model),
+        "latent factors and", zinbwave::nParams(model), "parameters", "\n\n")
+
+    default <- zinbwave::zinbModel()
+    for (category in names(pp)) {
+        parameters <- pp[[category]]
+        values <- lapply(parameters, function(x) {slot(model, x)})
+        short.values <- sapply(values, function(x) {
+            if ("matrix" %in% class(x)) {
+                if (nrow(x) == 1) {
+                    paste0(paste(head(x[1, ], n = 4), collapse = ", "), ",...")
+                } else if (ncol(x) == 1) {
+                    paste0(paste(head(x[, 1], n = 4), collapse = ", "), ",...")
+                } else {
+                    paste(nrow(x), "x", ncol(x), "matrix")
+                }
+            } else if (length(x) > 4) {
+                paste0(paste(head(x, n = 4), collapse = ", "), ",...")
+            } else {
+                paste(x, collapse = ", ")
+            }
+        })
+        values <- sapply(values, paste, collapse = ", ")
+        default.values <- lapply(parameters, function(x) {slot(default, x)})
+        default.values <- sapply(default.values, paste, collapse = ", ")
+        not.default <- values != default.values
+        names(short.values)[not.default] <- toupper(names(values[not.default]))
+        cat("Model", category, "\n")
+        print(noquote(short.values), print.gap = 2)
+        cat("\n")
+    }
+
+})
diff --git a/R/zinb-estimate.R b/R/zinb-estimate.R
new file mode 100644
index 0000000000000000000000000000000000000000..217283c8e91058b5199ee67462b5437b0ecbce0b
--- /dev/null
+++ b/R/zinb-estimate.R
@@ -0,0 +1,95 @@
+#' Estimate ZINB-WaVE simulation parameters
+#'
+#' Estimate simulation parameters for the ZINB-WaVE simulation from a real
+#' dataset.
+#'
+#' @param counts either a counts matrix or a SingleCellExperiment object
+#'        containing count data to estimate parameters from.
+#' @param design.samples design matrix of sample-level covariates.
+#' @param design.genes design matrix of gene-level covariates.
+#' @param common.disp logical. Whether or not a single dispersion for all
+#'        features is estimated.
+#' @param iter.init number of iterations to use for initalization.
+#' @param iter.opt number of iterations to use for optimization.
+#' @param stop.opt stopping criterion for optimization.
+#' @param params ZINBParams object to store estimated values in.
+#' @param verbose logical. Whether to print progress messages.
+#' @param BPPARAM A \code{\link[BiocParallel]{BiocParallelParam}} instance
+#'        giving the parallel back-end to be used. Default is
+#'        \code{\link[BiocParallel]{SerialParam}} which uses a single core.
+#' @param ... additional arguments passes to \code{\link[zinbwave]{zinbFit}}.
+#'
+#' @details
+#' The function is a wrapper around \code{\link[zinbwave]{zinbFit}} that takes
+#' the fitted model and inserts it into a \code{\link{ZINBParams}} object. See
+#' \code{\link{ZINBParams}} for more details on the parameters and
+#' \code{\link[zinbwave]{zinbFit}} for details of the estimation procedure.
+#'
+#' @return ZINBParams object containing the estimated parameters.
+#'
+#' @examples
+#' data("sc_example_counts")
+#' params <- zinbEstimate(sc_example_counts)
+#' params
+#' @importFrom BiocParallel SerialParam
+#' @export
+zinbEstimate <- function(counts, design.samples = NULL, design.genes = NULL,
+                         common.disp = TRUE, iter.init = 2, iter.opt = 25,
+                         stop.opt = 1e-04, params = newZINBParams(),
+                         verbose = TRUE, BPPARAM = SerialParam(), ...) {
+    UseMethod("zinbEstimate")
+}
+
+#' @rdname zinbEstimate
+#' @export
+zinbEstimate.SingleCellExperiment <- function(counts, design.samples = NULL,
+                                              design.genes = NULL,
+                                              common.disp = TRUE,
+                                              iter.init = 2, iter.opt = 25,
+                                              stop.opt = 1e-04,
+                                              params = newZINBParams(),
+                                              verbose = TRUE,
+                                              BPPARAM = SerialParam(), ...) {
+    counts <- BiocGenerics::counts(counts)
+    zinbEstimate(counts, design.samples, design.genes, common.disp,
+                 iter.init, iter.opt, stop.opt, params, verbose, BPPARAM, ...)
+}
+
+#' @rdname zinbEstimate
+#' @export
+zinbEstimate.matrix <- function(counts, design.samples = NULL,
+                                design.genes = NULL, common.disp = TRUE,
+                                iter.init = 2, iter.opt = 25, stop.opt = 1e-04,
+                                params = newZINBParams(), verbose = TRUE,
+                                BPPARAM = SerialParam(), ...) {
+
+    checkmate::assertClass(params, "ZINBParams")
+
+    if (verbose) {message("Removing all zero genes...")}
+    counts <- counts[rowSums(counts) > 0, ]
+
+    args.list <- list(Y = counts,
+                      commondispersion = common.disp,
+                      verbose = verbose,
+                      nb.repeat.initialize = iter.init,
+                      maxiter.optimize = iter.opt,
+                      stop.epsilon.optimize = stop.opt,
+                      BPPARAM = BPPARAM)
+
+    if (!is.null(design.samples)) {
+        args.list$X <- design.samples
+    }
+
+    if (!is.null(design.genes)) {
+        args.list$V <- design.genes
+    }
+
+    args.list <- c(args.list, list(...))
+
+    if (verbose) {message("Fitting model...")}
+    model <- do.call(zinbwave::zinbFit, args.list)
+
+    params <- setParams(params, model = model)
+
+    return(params)
+}
diff --git a/R/zinb-simulate.R b/R/zinb-simulate.R
new file mode 100644
index 0000000000000000000000000000000000000000..92e93c5aae1c4608da315c1dbecb44c546965804
--- /dev/null
+++ b/R/zinb-simulate.R
@@ -0,0 +1,74 @@
+#' ZINB-WaVE simulation
+#'
+#' Simulate counts using the ZINB-WaVE method.
+#'
+#' @param params ZINBParams object containing simulation parameters.
+#' @param verbose logical. Whether to print progress messages
+#' @param ... any additional parameter settings to override what is provided in
+#'        \code{params}.
+#'
+#' @details
+#' This function is just a wrapper around \code{\link[zinbwave]{zinbSim}} that
+#' takes a \code{\link{ZINBParams}}, runs the simulation then converts the
+#' output to a \code{\link[SingleCellExperiment]{SingleCellExperiment}} object.
+#' See \code{\link[zinbwave]{zinbSim}} and the ZINB-WaVE paper for
+#' more details about how the simulation works.
+#'
+#' @return SingleCellExperiment containing simulated counts
+#'
+#' @references
+#' Campbell K, Yau C. Uncovering genomic trajectories with heterogeneous genetic
+#' and environmental backgrounds across single-cells and populations. bioRxiv
+#' (2017).
+#'
+#' Risso D, Perraudeau F, Gribkova S, Dudoit S, Vert J-P. ZINB-WaVE: A general
+#' and flexible method for signal extraction from single-cell RNA-seq data
+#' bioRxiv (2017).
+#'
+#' Paper: \url{10.1101/125112}
+#'
+#' Code: \url{https://github.com/drisso/zinbwave}
+#'
+#' @examples
+#' sim <- zinbSimulate()
+#'
+#' @export
+#' @importFrom SingleCellExperiment SingleCellExperiment
+zinbSimulate <- function(params = newZINBParams(), verbose = TRUE, ...) {
+
+    checkmate::assertClass(params, "ZINBParams")
+    params <- setParams(params, ...)
+
+    # Get the parameters we are going to use
+    nCells <- getParam(params, "nCells")
+    nGenes <- getParam(params, "nGenes")
+    model <- getParam(params, "model")
+    seed <- getParam(params, "seed")
+
+    if (verbose) {message("Simulating counts...")}
+    zinb.sim <- zinbwave::zinbSim(model, seed)
+
+    if (verbose) {message("Creating final dataset...")}
+    cell.names <- paste0("Cell", seq_len(nCells))
+    gene.names <- paste0("Gene", seq_len(nGenes))
+
+    for (item in c("counts", "dataNB", "dataDropouts")) {
+        rownames(zinb.sim[[item]]) <- gene.names
+        colnames(zinb.sim[[item]]) <- cell.names
+    }
+
+    cells <- data.frame(Cell = cell.names)
+    rownames(cells) <- cell.names
+
+    features <- data.frame(Gene = gene.names)
+    rownames(features) <- gene.names
+
+    sim <- SingleCellExperiment(assays = list(counts = zinb.sim$counts,
+                                              TrueCounts = zinb.sim$dataNB,
+                                              Dropouts = zinb.sim$dataDropouts),
+                                rowData = features,
+                                colData = cells,
+                                metadata = list(params = params))
+
+    return(sim)
+}
diff --git a/man/ZINBParams.Rd b/man/ZINBParams.Rd
new file mode 100644
index 0000000000000000000000000000000000000000..e7fc22a1cdcf7a8627bfbb2c8edd5e0c095fb1b3
--- /dev/null
+++ b/man/ZINBParams.Rd
@@ -0,0 +1,32 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/AllClasses.R
+\docType{class}
+\name{ZINBParams}
+\alias{ZINBParams}
+\alias{ZINBParams-class}
+\title{The ZINBParams class}
+\description{
+S4 class that holds parameters for the ZINB-WaVE simulation.
+}
+\section{Parameters}{
+
+
+The ZINB-WaVE simulation uses the following parameters:
+
+\describe{
+    \item{\code{nGenes}}{The number of genes to simulate.}
+    \item{\code{nCells}}{The number of cells to simulate.}
+    \item{\code{[seed]}}{Seed to use for generating random numbers.}
+    \item{\code{model}}{Object describing a ZINB model.}
+}
+
+The majority of the parameters for this simulation are stored in a
+\code{\link[zinbwave]{ZinbModel}} object. Please refer to the documentation
+for this class and its constructor(\code{\link[zinbwave]{zinbModel}}) for
+details about all the parameters.
+
+The parameters not shown in brackets can be estimated from real data using
+\code{\link{zinbEstimate}}. For details of the ZINB-WaVE simulation
+see \code{\link{zinbSimulate}}.
+}
+
diff --git a/man/newParams.Rd b/man/newParams.Rd
index 106202ac305422fbac9a9b60c681fac9d0963685..bac473177a58a9f19db7e21ae2e75ce97a887d34 100644
--- a/man/newParams.Rd
+++ b/man/newParams.Rd
@@ -1,7 +1,8 @@
 % Generated by roxygen2: do not edit by hand
 % Please edit documentation in R/AllGenerics.R, R/Lun2Params-methods.R,
 %   R/LunParams-methods.R, R/MFAParams-methods.R, R/PhenoParams-methods.R,
-%   R/SCDDParams-methods.R, R/SimpleParams-methods.R, R/SplatParams-methods.R
+%   R/SCDDParams-methods.R, R/SimpleParams-methods.R, R/SplatParams-methods.R,
+%   R/ZINBParams-methods.R
 \name{newParams}
 \alias{newParams}
 \alias{newLun2Params}
@@ -11,6 +12,7 @@
 \alias{newSCDDParams}
 \alias{newSimpleParams}
 \alias{newSplatParams}
+\alias{newZINBParams}
 \title{New Params}
 \usage{
 newLun2Params(...)
@@ -26,6 +28,8 @@ newSCDDParams(...)
 newSimpleParams(...)
 
 newSplatParams(...)
+
+newZINBParams(...)
 }
 \arguments{
 \item{...}{additional parameters passed to \code{\link{setParams}}.}
diff --git a/man/setParam.Rd b/man/setParam.Rd
index ce8a32ca05e800826ca81b60da703c6c538966bd..3958a594224ee6318da204d9d4afa71715a82727 100644
--- a/man/setParam.Rd
+++ b/man/setParam.Rd
@@ -1,7 +1,7 @@
 % Generated by roxygen2: do not edit by hand
 % Please edit documentation in R/AllGenerics.R, R/Lun2Params-methods.R,
 %   R/LunParams-methods.R, R/Params-methods.R, R/PhenoParams-methods.R,
-%   R/SCDDParams-methods.R, R/SplatParams-methods.R
+%   R/SCDDParams-methods.R, R/SplatParams-methods.R, R/ZINBParams-methods.R
 \docType{methods}
 \name{setParam}
 \alias{setParam}
@@ -11,6 +11,7 @@
 \alias{setParam,PhenoParams-method}
 \alias{setParam,SCDDParams-method}
 \alias{setParam,SplatParams-method}
+\alias{setParam,ZINBParams-method}
 \title{Set a parameter}
 \usage{
 setParam(object, name, value)
@@ -26,6 +27,8 @@ setParam(object, name, value)
 \S4method{setParam}{SCDDParams}(object, name, value)
 
 \S4method{setParam}{SplatParams}(object, name, value)
+
+\S4method{setParam}{ZINBParams}(object, name, value)
 }
 \arguments{
 \item{object}{object to set parameter in.}
diff --git a/man/zinbEstimate.Rd b/man/zinbEstimate.Rd
new file mode 100644
index 0000000000000000000000000000000000000000..84956c9d532364f5793ff5e530328f17abcc2c7c
--- /dev/null
+++ b/man/zinbEstimate.Rd
@@ -0,0 +1,67 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/zinb-estimate.R
+\name{zinbEstimate}
+\alias{zinbEstimate}
+\alias{zinbEstimate.SingleCellExperiment}
+\alias{zinbEstimate.matrix}
+\title{Estimate ZINB-WaVE simulation parameters}
+\usage{
+zinbEstimate(counts, design.samples = NULL, design.genes = NULL,
+  common.disp = TRUE, iter.init = 2, iter.opt = 25, stop.opt = 1e-04,
+  params = newZINBParams(), verbose = TRUE, BPPARAM = SerialParam(), ...)
+
+\method{zinbEstimate}{SingleCellExperiment}(counts, design.samples = NULL,
+  design.genes = NULL, common.disp = TRUE, iter.init = 2, iter.opt = 25,
+  stop.opt = 1e-04, params = newZINBParams(), verbose = TRUE,
+  BPPARAM = SerialParam(), ...)
+
+\method{zinbEstimate}{matrix}(counts, design.samples = NULL,
+  design.genes = NULL, common.disp = TRUE, iter.init = 2, iter.opt = 25,
+  stop.opt = 1e-04, params = newZINBParams(), verbose = TRUE,
+  BPPARAM = SerialParam(), ...)
+}
+\arguments{
+\item{counts}{either a counts matrix or a SingleCellExperiment object
+containing count data to estimate parameters from.}
+
+\item{design.samples}{design matrix of sample-level covariates.}
+
+\item{design.genes}{design matrix of gene-level covariates.}
+
+\item{common.disp}{logical. Whether or not a single dispersion for all
+features is estimated.}
+
+\item{iter.init}{number of iterations to use for initalization.}
+
+\item{iter.opt}{number of iterations to use for optimization.}
+
+\item{stop.opt}{stopping criterion for optimization.}
+
+\item{params}{ZINBParams object to store estimated values in.}
+
+\item{verbose}{logical. Whether to print progress messages.}
+
+\item{BPPARAM}{A \code{\link[BiocParallel]{BiocParallelParam}} instance
+giving the parallel back-end to be used. Default is
+\code{\link[BiocParallel]{SerialParam}} which uses a single core.}
+
+\item{...}{additional arguments passes to \code{\link[zinbwave]{zinbFit}}.}
+}
+\value{
+ZINBParams object containing the estimated parameters.
+}
+\description{
+Estimate simulation parameters for the ZINB-WaVE simulation from a real
+dataset.
+}
+\details{
+The function is a wrapper around \code{\link[zinbwave]{zinbFit}} that takes
+the fitted model and inserts it into a \code{\link{ZINBParams}} object. See
+\code{\link{ZINBParams}} for more details on the parameters and
+\code{\link[zinbwave]{zinbFit}} for details of the estimation procedure.
+}
+\examples{
+data("sc_example_counts")
+params <- zinbEstimate(sc_example_counts)
+params
+}
diff --git a/man/zinbSimulate.Rd b/man/zinbSimulate.Rd
new file mode 100644
index 0000000000000000000000000000000000000000..17c6611ee19c0ce88a5412947f5945b2737f3738
--- /dev/null
+++ b/man/zinbSimulate.Rd
@@ -0,0 +1,46 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/zinb-simulate.R
+\name{zinbSimulate}
+\alias{zinbSimulate}
+\title{ZINB-WaVE simulation}
+\usage{
+zinbSimulate(params = newZINBParams(), verbose = TRUE, ...)
+}
+\arguments{
+\item{params}{ZINBParams object containing simulation parameters.}
+
+\item{verbose}{logical. Whether to print progress messages}
+
+\item{...}{any additional parameter settings to override what is provided in
+\code{params}.}
+}
+\value{
+SingleCellExperiment containing simulated counts
+}
+\description{
+Simulate counts using the ZINB-WaVE method.
+}
+\details{
+This function is just a wrapper around \code{\link[zinbwave]{zinbSim}} that
+takes a \code{\link{ZINBParams}}, runs the simulation then converts the
+output to a \code{\link[SingleCellExperiment]{SingleCellExperiment}} object.
+See \code{\link[zinbwave]{zinbSim}} and the ZINB-WaVE paper for
+more details about how the simulation works.
+}
+\examples{
+sim <- zinbSimulate()
+
+}
+\references{
+Campbell K, Yau C. Uncovering genomic trajectories with heterogeneous genetic
+and environmental backgrounds across single-cells and populations. bioRxiv
+(2017).
+
+Risso D, Perraudeau F, Gribkova S, Dudoit S, Vert J-P. ZINB-WaVE: A general
+and flexible method for signal extraction from single-cell RNA-seq data
+bioRxiv (2017).
+
+Paper: \url{10.1101/125112}
+
+Code: \url{https://github.com/drisso/zinbwave}
+}
diff --git a/tests/testthat/test-PhenoParams.R b/tests/testthat/test-PhenoParams.R
index 535724ad0b4ddf631c7e9e10a4430a9be9d092d1..d5a252daba9f5f300840bc949f8a199b79d6e7a1 100644
--- a/tests/testthat/test-PhenoParams.R
+++ b/tests/testthat/test-PhenoParams.R
@@ -1,7 +1,7 @@
 context("PhenoParams")
 
 test_that("constructor is valid", {
-    expect_true(validObject(newSCDDParams()))
+    expect_true(validObject(newPhenoParams()))
 })
 
 test_that("nGenes checks work", {
diff --git a/tests/testthat/test-ZINBParams.R b/tests/testthat/test-ZINBParams.R
new file mode 100644
index 0000000000000000000000000000000000000000..4a4cfdece7ffc70332b918b9e87bd0c1c885eb8c
--- /dev/null
+++ b/tests/testthat/test-ZINBParams.R
@@ -0,0 +1,13 @@
+context("ZINBParams")
+
+test_that("constructor is valid", {
+    expect_true(validObject(newZINBParams()))
+})
+
+test_that("nGenes checks work", {
+    params <- newZINBParams()
+    expect_error(setParam(params, "nGenes", 1),
+                 "nGenes cannot be set directly")
+    expect_error(setParam(params, "nCells", 1),
+                 "nGenes cannot be set directly")
+})
diff --git a/tests/testthat/test-zinb-simulate.R b/tests/testthat/test-zinb-simulate.R
new file mode 100644
index 0000000000000000000000000000000000000000..80b2b1c5bcd8062be039f7a2260458920353d6de
--- /dev/null
+++ b/tests/testthat/test-zinb-simulate.R
@@ -0,0 +1,5 @@
+context("ZINB-WaVE simulation")
+
+test_that("ZINB-WaVE simulation output is valid", {
+    expect_true(validObject(zinbSimulate()))
+})