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 496 487 490 497 527 498 454 528 480 547
#> gene_2 470 472 454 447 503 446 447 470 429 496
#> gene_3 476 460 500 451 458 475 460 432 468 492
#> gene_4 489 487 488 470 515 505 498 447 462 475
#> gene_5 485 546 481 511 484 479 514 493 468 504
#> gene_6 516 538 510 537 539 563 533 542 541 515
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  965.7487  956.6520  936.1020  987.3309 1014.2117  961.9135  993.1397
#> gene_2 1002.6537  931.2930  943.6789  948.4675  964.7358  994.7519  958.4227
#> gene_3 1050.5157 1080.9363 1180.0792  941.5597  997.7802 1031.9527 1012.2099
#> gene_4 1025.8436 1074.1624  968.6490 1010.1402 1017.1061 1009.8280  989.4136
#> gene_5  975.2759 1025.2550  911.5721  949.7127  965.2734  953.9865  925.8156
#> gene_6  914.0723  978.3376  997.4245  968.2488 1008.7732  992.1015  972.3097
#>                                     
#> gene_1  918.6621 1015.4525  940.6065
#> gene_2  961.5913  984.9547  872.7640
#> gene_3 1051.3892 1085.5960  945.6254
#> gene_4 1111.8101 1025.3861  904.4521
#> gene_5  924.1633  899.8734  971.3259
#> gene_6  911.2217  974.2427 1005.8994

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)

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 Under development (unstable) (2025-01-20 r87609)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              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: America/New_York
#> tzcode source: system (glibc)
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] SimBu_1.9.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] sass_0.4.9                  generics_0.1.3             
#>  [3] tidyr_1.3.1                 SparseArray_1.7.4          
#>  [5] lattice_0.22-6              digest_0.6.37              
#>  [7] magrittr_2.0.3              RColorBrewer_1.1-3         
#>  [9] evaluate_1.0.3              sparseMatrixStats_1.19.0   
#> [11] grid_4.5.0                  fastmap_1.2.0              
#> [13] jsonlite_1.8.9              Matrix_1.7-2               
#> [15] GenomeInfoDb_1.43.3         proxyC_0.4.1               
#> [17] httr_1.4.7                  purrr_1.0.2                
#> [19] scales_1.3.0                UCSC.utils_1.3.1           
#> [21] codetools_0.2-20            jquerylib_0.1.4            
#> [23] abind_1.4-8                 cli_3.6.3                  
#> [25] rlang_1.1.5                 crayon_1.5.3               
#> [27] XVector_0.47.2              Biobase_2.67.0             
#> [29] munsell_0.5.1               withr_3.0.2                
#> [31] cachem_1.1.0                DelayedArray_0.33.4        
#> [33] yaml_2.3.10                 S4Arrays_1.7.1             
#> [35] tools_4.5.0                 parallel_4.5.0             
#> [37] BiocParallel_1.41.0         dplyr_1.1.4                
#> [39] colorspace_2.1-1            ggplot2_3.5.1              
#> [41] GenomeInfoDbData_1.2.13     SummarizedExperiment_1.37.0
#> [43] BiocGenerics_0.53.5         vctrs_0.6.5                
#> [45] R6_2.5.1                    matrixStats_1.5.0          
#> [47] stats4_4.5.0                lifecycle_1.0.4            
#> [49] S4Vectors_0.45.2            IRanges_2.41.2             
#> [51] pkgconfig_2.0.3             gtable_0.3.6               
#> [53] bslib_0.8.0                 pillar_1.10.1              
#> [55] data.table_1.16.4           glue_1.8.0                 
#> [57] Rcpp_1.0.14                 tidyselect_1.2.1           
#> [59] xfun_0.50                   tibble_3.2.1               
#> [61] GenomicRanges_1.59.1        MatrixGenerics_1.19.1      
#> [63] knitr_1.49                  farver_2.1.2               
#> [65] htmltools_0.5.8.1           labeling_0.4.3             
#> [67] rmarkdown_2.29              compiler_4.5.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.