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.
#> Warning in BiocParallel::MulticoreParam(workers = 4): MulticoreParam() not
#> supported on Windows, use SnowParam()
#> 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 524 516 498 496 527 527 572 507 503
#> gene_2 498 502 492 491 499 519 506 489 502 478
#> gene_3 508 544 510 491 521 545 507 573 497 515
#> gene_4 497 468 451 453 392 505 489 508 472 472
#> gene_5 543 491 533 502 500 526 484 531 507 531
#> gene_6 514 515 531 503 517 502 524 553 544 530
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 1035.7999  981.5290  993.5668 1013.3958  988.6360  920.0082 1009.7324
#> gene_2  946.4048 1015.9853  971.1080  992.2168 1002.9359  945.5341 1002.7363
#> gene_3  967.2391  932.2821  965.1017  894.4934  975.2662  959.2517  995.2966
#> gene_4 1070.2462 1019.6009 1080.8562 1122.8259  972.9105 1010.5321 1086.0225
#> gene_5 1051.8988  961.3730 1013.3511 1020.8820 1116.8221 1085.8903  983.0459
#> gene_6  998.9803  974.9426 1052.0666 1008.5478  980.4540 1045.2421 1121.9053
#>                                    
#> gene_1 1016.4763  999.0915 994.4497
#> gene_2 1028.9139 1058.7327 919.8217
#> gene_3  945.0728 1009.0733 945.2907
#> gene_4 1018.2411 1013.4746 988.5993
#> gene_5 1088.2163  976.3887 930.6041
#> gene_6  975.5327  989.9842 940.4266

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.
#> Warning in BiocParallel::MulticoreParam(workers = 4): MulticoreParam() not
#> supported on Windows, use SnowParam()
#> 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.
#> Warning in BiocParallel::MulticoreParam(workers = 4): MulticoreParam() not
#> supported on Windows, use SnowParam()
#> 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.4.1 (2024-06-14 ucrt)
#> Platform: x86_64-w64-mingw32/x64
#> Running under: Windows Server 2022 x64 (build 20348)
#> 
#> Matrix products: default
#> 
#> 
#> locale:
#> [1] LC_COLLATE=C                          
#> [2] LC_CTYPE=English_United States.utf8   
#> [3] LC_MONETARY=English_United States.utf8
#> [4] LC_NUMERIC=C                          
#> [5] LC_TIME=English_United States.utf8    
#> 
#> 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.7.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] SummarizedExperiment_1.35.4 gtable_0.3.5               
#>  [3] xfun_0.48                   bslib_0.8.0                
#>  [5] ggplot2_3.5.1               Biobase_2.65.1             
#>  [7] lattice_0.22-6              vctrs_0.6.5                
#>  [9] tools_4.4.1                 generics_0.1.3             
#> [11] stats4_4.4.1                parallel_4.4.1             
#> [13] tibble_3.2.1                fansi_1.0.6                
#> [15] highr_0.11                  pkgconfig_2.0.3            
#> [17] Matrix_1.7-1                data.table_1.16.2          
#> [19] RColorBrewer_1.1-3          S4Vectors_0.43.2           
#> [21] sparseMatrixStats_1.17.2    lifecycle_1.0.4            
#> [23] GenomeInfoDbData_1.2.13     compiler_4.4.1             
#> [25] farver_2.1.2                munsell_0.5.1              
#> [27] codetools_0.2-20            GenomeInfoDb_1.41.2        
#> [29] htmltools_0.5.8.1           sass_0.4.9                 
#> [31] yaml_2.3.10                 pillar_1.9.0               
#> [33] crayon_1.5.3                jquerylib_0.1.4            
#> [35] tidyr_1.3.1                 BiocParallel_1.39.0        
#> [37] DelayedArray_0.31.14        cachem_1.1.0               
#> [39] abind_1.4-8                 tidyselect_1.2.1           
#> [41] digest_0.6.37               dplyr_1.1.4                
#> [43] purrr_1.0.2                 labeling_0.4.3             
#> [45] fastmap_1.2.0               grid_4.4.1                 
#> [47] colorspace_2.1-1            cli_3.6.3                  
#> [49] SparseArray_1.5.45          magrittr_2.0.3             
#> [51] S4Arrays_1.5.11             utf8_1.2.4                 
#> [53] withr_3.0.1                 UCSC.utils_1.1.0           
#> [55] scales_1.3.0                rmarkdown_2.28             
#> [57] XVector_0.45.0              httr_1.4.7                 
#> [59] matrixStats_1.4.1           proxyC_0.4.1               
#> [61] evaluate_1.0.1              knitr_1.48                 
#> [63] GenomicRanges_1.57.2        IRanges_2.39.2             
#> [65] rlang_1.1.4                 Rcpp_1.0.13                
#> [67] glue_1.8.0                  BiocGenerics_0.51.3        
#> [69] jsonlite_1.8.9              R6_2.5.1                   
#> [71] MatrixGenerics_1.17.0       zlibbioc_1.51.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.