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 489 565 513 518 470 486 514 508 484 529
#> gene_2 506 534 515 528 474 464 550 515 533 475
#> gene_3 510 540 514 543 512 520 547 525 531 485
#> gene_4 494 550 483 530 483 532 488 546 476 496
#> gene_5 487 490 540 502 529 538 482 496 497 546
#> gene_6 519 501 526 522 526 544 514 564 585 529
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  925.1227 1015.9095 1049.7771  912.9366 1032.6119 1044.1196 1059.9838
#> gene_2 1013.2067  908.8664  944.5299  931.5420  999.0993 1009.6636  966.5923
#> gene_3 1046.2549 1007.9476 1014.2529  955.0390  982.3698  928.7755  967.3107
#> gene_4 1012.1362 1184.4737 1038.2242 1028.1446 1059.5192 1129.0946 1053.0866
#> gene_5 1022.5235 1042.6725  906.2905  943.6566 1009.2484  983.1068  993.6903
#> gene_6  973.7694  930.7787  983.4536  967.6135 1074.7520 1040.9711  971.6706
#>                                     
#> gene_1  987.4353 1005.5591 1070.3368
#> gene_2  952.0123  946.6407  982.6556
#> gene_3  988.0813  967.2723  986.7432
#> gene_4 1068.3893 1073.8926  945.6184
#> gene_5  995.9577  876.4740  969.1954
#> gene_6 1070.4652 1014.4083 1023.6900

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 every 8 hours.
#> 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.5.1 Patched (2025-09-10 r88807)
#> Platform: x86_64-apple-darwin20
#> Running under: macOS Monterey 12.7.6
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRblas.0.dylib 
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.5-x86_64/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.12.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] sass_0.4.10                 generics_0.1.4             
#>  [3] tidyr_1.3.1                 SparseArray_1.10.0         
#>  [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.22.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.50.0             
#> [25] Biobase_2.70.0              withr_3.0.2                
#> [27] cachem_1.1.0                DelayedArray_0.36.0        
#> [29] yaml_2.3.10                 S4Arrays_1.10.0            
#> [31] tools_4.5.1                 parallel_4.5.1             
#> [33] BiocParallel_1.44.0         dplyr_1.1.4                
#> [35] ggplot2_4.0.0               SummarizedExperiment_1.40.0
#> [37] BiocGenerics_0.56.0         vctrs_0.6.5                
#> [39] R6_2.6.1                    matrixStats_1.5.0          
#> [41] stats4_4.5.1                lifecycle_1.0.4            
#> [43] Seqinfo_1.0.0               S4Vectors_0.48.0           
#> [45] IRanges_2.44.0              pkgconfig_2.0.3            
#> [47] gtable_0.3.6                bslib_0.9.0                
#> [49] pillar_1.11.1               data.table_1.17.8          
#> [51] glue_1.8.0                  Rcpp_1.1.0                 
#> [53] xfun_0.53                   tibble_3.3.0               
#> [55] GenomicRanges_1.62.0        tidyselect_1.2.1           
#> [57] dichromat_2.0-0.1           MatrixGenerics_1.22.0      
#> [59] knitr_1.50                  farver_2.1.2               
#> [61] htmltools_0.5.8.1           labeling_0.4.3             
#> [63] rmarkdown_2.30              compiler_4.5.1             
#> [65] S7_0.2.0

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.