diff --git a/NAMESPACE b/NAMESPACE
index 5017918a15815477a8c8c0152abb2c45d72d8349..b97422b694085a5fc997f3dbd4527202931578d2 100644
--- a/NAMESPACE
+++ b/NAMESPACE
@@ -31,9 +31,11 @@ export(mfaSimulate)
 export(newLun2Params)
 export(newLunParams)
 export(newMFAParams)
+export(newPhenoParams)
 export(newSCDDParams)
 export(newSimpleParams)
 export(newSplatParams)
+export(phenoSimulate)
 export(scDDEstimate)
 export(scDDSimulate)
 export(setParam)
@@ -49,6 +51,7 @@ export(summariseDiff)
 exportClasses(Lun2Params)
 exportClasses(LunParams)
 exportClasses(MFAParams)
+exportClasses(PhenoParams)
 exportClasses(SCDDParams)
 exportClasses(SimpleParams)
 exportClasses(SplatParams)
diff --git a/R/pheno-simulate.R b/R/pheno-simulate.R
new file mode 100644
index 0000000000000000000000000000000000000000..70042e49966ddba3d13947cce5aaeb4e16e94b37
--- /dev/null
+++ b/R/pheno-simulate.R
@@ -0,0 +1,84 @@
+#' PhenoPath simulation
+#'
+#' Simulate counts from a pseudotime trajectory using the PhenoPath method.
+#'
+#' @param params PhenoParams 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[phenopath]{simulate_phenopath}} that takes a
+#' \code{\link{PhenoParams}}, runs the simulation then converts the
+#' output to a \code{\link[SingleCellExperiment]{SingleCellExperiment}} object.
+#' See \code{\link[phenopath]{simulate_phenopath}} and the PhenoPath 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).
+#'
+#' Paper: \url{10.1101/159913}
+#'
+#' Code: \url{https://github.com/kieranrcampbell/phenopath}
+#'
+#' @examples
+#' sim <- phenoSimulate()
+#'
+#' @export
+#' @importFrom SingleCellExperiment SingleCellExperiment
+phenoSimulate <- function(params = newPhenoParams(), verbose = TRUE, ...) {
+
+    checkmate::assertClass(params, "PhenoParams")
+    params <- setParams(params, ...)
+
+    # Set random seed
+    seed <- getParam(params, "seed")
+    set.seed(seed)
+
+    # Get the parameters we are going to use
+    nCells <- getParam(params, "nCells")
+    nGenes <- getParam(params, "nGenes")
+    n.de <- getParam(params, "n.de")
+    n.pst <- getParam(params, "n.pst")
+    n.pst.beta <- getParam(params, "n.pst.beta")
+    n.de.pst.beta <- getParam(params, "n.de.pst.beta")
+
+    if (verbose) {message("Simulating counts...")}
+    pheno.sim <- phenopath::simulate_phenopath(N = nCells,
+                                               G_de = n.de,
+                                               G_pst = n.pst,
+                                               G_pst_beta = n.pst.beta,
+                                               G_de_pst_beta = n.de.pst.beta)
+
+    if (verbose) {message("Creating final dataset...")}
+    cell.names <- paste0("Cell", seq_len(nCells))
+    gene.names <- paste0("Gene", seq_len(nGenes))
+
+    counts <- t(pheno.sim$y)
+    rownames(counts) <- gene.names
+    colnames(counts) <- cell.names
+
+    cells <- data.frame(Cell = cell.names,
+                        Covariate = pheno.sim$x,
+                        Pseudotime = pheno.sim$z)
+    rownames(cells) <- cell.names
+
+    features <- data.frame(Gene = gene.names,
+                           Alpha = pheno.sim$parameters$alpha,
+                           Lambda = pheno.sim$parameters$lambda,
+                           Beta = pheno.sim$parameters$beta,
+                           Regime = pheno.sim$parameters$regime)
+    rownames(features) <- gene.names
+
+    sim <- SingleCellExperiment(assays = list(counts = counts),
+                                rowData = features,
+                                colData = cells,
+                                metadata = list(params = params))
+
+    return(sim)
+}
diff --git a/man/PhenoParams.Rd b/man/PhenoParams.Rd
new file mode 100644
index 0000000000000000000000000000000000000000..b1e135f28ffa19cc6fa11e950bf9295a2d443f6a
--- /dev/null
+++ b/man/PhenoParams.Rd
@@ -0,0 +1,34 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/AllClasses.R
+\docType{class}
+\name{PhenoParams}
+\alias{PhenoParams}
+\alias{PhenoParams-class}
+\title{The PhenoParams class}
+\description{
+S4 class that holds parameters for the PhenoPath simulation.
+}
+\section{Parameters}{
+
+
+The PhenoPath 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{[n.de]}}{Number of genes to simulate from the differential
+    expression regime}
+    \item{\code{[n.pst]}}{Number of genes to simulate from the pseudotime
+    regime}
+    \item{\code{[n.pst.beta]}}{Number of genes to simulate from the
+    pseudotime + beta interactions regime}
+    \item{\code{[n.de.pst.beta]}}{Number of genes to simulate from the
+    differential expression + pseudotime + interactions regime}
+}
+
+The parameters not shown in brackets can be estimated from real data using
+\code{\link{phenoEstimate}}. For details of the PhenoPath simulation
+see \code{\link{phenoSimulate}}.
+}
+
diff --git a/man/newParams.Rd b/man/newParams.Rd
index de3e6b85de30198d6b9b53158119455e52e40d50..106202ac305422fbac9a9b60c681fac9d0963685 100644
--- a/man/newParams.Rd
+++ b/man/newParams.Rd
@@ -1,12 +1,13 @@
 % 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/SCDDParams-methods.R,
