Welcome to Splatter! Splatter is an R package for the simple simulation of single-cell RNA sequencing data. This vignette gives an overview and introduction to Splatter’s functionality.
Splatter can be installed from Bioconductor:
source("https://bioconductor.org/biocLite.R")
biocLite("splatter")
To install the most recent development version from Github use:
biocLite("Oshlack/splatter", dependencies = TRUE,
build_vignettes = TRUE)
Assuming you already have a matrix of count data similar to that you wish to simulate there are two simple steps to creating a simulated data set with Splatter. Here is an example using the example dataset in the scater
package:
# Load package
library(splatter)
## Loading required package: scater
## Loading required package: Biobase
## Loading required package: BiocGenerics
## Loading required package: parallel
##
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:parallel':
##
## clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
## clusterExport, clusterMap, parApply, parCapply, parLapply,
## parLapplyLB, parRapply, parSapply, parSapplyLB
## The following objects are masked from 'package:stats':
##
## IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
##
## Filter, Find, Map, Position, Reduce, anyDuplicated, append,
## as.data.frame, cbind, colMeans, colSums, colnames, do.call,
## duplicated, eval, evalq, get, grep, grepl, intersect,
## is.unsorted, lapply, lengths, mapply, match, mget, order,
## paste, pmax, pmax.int, pmin, pmin.int, rank, rbind, rowMeans,
## rowSums, rownames, sapply, setdiff, sort, table, tapply,
## union, unique, unsplit, which, which.max, which.min
## Welcome to Bioconductor
##
## Vignettes contain introductory material; view with
## 'browseVignettes()'. To cite Bioconductor, see
## 'citation("Biobase")', and for packages 'citation("pkgname")'.
## Loading required package: ggplot2
##
## Attaching package: 'scater'
## The following object is masked from 'package:stats':
##
## filter
# Load example data
data("sc_example_counts")
# Estimate parameters from example data
params <- splatEstimate(sc_example_counts)
# Simulate data using estimated parameters
sim <- splatSimulate(params, dropout.present = FALSE)
## Getting parameters...
## Creating simulation object...
## Simulating library sizes...
## Simulating gene means...
## Simulating BCV...
## Simulating counts..
## Simulating dropout (if needed)...
## Creating final SCESet...
## Done!
These steps will be explained in detail in the following sections but briefly the first step takes a dataset and estimates simulation parameters from it and the second step takes those parameters and simulates a new dataset.
Before we look at how we estimate parameters let’s first look at how Splatter simulates data and what those parameters are. We use the term ‘Splat’ to refer to the Splatter’s own simulation and differentiate it from the package itself. The core of the Splat model is a gamma-Poisson distribution used to generate a gene by cell matrix of counts. Mean expression levels for each gene are simulated from a gamma distribution and the Biological Coefficient of Variation is used to enforce a mean-variance trend before counts are simulated from a Poisson distribution. Splat also allows you to simulate expression outlier genes (genes with mean expression outside the gamma distribution) and dropout (random knock out of counts based on mean expression). Each cell is given an expected library size (simulated from a log-normal distribution) that makes it easier to match to a given dataset.
Splat can also simulate differential expression between groups of different types of cells or differentiation paths between different cells types where expression changes in a continuous way. These are described further in the simulating counts section.
The parameters required for the Splat simulation are briefly described here:
nGenes
- The number of genes to simulate.nCells
- The number of cells to simulate.nGroups
- The number of groups or paths to simulate.groupCells
- The number of cells in each group/path.seed
- Seed to use for generating random numbers.mean.shape
- Shape parameter for the mean gamma distribution.mean.rate
- Rate parameter for the mean gamma distribution.lib.loc
- Location (meanlog) parameter for the library size log-normal distribution.lib.scale
- Scale (sdlog) parameter for the library size log-normal distribution.out.prob
- Probability that a gene is an expression outlier.out.facLoc
- Location (meanlog) parameter for the expression outlier factor log-normal distribution.out.facScale
- Scale (sdlog) parameter for the expression outlier factor log-normal distribution.de.prob
- Probability that a gene is differentially expressed in each group or path.de.loProb
- Probability that a differentially expressed gene is down-regulated.de.facLoc
- Location (meanlog) parameter for the differential expression factor log-normal distribution.de.facScale
- Scale (sdlog) parameter for the differential expression factor log-normal distribution.bcv.common
- Underlying common dispersion across all genes.bcv.df
- Degrees of Freedom for the BCV inverse chi-squared distribution.dropout.present
- Logical. Whether to simulate dropout.dropout.mid
- Midpoint parameter for the dropout logistic function.dropout.shape
- Shape parameter for the dropout logistic function.path.from
- Vector giving the originating point of each path.path.length
- Vector giving the number of steps to simulate along each path.path.skew
- Vector giving the skew of each path.path.nonlinearProb
- Probability that a gene changes expression in a non-linear way along the differentiation path.path.sigmaFac
- Sigma factor for non-linear gene paths.While this may look like a lot of parameters Splatter attempts to make it easy for the user, both by providing sensible defaults and making it easy to estimate many of the parameters from real data. For more details on the parameters see ?SplatParams
.
SplatParams
objectAll the parameters for the Splat simulation are stored in a SplatParams
object. Let’s create a new one and see what it looks like.
params <- newSplatParams()
params
## A Params object of class SplatParams
## Parameters can be (estimable) or [not estimable], 'Default' or 'NOT DEFAULT'.
##
## Global:
## (Genes) (Cells) [Seed]
## 10000 100 945700
##
## 23 additional parameters
##
## Groups:
## [Groups] [Group Cells]
## 1 100
##
## Mean:
## (Rate) (Shape)
## 0.3 0.6
##
## Library size:
## (Location) (Scale)
## 11 0.2
##
## Exprs outliers:
## (Probability) (Location) (Scale)
## 0.05 4 0.5
##
## Diff expr:
## [Probability] [Down Prob] [Location] [Scale]
## 0.1 0.5 0.1 0.4
##
## BCV:
## (Common Disp) (DoF)
## 0.1 60
##
## Dropout:
## [Present] (Midpoint) (Shape)
## FALSE 0 -1
##
## Paths:
## [From] [Length] [Skew] [Non-linear]
## 0 100 0.5 0.1
## [Sigma Factor]
## 0.8
As well as telling us what type of object we have (“A Params
object of class SplatParams
”) and showing us the values of the parameter this output gives us some extra information. We can see which parameters can be estimated by the splatEstimate
function (those in parentheses), which can’t be estimated (those in brackets) and which have been changed from their default values (those in ALL CAPS).
If we want to look at a particular parameter, for example the number of genes to simulate, we can extract it using the getParam
function:
getParam(params, "nGenes")
## [1] 10000
Alternatively, to give a parameter a new value we can use the setParam
function:
params <- setParam(params, "nGenes", 5000)
getParam(params, "nGenes")
## [1] 5000
If we want to extract multiple parameters (as a list) or set multiple parameters we can use the getParams
or setParams
functions:
# Set multiple parameters at once (using a list)
params <- setParams(params, update = list(nGenes = 8000, mean.rate = 0.5))
# Extract multiple parameters as a list
getParams(params, c("nGenes", "mean.rate", "mean.shape"))
## $nGenes
## [1] 8000
##
## $mean.rate
## [1] 0.5
##
## $mean.shape
## [1] 0.6
# Set multiple parameters at once (using additional arguments)
params <- setParams(params, mean.shape = 0.5, de.prob = 0.2)
params
## A Params object of class SplatParams
## Parameters can be (estimable) or [not estimable], 'Default' or 'NOT DEFAULT'.
##
## Global:
## (GENES) (Cells) [Seed]
## 8000 100 945700
##
## 23 additional parameters
##
## Groups:
## [Groups] [Group Cells]
## 1 100
##
## Mean:
## (RATE) (SHAPE)
## 0.5 0.5
##
## Library size:
## (Location) (Scale)
## 11 0.2
##
## Exprs outliers:
## (Probability) (Location) (Scale)
## 0.05 4 0.5
##
## Diff expr:
## [PROBABILITY] [Down Prob] [Location] [Scale]
## 0.2 0.5 0.1 0.4
##
## BCV:
## (Common Disp) (DoF)
## 0.1 60
##
## Dropout:
## [Present] (Midpoint) (Shape)
## FALSE 0 -1
##
## Paths:
## [From] [Length] [Skew] [Non-linear]
## 0 100 0.5 0.1
## [Sigma Factor]
## 0.8
The parameters with have changed are now shown in ALL CAPS to indicate that they been changed form the default.
We can also set parameters directly when we call newSplatParams
:
params <- newSplatParams(lib.loc = 12, lib.scale = 0.6)
getParams(params, c("lib.loc", "lib.scale"))
## $lib.loc
## [1] 12
##
## $lib.scale
## [1] 0.6
Splat allows you to estimate many of it’s parameters from a data set containing counts using the splatEstimate
function.
# Check that sc_example counts is an integer matrix
class(sc_example_counts)
## [1] "matrix"
typeof(sc_example_counts)
## [1] "integer"
# Check the dimensions, each row is a gene, each column is a cell
dim(sc_example_counts)
## [1] 2000 40
# Show the first few entries
sc_example_counts[1:5, 1:5]
## Cell_001 Cell_002 Cell_003 Cell_004 Cell_005
## Gene_0001 0 123 2 0 0
## Gene_0002 575 65 3 1561 2311
## Gene_0003 0 0 0 0 1213
## Gene_0004 0 1 0 0 0
## Gene_0005 0 0 11 0 0
params <- splatEstimate(sc_example_counts)
Here we estimated parameters from a counts matrix but splatEstimate
can also take an SCESet
object from the scater
package. The estimation process has the following steps:
estimateDisp
function from the edgeR
package.For more details of the estimation procedures see ?splatEstimate
.
Once we have a set of parameters we are happy with we can use splatSimulate
to simulate counts. If we want to make small adjustments to the parameters we can provide them as additional arguments, alternatively if we don’t supply any parameters the defaults will be used:
sim <- splatSimulate(params, nGenes = 1000, dropout.present = FALSE)
## Getting parameters...
## Creating simulation object...
## Simulating library sizes...
## Simulating gene means...
## Simulating BCV...
## Simulating counts..
## Simulating dropout (if needed)...
## Creating final SCESet...
## Done!
sim
## SCESet (storageMode: lockedEnvironment)
## assayData: 1000 features, 40 samples
## element names: BCV, BaseCellMeans, CellMeans, TrueCounts, counts, exprs
## protocolData: none
## phenoData
## sampleNames: Cell1 Cell2 ... Cell40 (40 total)
## varLabels: Cell ExpLibSize
## varMetadata: labelDescription
## featureData
## featureNames: Gene1 Gene2 ... Gene1000 (1000 total)
## fvarLabels: Gene BaseGeneMean OutlierFactor GeneMean
## fvarMetadata: labelDescription
## experimentData: use 'experimentData(object)'
## Annotation:
Looking at the output of splatSimulate
we can see that sim
is an SCESet
object with 1000 features (genes) and 40 samples (cells). The main part of this object is a features by samples matrix containing the simulated counts (accessed using counts
), although it can also hold other expression measures such as FPKM or TPM. Additionaly an SCESet
contains phenotype information about each cell (accessed using pData
) and feature information about each gene (accessed using fData
). Splatter uses these slots to store information about the intermediate values of the simulation.
# Access the counts
counts(sim)[1:5, 1:5]
## Cell1 Cell2 Cell3 Cell4 Cell5
## Gene1 0 2 0 5 29
## Gene2 0 1421 0 350 1
## Gene3 5 0 163 0 8515
## Gene4 0 3 6 0 0
## Gene5 0 19 0 0 0
# Information about genes
head(fData(sim))
## Gene BaseGeneMean OutlierFactor GeneMean
## Gene1 Gene1 1.126824 1 1.126824
## Gene2 Gene2 81.914008 1 81.914008
## Gene3 Gene3 356.083786 1 356.083786
## Gene4 Gene4 1.223045 1 1.223045
## Gene5 Gene5 3.570305 1 3.570305
## Gene6 Gene6 1.376143 1 1.376143
# Information about cells
head(pData(sim))
## Cell ExpLibSize
## Cell1 Cell1 53647.02
## Cell2 Cell2 229422.33
## Cell3 Cell3 237554.11
## Cell4 Cell4 868648.65
## Cell5 Cell5 452927.73
## Cell6 Cell6 400472.04
# Gene by cell matrices
names(assayData(sim))
## [1] "TrueCounts" "BaseCellMeans" "counts" "BCV"
## [5] "CellMeans" "exprs"
# Example of cell means matrix
get_exprs(sim, "CellMeans")[1:5, 1:5]
## Cell1 Cell2 Cell3 Cell4 Cell5
## Gene1 2.461662e-38 1.952075e+00 2.936602e-10 3.668324e+00 2.845669e+01
## Gene2 1.424940e-08 1.454941e+03 3.464020e-02 3.411833e+02 1.108474e-01
## Gene3 6.089158e+00 1.294014e-02 1.733477e+02 3.618921e-17 8.583661e+03
## Gene4 2.375904e-34 3.879437e+00 4.183261e+00 8.752412e-16 1.891090e-11
## Gene5 2.065816e-02 1.417084e+01 3.082099e-16 1.773306e-06 7.721918e-27
An additional (big) advantage of outputting an SCESet
is that we get immediate access to all of the scater
functions. For example we can make a PCA plot:
plotPCA(sim)
(NOTE: Your values and plots may look different as the simulation is random and produces different results each time it is run.)
For more details of the SCESet
and what you can do with scater
refer to the scater
documentation and vignette.
The splatSimulate
function outputs the following additional information about the simulation:
pData
)
Cell
- Unique cell identifier.Group
- The group or path the cell belongs to.ExpLibSize
- The expected library size for that cell.Step
(paths only) - How far along the path each cell is.fData
)
Gene
- Unique gene identifier.BaseGeneMean
- The base expression level for that gene.OutlierFactor
- Expression outlier factor for that gene (1 is not an outlier).GeneMean
- Expression level after applying outlier factors.DEFac[Group]
- The differential expression factor for each gene in a particular group (1 is not differentially expressed).GeneMean[Group]
- Expression level of a gene in a particular group after applying differential expression factors.assayData
)
BaseCellMeans
- The expression of genes in each cell adjusted for expected library size.BCV
- The Biological Coefficient of Variation for each gene in each cell.CellMeans
- The expression level of genes in each cell adjusted for BCV.TrueCounts
- The simulated counts before dropout.Dropout
- Logical matrix showing which counts have been dropped in which cells.Values that have been added by Splatter are named using UpperCamelCase
to separate them from the underscore_naming
used by scater
. For more information on the simulation see ?splatSimulate
.
So far we have only simulated a single population of cells but often we are interested in investigating a mixed population of cells and looking to see what cell types are present or what differences there are between them. Splatter is able to simulate these situations by changing the method
argument Here we are going to simulate two groups, each with 50 cells, by specifying the groupCells
parameter and setting the method
parameter to "groups"
:
(NOTE: We have also set the verbose
argument to FALSE
to stop Splatter printing progress messages.)
sim.groups <- splatSimulate(groupCells = c(100, 100), method = "groups",
verbose = FALSE)
plotPCA(sim.groups, colour_by = "Group")
The other situation that is often of interest is a differentiation process where one cell type is changing into another. Splatter approximates this process by simulating a series of steps between two groups and randomly assigning each cell to a step. We can create this kind of simulation using the "paths"
method.
sim.paths <- splatSimulate(method = "paths", verbose = FALSE)
plotPCA(sim.paths, colour_by = "Step")
Here the colours represent the “step” of each cell or how far along the differentiation path it is. We can see that the cells with dark colours are more similar to the originating cell type and the light coloured cells are closer to the final, differentiated, cell type. By setting additional parameters it is possible to simulate more complex process (for example multiple mature cell types from a single progenitor).
Each of the Splatter simulation methods has it’s own convenience function. To simulate a single population use splatSimulateSingle()
(equivalent to splatSimulate(method = "single")
), to simulate grops use splatSimulateGroups()
(equivalent to splatSimulate(method = "groups")
) or to simulate paths use splatSimulatePaths()
(equivalent to splatSimulate(method = "paths")
).
As well as it’s own simulation method the Splatter package contains implementations of other single-cell RNA-seq simulations that have been published or wrappers around simulations included in other packages. To see all the available simulations run the listSims()
function:
listSims()
## Splatter currently contains 8 simulations
##
## Splat (splat)
## DOI: Github:
## The Splat simulation generates means from a gamma distribution, adjusts them for BCV and generates counts from a gamma-poisson. Dropout can be optionally added.
##
## Splat Single (splatSingle)
## DOI: Github:
## The Splat simulation with a single population.
##
## Splat Groups (splatGroups)
## DOI: Github:
## The Splat simulation with multiple groups. Each group can have it's own differential expression probability and fold change distribution.
##
## Splat Paths (splatPaths)
## DOI: Github:
## The Splat simulation with differentiation paths. Each path can have it's own length, skew and probability. Genes can change in non-linear ways.
##
## Simple (simple)
## DOI: Github:
## A simple simulation with gamma means and negative binomial counts.
##
## Lun (lun)
## DOI: 10.1186/s13059-016-0947-7 Github: MarioniLab/Deconvolution2016
## Gamma distributed means and negative binomial counts. Cells are given a size factor and differential expression can be simulated with fixed fold changes.
##
## Lun 2 (lun2)
## DOI: 10.1101/073973 Github: MarioniLab/PlateEffects2016
## Negative binomial counts where the means and dispersions have been sampled from a real dataset. The core feature of the Lun 2 simulation is the addition of plate effects. Differential expression can be added between two groups of plates and optionally a zero-inflated negative-binomial can be used.
##
## scDD (scDD)
## DOI: 10.1186/s13059-016-1077-y Github: kdkorthauer/scDD
## The scDD simulation samples a given dataset and can simulate differentially expressed and differentially distributed genes between two conditions.
(or more conveniently for the vignette as a table)
knitr::kable(listSims(print = FALSE))
Name | Prefix | DOI | Github | Description |
---|---|---|---|---|
Splat | splat | The Splat simulation generates means from a gamma distribution, adjusts them for BCV and generates counts from a gamma-poisson. Dropout can be optionally added. | ||
Splat Single | splatSingle | The Splat simulation with a single population. | ||
Splat Groups | splatGroups | The Splat simulation with multiple groups. Each group can have it’s own differential expression probability and fold change distribution. | ||
Splat Paths | splatPaths | The Splat simulation with differentiation paths. Each path can have it’s own length, skew and probability. Genes can change in non-linear ways. | ||
Simple | simple | A simple simulation with gamma means and negative binomial counts. | ||
Lun | lun | 10.1186/s13059-016-0947-7 | MarioniLab/Deconvolution2016 | Gamma distributed means and negative binomial counts. Cells are given a size factor and differential expression can be simulated with fixed fold changes. |
Lun 2 | lun2 | 10.1101/073973 | MarioniLab/PlateEffects2016 | Negative binomial counts where the means and dispersions have been sampled from a real dataset. The core feature of the Lun 2 simulation is the addition of plate effects. Differential expression can be added between two groups of plates and optionally a zero-inflated negative-binomial can be used. |
scDD | scDD | 10.1186/s13059-016-1077-y | kdkorthauer/scDD | The scDD simulation samples a given dataset and can simulate differentially expressed and differentially distributed genes between two conditions. |
Each simulation has it’s own prefix which gives the name of the functions associated with that simulation. For example the prefix for the simple simulation is simple
so it would store it’s parameters in a SimpleParams
object that can be created using newSimpleParams()
or estimated from real data using simpleEstimate()
. To simulate data using that simulation you would use simpleSimulate()
. Each simulation returns an SCESet
object with intermediate values similar to that returned by splatSimulate()
. For more detailed information on each simulation see the appropriate help page (eg. ?simpleSimulate
for information on how the simple simulation works or ?lun2Estimate
for details of how the Lun 2 simulation estimates parameters) or refer to the appropriate paper or package.
Splatter is designed to simulate count data but some analysis methods expect other expression values, particularly length-normalised values such as TPM or FPKM. The scater
package has functions for adding these values to an SCESet
object but they require a length for each gene. The addGeneLengths
can be used to simulate these lengths:
sim <- simpleSimulate(verbose = FALSE)
sim <- addGeneLengths(sim)
head(fData(sim))
## Gene GeneMean Length
## Gene1 Gene1 1.603842e-07 1472
## Gene2 Gene2 1.770719e-01 1888
## Gene3 Gene3 6.478439e+00 2036
## Gene4 Gene4 1.478556e+00 1610
## Gene5 Gene5 1.859077e+00 1306
## Gene6 Gene6 2.803226e-01 3151
We can then use scater
to calculate TPM:
tpm(sim) <- calculateTPM(sim, fData(sim)$Length)
tpm(sim)[1:5, 1:5]
## Cell1 Cell2 Cell3 Cell4 Cell5
## Gene1 0.0000 0.00000 0.0000 0.0000 0.0000
## Gene2 0.0000 0.00000 0.0000 0.0000 0.0000
## Gene3 705.6530 77.17703 470.1049 618.5920 641.5033
## Gene4 198.3036 195.19558 0.0000 0.0000 0.0000
## Gene5 488.9262 0.00000 366.4371 241.0899 125.0096
The default method used by addGeneLengths
to simulate lengths is to generate values from a log-normal distribution which are then rounded to give an integer length. The parameters for this distribution are based on human coding genes but can be adjusted if needed (for example for other species). Alternatively lengths can be sampled from a provided vector (see ?addGeneLengths
for details and an example).
One thing you might like to do after simulating data is to compare it to a real dataset, or compare simulations with different parameters or models. Splatter provides a function compareSCESets
that aims to make these comparisons easier. As the name suggests this function takes a list of SCESet
objects, combines the datasets and produces some plots comparing them. Let’s make two small simulations and see how they compare.
sim1 <- splatSimulate(nGenes = 1000, groupCells = 20, verbose = FALSE)
sim2 <- simpleSimulate(nGenes = 1000, nCells = 20, verbose = FALSE)
comparison <- compareSCESets(list(Splat = sim1, Simple = sim2))
names(comparison)
## [1] "FeatureData" "PhenoData" "Plots"
names(comparison$Plots)
## [1] "Means" "Variances" "MeanVar" "LibrarySizes"
## [5] "ZerosGene" "ZerosCell" "MeanZeros"
The returned list has three items. The first two are the combined datasets by gene (FeatureData
) and by cell (PhenoData
) and the third contains some comparison plots (produced using ggplot2
), for example a plot of the distribution of means:
comparison$Plots$Means
These are only a few of the plots you might want to consider but it should be easy to make more using the returned data. For example, we could plot the number of expressed genes against the library size:
library("ggplot2")
ggplot(comparison$PhenoData,
aes(x = total_counts, y = total_features, colour = Dataset)) +
geom_point()
Sometimes instead of visually comparing datasets it may be more interesting to look at the differences between them. We can do this using the diffSCESets
function. Similar to compareSCESets
this function takes a list of SCESet
objects but now we also specify one to be a reference. A series of similar plots are returned but instead of showing the overall distributions they demonstrate differences from the reference.
difference <- diffSCESets(list(Splat = sim1, Simple = sim2), ref = "Simple")
difference$Plots$Means
We also get a series of Quantile-Quantile plot that can be used to compare distributions.
difference$QQPlots$Means
Each of these comparisons makes several plots which can be a lot to look at. To make this easier, or to produce figures for publications, you can make use of the functions makeCompPanel
, makeDiffPanel
and makeOverallPanel
.
These functions combine the plots into a single panel using the cowplot
package. The panels can be quite large and hard to view (for example in RStudio’s plot viewer) so it can be better to output the panels and view them separately. Luckily cowplot
provides a convenient function for saving the images. Here are some suggested parameters for outputting each of the panels:
# This code is just an example and is not run
panel <- makeCompPanel(comparison)
cowplot::save_plot("comp_panel.png", panel, nrow = 4, ncol = 3)
panel <- makeDiffPanel(difference)
cowplot::save_plot("diff_panel.png", panel, nrow = 3, ncol = 5)
panel <- makeOverallPanel(comparison, difference)
cowplot::save_plot("overall_panel.png", panel, ncol = 4, nrow = 7)
If you use Splatter in your work please cite our paper:
citation("splatter")
##
## Zappia L, Phipson B, Oshlack A. Splatter: Simulation Of
## Single-Cell RNA Sequencing Data. bioRxiv. 2017;
## doi:10.1101/133173
##
## A BibTeX entry for LaTeX users is
##
## @Article{,
## author = {Luke Zappia and Belinda Phipson and Alicia Oshlack},
## title = {Splatter: Simulation Of Single-Cell RNA Sequencing Data},
## journal = {bioRxiv},
## year = {2017},
## url = {http://dx.doi.org/10.1101/133173},
## doi = {10.1101/133173},
## }
sessionInfo()
## R version 3.4.0 (2017-04-21)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.5-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.5-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] parallel stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] splatter_1.0.3 scater_1.4.0 ggplot2_2.2.1
## [4] Biobase_2.36.2 BiocGenerics_0.22.0 BiocStyle_2.4.0
##
## loaded via a namespace (and not attached):
## [1] viridis_0.4.0 edgeR_3.18.1
## [3] splines_3.4.0 viridisLite_0.2.0
## [5] shiny_1.0.3 assertthat_0.2.0
## [7] highr_0.6 sp_1.2-4
## [9] stats4_3.4.0 GenomeInfoDbData_0.99.0
## [11] vipor_0.4.5 yaml_2.1.14
## [13] RSQLite_1.1-2 backports_1.1.0
## [15] lattice_0.20-35 limma_3.32.2
## [17] digest_0.6.12 GenomicRanges_1.28.3
## [19] XVector_0.16.0 checkmate_1.8.2
## [21] colorspace_1.3-2 cowplot_0.7.0
## [23] htmltools_0.3.6 httpuv_1.3.3
## [25] Matrix_1.2-10 plyr_1.8.4
## [27] XML_3.98-1.7 biomaRt_2.32.0
## [29] zlibbioc_1.22.0 xtable_1.8-2
## [31] scales_0.4.1 BiocParallel_1.10.1
## [33] tibble_1.3.1 IRanges_2.10.2
## [35] SummarizedExperiment_1.6.1 lazyeval_0.2.0
## [37] survival_2.41-3 magrittr_1.5
## [39] mime_0.5 memoise_1.1.0
## [41] evaluate_0.10 MASS_7.3-47
## [43] beeswarm_0.2.3 shinydashboard_0.6.0
## [45] tools_3.4.0 fitdistrplus_1.0-9
## [47] data.table_1.10.4 matrixStats_0.52.2
## [49] stringr_1.2.0 S4Vectors_0.14.2
## [51] munsell_0.4.3 locfit_1.5-9.1
## [53] DelayedArray_0.2.4 AnnotationDbi_1.38.0
## [55] akima_0.6-2 compiler_3.4.0
## [57] GenomeInfoDb_1.12.0 rlang_0.1.1
## [59] rhdf5_2.20.0 grid_3.4.0
## [61] RCurl_1.95-4.8 tximport_1.4.0
## [63] rjson_0.2.15 labeling_0.3
## [65] bitops_1.0-6 rmarkdown_1.5
## [67] gtable_0.2.0 DBI_0.6-1
## [69] reshape2_1.4.2 R6_2.2.1
## [71] gridExtra_2.2.1 knitr_1.16
## [73] dplyr_0.5.0 rprojroot_1.2
## [75] stringi_1.1.5 ggbeeswarm_0.5.3
## [77] Rcpp_0.12.11