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 447 517 489 500 466 523 426 466 509 495
#> gene_2 531 506 479 532 527 584 498 534 546 506
#> gene_3 472 445 491 556 436 471 455 474 447 528
#> gene_4 476 476 481 463 503 506 491 476 508 493
#> gene_5 514 412 492 512 484 519 455 508 543 514
#> gene_6 520 569 492 551 546 511 530 544 542 552
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  866.1874 981.0445  976.9384 1025.2350  970.3928 957.6776  907.7802
#> gene_2  966.3444 980.0366  932.7332 1111.4390  968.4383 992.0178 1075.9928
#> gene_3  923.6296 906.5138  990.6548  977.4324  901.9055 978.0750  919.1477
#> gene_4  944.5888 870.1909  947.4791  837.6072  910.6976 958.9851  918.8770
#> gene_5 1099.7801 977.2179  973.4431 1056.1359 1006.4795 976.4724 1063.7510
#> gene_6  982.0899 959.1222 1054.4973  994.7651 1003.7294 885.7927  986.0464
#>                                     
#> gene_1  918.6967  892.4890  928.4991
#> gene_2 1015.8794  985.6254  907.4248
#> gene_3  947.1294 1054.5231 1062.5693
#> gene_4  909.7051  831.4651  832.5745
#> gene_5 1046.5421  966.8061 1023.6006
#> gene_6  957.2819  982.7361 1113.7042

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 version 4.4.0 alpha (2024-03-27 r86216)
#> Platform: aarch64-apple-darwin20
#> Running under: macOS Ventura 13.6.5
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib 
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.0
#> 
#> locale:
#> [1] en_US.UTF-8/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.6.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] SummarizedExperiment_1.34.0 gtable_0.3.4               
#>  [3] xfun_0.43                   bslib_0.6.2                
#>  [5] ggplot2_3.5.0               Biobase_2.64.0             
#>  [7] lattice_0.22-6              vctrs_0.6.5                
#>  [9] tools_4.4.0                 generics_0.1.3             
#> [11] stats4_4.4.0                parallel_4.4.0             
#> [13] tibble_3.2.1                fansi_1.0.6                
#> [15] highr_0.10                  pkgconfig_2.0.3            
#> [17] Matrix_1.7-0                data.table_1.15.4          
#> [19] RColorBrewer_1.1-3          S4Vectors_0.42.0           
#> [21] sparseMatrixStats_1.16.0    RcppParallel_5.1.7         
#> [23] lifecycle_1.0.4             GenomeInfoDbData_1.2.12    
#> [25] farver_2.1.1                compiler_4.4.0             
#> [27] munsell_0.5.0               codetools_0.2-19           
#> [29] GenomeInfoDb_1.40.0         htmltools_0.5.8            
#> [31] sass_0.4.9                  yaml_2.3.8                 
#> [33] pillar_1.9.0                crayon_1.5.2               
#> [35] jquerylib_0.1.4             tidyr_1.3.1                
#> [37] BiocParallel_1.38.0         DelayedArray_0.30.0        
#> [39] cachem_1.0.8                abind_1.4-5                
#> [41] tidyselect_1.2.1            digest_0.6.35              
#> [43] dplyr_1.1.4                 purrr_1.0.2                
#> [45] labeling_0.4.3              fastmap_1.1.1              
#> [47] grid_4.4.0                  colorspace_2.1-0           
#> [49] cli_3.6.2                   SparseArray_1.4.0          
#> [51] magrittr_2.0.3              S4Arrays_1.4.0             
#> [53] utf8_1.2.4                  withr_3.0.0                
#> [55] UCSC.utils_1.0.0            scales_1.3.0               
#> [57] rmarkdown_2.26              XVector_0.44.0             
#> [59] httr_1.4.7                  matrixStats_1.2.0          
#> [61] proxyC_0.3.4                evaluate_0.23              
#> [63] knitr_1.45                  GenomicRanges_1.56.0       
#> [65] IRanges_2.38.0              rlang_1.1.3                
#> [67] Rcpp_1.0.12                 glue_1.7.0                 
#> [69] BiocGenerics_0.50.0         jsonlite_1.8.8             
#> [71] R6_2.5.1                    MatrixGenerics_1.16.0      
#> [73] zlibbioc_1.50.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.