-%   R/SimpleParams-methods.R, R/SplatParams-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
 \name{newParams}
 \alias{newParams}
 \alias{newLun2Params}
 \alias{newLunParams}
 \alias{newMFAParams}
+\alias{newPhenoParams}
 \alias{newSCDDParams}
 \alias{newSimpleParams}
 \alias{newSplatParams}
@@ -18,6 +19,8 @@ newLunParams(...)
 
 newMFAParams(...)
 
+newPhenoParams(...)
+
 newSCDDParams(...)
 
 newSimpleParams(...)
diff --git a/man/phenoSimulate.Rd b/man/phenoSimulate.Rd
new file mode 100644
index 0000000000000000000000000000000000000000..2effe2029d624eccd684c90b8685b215704e7aa8
--- /dev/null
+++ b/man/phenoSimulate.Rd
@@ -0,0 +1,43 @@
+% Generated by roxygen2: do not edit by hand
+% Please edit documentation in R/pheno-simulate.R
+\name{phenoSimulate}
+\alias{phenoSimulate}
+\title{PhenoPath simulation}
+\usage{
+phenoSimulate(params = newPhenoParams(), verbose = TRUE, ...)
+}
+\arguments{
+\item{params}{PhenoParams 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 from a pseudotime trajectory using the PhenoPath method.
+}
+\details{
+This function is just a wrapper around
+\code{\link[phenopath]{simulate_phenopath}} that takes a
+\code{\link{PhenoParams}}, runs the simulation then converts the
+output to a \code{\link[SingleCellExperiment]{SingleCellExperiment}} object.
+See \code{\link[phenopath]{simulate_phenopath}} and the PhenoPath paper for
+more details about how the simulation works.
+}
+\examples{
+sim <- phenoSimulate()
+
+}
+\references{
+Campbell K, Yau C. Uncovering genomic trajectories with heterogeneous genetic
+and environmental backgrounds across single-cells and populations. bioRxiv
+(2017).
+
+Paper: \url{10.1101/159913}
+
+Code: \url{https://github.com/kieranrcampbell/phenopath}
+}
diff --git a/man/setParam.Rd b/man/setParam.Rd
index 9c1ba5341613eb970da5a70d097a52fc521d0f60..ce8a32ca05e800826ca81b60da703c6c538966bd 100644
--- a/man/setParam.Rd
+++ b/man/setParam.Rd
@@ -1,13 +1,14 @@
 % 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/SCDDParams-methods.R,
-%   R/SplatParams-methods.R
+%   R/LunParams-methods.R, R/Params-methods.R, R/PhenoParams-methods.R,
+%   R/SCDDParams-methods.R, R/SplatParams-methods.R
 \docType{methods}
 \name{setParam}
 \alias{setParam}
 \alias{setParam,Lun2Params-method}
 \alias{setParam,LunParams-method}
 \alias{setParam,Params-method}
+\alias{setParam,PhenoParams-method}
 \alias{setParam,SCDDParams-method}
 \alias{setParam,SplatParams-method}
 \title{Set a parameter}
@@ -20,6 +21,8 @@ setParam(object, name, value)
 
 \S4method{setParam}{Params}(object, name, value)
 
+\S4method{setParam}{PhenoParams}(object, name, value)
+
 \S4method{setParam}{SCDDParams}(object, name, value)
 
 \S4method{setParam}{SplatParams}(object, name, value)
diff --git a/tests/testthat/test-pheno-simulate.R b/tests/testthat/test-pheno-simulate.R
new file mode 100644
index 0000000000000000000000000000000000000000..f67a25e72e0ca8cd8ef82537f83a752ae23a429a
--- /dev/null
+++ b/tests/testthat/test-pheno-simulate.R
@@ -0,0 +1,5 @@
+context("PhenoPath simulation")
+
+test_that("PhenoPath simulation output is valid", {
+    expect_true(validObject(phenoSimulate()))
+})