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 510 472 536 509 521 513 504 484 479 508
#> gene_2 453 468 486 517 429 469 481 449 472 481
#> gene_3 495 525 516 462 486 477 522 533 485 491
#> gene_4 529 536 500 494 510 491 485 490 489 488
#> gene_5 498 488 548 499 542 482 535 489 536 493
#> gene_6 529 472 520 476 529 472 496 516 517 505
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 1001.4480 1019.0515  986.3953 1043.6026 1054.4960 1046.4670  976.4499
#> gene_2  999.9848  968.6849 1045.7449 1052.5533 1027.5327 1001.6675  906.5013
#> gene_3  967.5137 1067.0218 1145.3987  988.4225 1046.3039  992.8773 1206.9548
#> gene_4  923.7367  955.2821  935.6155 1096.6408  998.2858  997.7081  973.6348
#> gene_5  857.6562  931.5557 1025.4675  968.5354 1036.4539  975.1879 1085.3867
#> gene_6 1108.1948 1071.6095 1041.1932 1029.2884 1021.6934 1034.1712  906.3912
#>                                     
#> gene_1  983.5070 1005.4866 1054.1475
#> gene_2  999.1710 1057.0826 1070.5088
#> gene_3 1053.7727 1165.0212 1049.0414
#> gene_4 1043.6898 1059.3759  967.3083
#> gene_5  963.4586  958.5773  950.4356
#> gene_6 1119.3173 1012.4783 1002.4219

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.3.2 Patched (2023-11-01 r85457)
#> Platform: x86_64-apple-darwin20 (64-bit)
#> Running under: macOS Monterey 12.7.1
#> 
#> Matrix products: default
#> BLAS:   /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRblas.0.dylib 
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRlapack.dylib;  LAPACK version 3.11.0
#> 
#> 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.4.3
#> 
#> loaded via a namespace (and not attached):
#>  [1] SummarizedExperiment_1.32.0 gtable_0.3.4               
#>  [3] xfun_0.41                   bslib_0.6.1                
#>  [5] ggplot2_3.4.4               Biobase_2.62.0             
#>  [7] lattice_0.22-5              vctrs_0.6.5                
#>  [9] tools_4.3.2                 bitops_1.0-7               
#> [11] generics_0.1.3              stats4_4.3.2               
#> [13] parallel_4.3.2              tibble_3.2.1               
#> [15] fansi_1.0.6                 highr_0.10                 
#> [17] pkgconfig_2.0.3             Matrix_1.6-4               
#> [19] data.table_1.14.10          RColorBrewer_1.1-3         
#> [21] S4Vectors_0.40.2            sparseMatrixStats_1.14.0   
#> [23] RcppParallel_5.1.7          lifecycle_1.0.4            
#> [25] GenomeInfoDbData_1.2.11     farver_2.1.1               
#> [27] compiler_4.3.2              munsell_0.5.0              
#> [29] codetools_0.2-19            GenomeInfoDb_1.38.3        
#> [31] htmltools_0.5.7             sass_0.4.8                 
#> [33] RCurl_1.98-1.13             yaml_2.3.8                 
#> [35] pillar_1.9.0                crayon_1.5.2               
#> [37] jquerylib_0.1.4             tidyr_1.3.0                
#> [39] BiocParallel_1.36.0         DelayedArray_0.28.0        
#> [41] cachem_1.0.8                abind_1.4-5                
#> [43] tidyselect_1.2.0            digest_0.6.33              
#> [45] dplyr_1.1.4                 purrr_1.0.2                
#> [47] labeling_0.4.3              fastmap_1.1.1              
#> [49] grid_4.3.2                  colorspace_2.1-0           
#> [51] cli_3.6.2                   SparseArray_1.2.2          
#> [53] magrittr_2.0.3              S4Arrays_1.2.0             
#> [55] utf8_1.2.4                  withr_2.5.2                
#> [57] scales_1.3.0                rmarkdown_2.25             
#> [59] XVector_0.42.0              matrixStats_1.2.0          
#> [61] proxyC_0.3.4                evaluate_0.23              
#> [63] knitr_1.45                  GenomicRanges_1.54.1       
#> [65] IRanges_2.36.0              rlang_1.1.2                
#> [67] Rcpp_1.0.11                 glue_1.6.2                 
#> [69] BiocGenerics_0.48.1         jsonlite_1.8.8             
#> [71] R6_2.5.1                    MatrixGenerics_1.14.0      
#> [73] zlibbioc_1.48.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.