To install the developmental version of the package, run:
To install from Bioconductor:
As complex tissues are typically composed of various cell types,
deconvolution tools have been developed to computationally infer their
cellular composition from bulk RNA sequencing (RNA-seq) data. To
comprehensively assess deconvolution performance, gold-standard datasets
are indispensable. Gold-standard, experimental techniques like flow
cytometry or immunohistochemistry are resource-intensive and cannot be
systematically applied to the numerous cell types and tissues profiled
with high-throughput transcriptomics. The simulation of ‘pseudo-bulk’
data, generated by aggregating single-cell RNA-seq (scRNA-seq)
expression profiles in pre-defined proportions, offers a scalable and
cost-effective alternative. This makes it feasible to create in silico
gold standards that allow fine-grained control of cell-type fractions
not conceivable in an experimental setup. However, at present, no
simulation software for generating pseudo-bulk RNA-seq data
exists.
SimBu was developed to simulate pseudo-bulk samples based on various
simulation scenarios, designed to test specific features of
deconvolution methods. A unique feature of SimBu is the modelling of
cell-type-specific mRNA bias using experimentally-derived or data-driven
scaling factors. Here, we show that SimBu can generate realistic
pseudo-bulk data, recapitulating the biological and statistical features
of real RNA-seq data. Finally, we illustrate the impact of mRNA bias on
the evaluation of deconvolution tools and provide recommendations for
the selection of suitable methods for estimating mRNA content.
This chapter covers all you need to know to quickly simulate some
pseudo-bulk samples!
This package can simulate samples from local or public data. This
vignette will work with artificially generated data as it serves as an
overview for the features implemented in SimBu. For the public data
integration using sfaira (Fischer et al. 2020), please refer to
the “Public Data Integration”
vignette.
We will create some toy data to use for our simulations; two matrices with 300 cells each and 1000 genes/features. One represents raw count data, while the other matrix represents scaled TPM-like data. We will assign these cells to some immune cell types.
counts <- Matrix::Matrix(matrix(stats::rpois(3e5, 5), ncol = 300), sparse = TRUE)
tpm <- Matrix::Matrix(matrix(stats::rpois(3e5, 5), ncol = 300), sparse = TRUE)
tpm <- Matrix::t(1e6 * Matrix::t(tpm) / Matrix::colSums(tpm))
colnames(counts) <- paste0("cell_", rep(1:300))
colnames(tpm) <- paste0("cell_", rep(1:300))
rownames(counts) <- paste0("gene_", rep(1:1000))
rownames(tpm) <- paste0("gene_", rep(1:1000))
annotation <- data.frame(
"ID" = paste0("cell_", rep(1:300)),
"cell_type" = c(
rep("T cells CD4", 50),
rep("T cells CD8", 50),
rep("Macrophages", 100),
rep("NK cells", 10),
rep("B cells", 70),
rep("Monocytes", 20)
)
)SimBu uses the SummarizedExperiment
class as storage for count data as well as annotation data. Currently it
is possible to store two matrices at the same time: raw counts and
TPM-like data (this can also be some other scaled count matrix, such as
RPKM, but we recommend to use TPMs). These two matrices have to have the
same dimensions and have to contain the same genes and cells. Providing
the raw count data is mandatory!
SimBu scales the matrix that is added via the tpm_matrix
slot by default to 1e6 per cell, if you do not want this, you can switch
it off by setting the scale_tpm parameter to
FALSE. Additionally, the cell type annotation of the cells
has to be given in a dataframe, which has to include the two columns
ID and cell_type. If additional columns from
this annotation should be transferred to the dataset, simply give the
names of them in the additional_cols parameter.
To generate a dataset that can be used in SimBu, you can use the
dataset() method; other methods exist as well, which are
covered in the “Inputs &
Outputs” vignette.
ds <- SimBu::dataset(
annotation = annotation,
count_matrix = counts,
tpm_matrix = tpm,
name = "test_dataset"
)
#> Filtering genes...
#> Created dataset.SimBu offers basic filtering options for your dataset, which you can
apply during dataset generation:
filter_genes: if TRUE, genes which have expression
values of 0 in all cells will be removed.
variance_cutoff: remove all genes with a expression
variance below the chosen cutoff.
type_abundance_cutoff: remove all cells, which belong to
a cell type that appears less the the given amount.
We are now ready to simulate the first pseudo bulk samples with the created dataset:
simulation <- SimBu::simulate_bulk(
data = ds,
scenario = "random",
scaling_factor = "NONE",
ncells = 100,
nsamples = 10,
BPPARAM = BiocParallel::MulticoreParam(workers = 4), # this will use 4 threads to run the simulation
run_parallel = TRUE
) # multi-threading to TRUE
#> Using parallel generation of simulations.
#> Finished simulation.ncells sets the number of cells in each sample, while
nsamples sets the total amount of simulated samples.
If you want to simulate a specific sequencing depth in your simulations,
you can use the total_read_counts parameter to do so. Note
that this parameter is only applied on the counts matrix (if supplied),
as TPMs will be scaled to 1e6 by default.
SimBu can add mRNA bias by using different scaling factors to the
simulations using the scaling_factor parameter. A detailed
explanation can be found in the “Scaling factor”
vignette.
Currently there are 6 scenarios implemented in the
package:
even: this creates samples, where all existing
cell-types in the dataset appear in the same proportions. So using a
dataset with 3 cell-types, this will simulate samples, where all
cell-type fractions are 1/3. In order to still have a slight variation
between cell type fractions, you can increase the
balance_uniform_mirror_scenario parameter (default to
0.01). Setting it to 0 will generate simulations with exactly the same
cell type fractions.
random: this scenario will create random cell type
fractions using all present types for each sample. The random sampling
is based on the uniform distribution.
mirror_db: this scenario will mirror the exact fractions
of cell types which are present in the provided dataset. If it consists
of 20% T cells, 30% B cells and 50% NK cells, all simulated samples will
mirror these fractions. Similar to the uniform scenario, you can add a
small variation to these fractions with the
balance_uniform_mirror_scenario parameter.
weighted: here you need to set two additional parameters
for the simulate_bulk() function:
weighted_cell_type sets the cell-type you want to be
over-representing and weighted_amount sets the fraction of
this cell-type. You could for example use B-cell and
0.5 to create samples, where 50% are B-cells and the rest
is filled randomly with other cell-types.
pure: this creates simulations of only one single
cell-type. You have to provide the name of this cell-type with the
pure_cell_type parameter.
custom: here you are able to create your own set of
cell-type fractions. When using this scenario, you additionally need to
provide a dataframe in the custom_scenario_data parameter,
where each row represents one sample (therefore the number of rows need
to match the nsamples parameter). Each column has to
represent one cell-type, which also occurs in the dataset and describes
the fraction of this cell-type in a sample. The fractions per sample
need to sum up to 1. An example can be seen here:
pure_scenario_dataframe <- data.frame(
"B cells" = c(0.2, 0.1, 0.5, 0.3),
"T cells" = c(0.3, 0.8, 0.2, 0.5),
"NK cells" = c(0.5, 0.1, 0.3, 0.2),
row.names = c("sample1", "sample2", "sample3", "sample4")
)
pure_scenario_dataframe
#> B.cells T.cells NK.cells
#> sample1 0.2 0.3 0.5
#> sample2 0.1 0.8 0.1
#> sample3 0.5 0.2 0.3
#> sample4 0.3 0.5 0.2The simulation object contains three named
entries:
bulk: a SummarizedExperiment object with the
pseudo-bulk dataset(s) stored in the assays. They can be
accessed like this:utils::head(SummarizedExperiment::assays(simulation$bulk)[["bulk_counts"]])
#> 6 x 10 sparse Matrix of class "dgCMatrix"
#> [[ suppressing 10 column names 'random_sample1', 'random_sample2', 'random_sample3' ... ]]
#>
#> gene_1 518 477 472 491 459 494 470 473 506 474
#> gene_2 481 509 520 487 510 524 509 469 453 545
#> gene_3 499 488 495 465 499 493 509 464 524 496
#> gene_4 493 476 472 470 535 495 527 489 496 487
#> gene_5 490 478 506 502 517 526 490 492 481 497
#> gene_6 522 522 519 527 528 496 539 545 532 464
utils::head(SummarizedExperiment::assays(simulation$bulk)[["bulk_tpm"]])
#> 6 x 10 sparse Matrix of class "dgCMatrix"
#> [[ suppressing 10 column names 'random_sample1', 'random_sample2', 'random_sample3' ... ]]
#>
#> gene_1 1022.9054 1087.6912 1086.4826 1083.5372 1034.842 1067.3905 1077.5534
#> gene_2 979.0954 991.9996 1086.1783 1020.2198 1064.737 952.9807 1020.1111
#> gene_3 1124.0538 1034.5137 1034.1799 943.3476 1068.523 1000.3675 1020.8917
#> gene_4 962.2917 968.1803 1009.7378 974.2652 1010.294 959.1355 998.0130
#> gene_5 963.2811 929.7096 881.7054 932.1803 941.984 909.6101 1022.7713
#> gene_6 898.9772 1032.6962 1037.1859 968.9235 1011.939 929.3633 933.5433
#>
#> gene_1 1003.0721 1116.7824 1051.7913
#> gene_2 1008.0346 1021.4528 936.5431
#> gene_3 933.3899 955.3950 1086.7755
#> gene_4 1009.8790 881.7924 1017.3805
#> gene_5 990.5483 988.6337 931.5012
#> gene_6 1021.9958 962.1853 921.4133If only a single matrix was given to the dataset initially, only one assay is filled.
cell_fractions: a table where rows represent the
simulated samples and columns represent the different simulated
cell-types. The entries in the table store the specific cell-type
fraction per sample.
scaling_vector: a named list, with the used scaling
value for each cell from the single cell dataset.
It is also possible to merge simulations:
simulation2 <- SimBu::simulate_bulk(
data = ds,
scenario = "even",
scaling_factor = "NONE",
ncells = 1000,
nsamples = 10,
BPPARAM = BiocParallel::MulticoreParam(workers = 4),
run_parallel = TRUE
)
#> Using parallel generation of simulations.
#> Finished simulation.
merged_simulations <- SimBu::merge_simulations(list(simulation, simulation2))Finally here is a barplot of the resulting simulation:
SimBu::plot_simulation(simulation = merged_simulations)
#> Warning: `aes_string()` was deprecated in ggplot2 3.0.0.
#> ℹ Please use tidy evaluation idioms with `aes()`.
#> ℹ See also `vignette("ggplot2-in-packages")` for more information.
#> ℹ The deprecated feature was likely used in the SimBu package.
#> Please report the issue at <https://github.com/omnideconv/SimBu/issues>.
#> This warning is displayed once every 8 hours.
#> Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
#> generated.Sometimes, you are only interested in specific cell-types (for
example T cells), but the dataset you are using has too many other
cell-types; you can handle this issue during simulation using the
whitelist parameter:
simulation <- SimBu::simulate_bulk(
data = ds,
scenario = "random",
scaling_factor = "NONE",
ncells = 1000,
nsamples = 20,
BPPARAM = BiocParallel::MulticoreParam(workers = 4),
run_parallel = TRUE,
whitelist = c("T cells CD4", "T cells CD8")
)
#> Using parallel generation of simulations.
#> Finished simulation.
SimBu::plot_simulation(simulation = simulation)In the same way, you can also provide a blacklist
parameter, where you name the cell-types you don’t want
to be included in your simulation.
utils::sessionInfo()
#> R version 4.5.1 (2025-06-13)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> 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
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] SimBu_1.13.0 rmarkdown_2.30
#>
#> loaded via a namespace (and not attached):
#> [1] sass_0.4.10 generics_0.1.4
#> [3] tidyr_1.3.1 SparseArray_1.9.2
#> [5] lattice_0.22-7 digest_0.6.37
#> [7] magrittr_2.0.4 RColorBrewer_1.1-3
#> [9] evaluate_1.0.5 sparseMatrixStats_1.21.0
#> [11] grid_4.5.1 fastmap_1.2.0
#> [13] jsonlite_2.0.0 Matrix_1.7-4
#> [15] proxyC_0.5.2 purrr_1.1.0
#> [17] scales_1.4.0 codetools_0.2-20
#> [19] jquerylib_0.1.4 abind_1.4-8
#> [21] cli_3.6.5 crayon_1.5.3
#> [23] rlang_1.1.6 XVector_0.49.3
#> [25] Biobase_2.69.1 withr_3.0.2
#> [27] cachem_1.1.0 DelayedArray_0.37.0
#> [29] yaml_2.3.10 S4Arrays_1.11.0
#> [31] tools_4.5.1 parallel_4.5.1
#> [33] BiocParallel_1.43.4 dplyr_1.1.4
#> [35] ggplot2_4.0.0 SummarizedExperiment_1.39.2
#> [37] BiocGenerics_0.55.4 buildtools_1.0.0
#> [39] vctrs_0.6.5 R6_2.6.1
#> [41] matrixStats_1.5.0 stats4_4.5.1
#> [43] lifecycle_1.0.4 Seqinfo_0.99.4
#> [45] S4Vectors_0.49.0 IRanges_2.45.0
#> [47] pkgconfig_2.0.3 gtable_0.3.6
#> [49] bslib_0.9.0 pillar_1.11.1
#> [51] data.table_1.17.8 glue_1.8.0
#> [53] Rcpp_1.1.0 tidyselect_1.2.1
#> [55] xfun_0.53 tibble_3.3.0
#> [57] GenomicRanges_1.63.0 sys_3.4.3
#> [59] MatrixGenerics_1.21.0 knitr_1.50
#> [61] farver_2.1.2 htmltools_0.5.8.1
#> [63] labeling_0.4.3 maketools_1.3.2
#> [65] compiler_4.5.1 S7_0.2.0