Getting started with SimBu

Alexander Dietrich

Installation

To install the developmental version of the package, run:

install.packages("devtools")
devtools::install_github("omnideconv/SimBu")

To install from Bioconductor:

if (!require("BiocManager", quietly = TRUE)) {
  install.packages("BiocManager")
}

BiocManager::install("SimBu")
library(SimBu)

Introduction

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.

Getting started

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)
  )
)

Creating a dataset

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:

Simulate pseudo bulk datasets

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:

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.2

Results

The simulation object contains three named entries:

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 522 552 524 480 506 524 550 491 479 539
#> gene_2 541 467 466 465 495 487 508 503 486 509
#> gene_3 502 516 522 503 486 500 461 466 507 507
#> gene_4 514 500 513 560 492 532 555 558 527 523
#> gene_5 473 480 504 514 480 467 506 484 523 457
#> gene_6 470 485 481 508 493 507 466 466 471 520
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 1009.5311 1003.9449 1008.4060  964.5361  980.8886  990.7845 1019.4537
#> gene_2 1009.2044 1084.4546  906.1917  910.6995  954.7952  971.8023  962.7392
#> gene_3  967.6697  984.7055 1049.1025 1041.4782 1031.3299  956.9811 1003.5883
#> gene_4  948.0173  948.1869  982.1590  911.9107  913.9551 1004.9987  890.4766
#> gene_5 1019.3639 1018.6326 1072.9398  999.9606  996.8291 1035.4879  976.2823
#> gene_6 1019.6525 1002.2042  996.3735 1037.6534  943.5846  988.6765  973.8709
#>                                     
#> gene_1  912.0217  952.4540  986.4352
#> gene_2  921.9204  987.9012 1011.4397
#> gene_3 1090.0348 1025.7881 1041.5170
#> gene_4  979.6003  887.4653  906.2881
#> gene_5  939.3864 1013.9636 1008.4729
#> gene_6 1056.9306 1030.5453  982.9749

If only a single matrix was given to the dataset initially, only one assay is filled.

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 per session.
#> Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
#> generated.

More features

Simulate using a whitelist (and blacklist) of cell-types

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.6.0 Patched (2026-04-24 r89963)
#> Platform: aarch64-apple-darwin23
#> Running under: macOS Tahoe 26.3.1
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.6/Resources/lib/libRblas.0.dylib 
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.6/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.1
#> 
#> locale:
#> [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#> 
#> time zone: America/New_York
#> tzcode source: internal
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] SimBu_1.14.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] sass_0.4.10                 generics_0.1.4             
#>  [3] tidyr_1.3.2                 SparseArray_1.12.0         
#>  [5] lattice_0.22-9              digest_0.6.39              
#>  [7] magrittr_2.0.5              RColorBrewer_1.1-3         
#>  [9] evaluate_1.0.5              sparseMatrixStats_1.24.0   
#> [11] grid_4.6.0                  fastmap_1.2.0              
#> [13] jsonlite_2.0.0              Matrix_1.7-5               
#> [15] proxyC_0.5.2                purrr_1.2.2                
#> [17] scales_1.4.0                codetools_0.2-20           
#> [19] jquerylib_0.1.4             abind_1.4-8                
#> [21] cli_3.6.6                   rlang_1.2.0                
#> [23] XVector_0.52.0              Biobase_2.72.0             
#> [25] withr_3.0.2                 cachem_1.1.0               
#> [27] DelayedArray_0.38.0         yaml_2.3.12                
#> [29] otel_0.2.0                  S4Arrays_1.12.0            
#> [31] tools_4.6.0                 parallel_4.6.0             
#> [33] BiocParallel_1.46.0         dplyr_1.2.1                
#> [35] ggplot2_4.0.3               SummarizedExperiment_1.42.0
#> [37] BiocGenerics_0.58.0         vctrs_0.7.3                
#> [39] R6_2.6.1                    matrixStats_1.5.0          
#> [41] stats4_4.6.0                lifecycle_1.0.5            
#> [43] Seqinfo_1.2.0               S4Vectors_0.50.0           
#> [45] IRanges_2.46.0              pkgconfig_2.0.3            
#> [47] gtable_0.3.6                bslib_0.10.0               
#> [49] pillar_1.11.1               data.table_1.18.2.1        
#> [51] glue_1.8.1                  Rcpp_1.1.1-1.1             
#> [53] xfun_0.57                   tibble_3.3.1               
#> [55] GenomicRanges_1.64.0        tidyselect_1.2.1           
#> [57] dichromat_2.0-0.1           MatrixGenerics_1.24.0      
#> [59] knitr_1.51                  farver_2.1.2               
#> [61] htmltools_0.5.9             labeling_0.4.3             
#> [63] rmarkdown_2.31              compiler_4.6.0             
#> [65] S7_0.2.2

References

Fischer, David S., Leander Dony, Martin König, Abdul Moeed, Luke Zappia, Sophie Tritschler, Olle Holmberg, Hananeh Aliee, and Fabian J. Theis. 2020. “Sfaira Accelerates Data and Model Reuse in Single Cell Genomics.” bioRxiv. https://doi.org/10.1101/2020.12.16.419036.