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 480 515 479 470 472 424 501 479 508 486
#> gene_2 520 504 429 479 467 477 475 452 520 461
#> gene_3 508 473 489 520 495 507 482 439 472 508
#> gene_4 493 526 518 477 480 549 487 521 505 511
#> gene_5 488 442 536 573 494 532 521 486 509 528
#> gene_6 494 490 502 502 496 524 500 498 460 455
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  851.5997  961.0056  982.7224 1042.0401 1033.6675 1031.8659  980.7262
#> gene_2 1017.9247 1073.0034  996.9593 1049.7167  976.5608 1113.4651 1117.3500
#> gene_3  910.8170 1011.2159  988.8510  924.0470 1004.9422  855.6023  886.6514
#> gene_4  960.9436  978.2517 1013.9056 1051.1849 1024.5927  986.0201 1002.8083
#> gene_5  926.9408  988.1291  918.0845  923.4818 1042.0072 1032.4121  959.1604
#> gene_6  949.3165  930.2015 1012.0027  881.1551  946.0135  985.1899  878.7972
#>                                     
#> gene_1  968.7519 1012.5160  995.4244
#> gene_2 1051.8021 1008.0397 1042.6429
#> gene_3  896.7873 1065.5241  819.3697
#> gene_4  945.1486 1032.6643 1004.3060
#> gene_5 1005.8687  944.7735 1006.7968
#> gene_6  955.7845  898.5142  980.0798

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.1 (2024-06-14)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.20-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.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.8.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] SummarizedExperiment_1.36.0 gtable_0.3.6               
#>  [3] xfun_0.48                   bslib_0.8.0                
#>  [5] ggplot2_3.5.1               Biobase_2.66.0             
#>  [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.44.0           
#> [21] sparseMatrixStats_1.18.0    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.42.0        
#> [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.40.0        
#> [37] DelayedArray_0.32.0         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.6.0           magrittr_2.0.3             
#> [51] S4Arrays_1.6.0              utf8_1.2.4                 
#> [53] withr_3.0.2                 UCSC.utils_1.2.0           
#> [55] scales_1.3.0                rmarkdown_2.28             
#> [57] XVector_0.46.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.58.0        IRanges_2.40.0             
#> [65] rlang_1.1.4                 Rcpp_1.0.13                
#> [67] glue_1.8.0                  BiocGenerics_0.52.0        
#> [69] jsonlite_1.8.9              R6_2.5.1                   
#> [71] MatrixGenerics_1.18.0       zlibbioc_1.52.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.