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(rpois(3e5, 5), ncol = 300), sparse = TRUE)
tpm <- Matrix::Matrix(matrix(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:

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 551 522 519 524 542 521 580 536 518 563
#> gene_2 461 490 459 489 480 488 503 492 482 499
#> gene_3 471 509 514 485 468 513 495 479 453 518
#> gene_4 518 483 475 518 476 502 519 520 492 488
#> gene_5 469 458 487 479 461 468 513 488 495 485
#> gene_6 528 506 502 519 480 530 490 530 493 536
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 1021.2603  852.7526  868.7378  879.7911  956.0529  914.8028  896.2281
#> gene_2  998.0065  938.8738  937.9736 1029.0375  952.5837 1060.5080 1000.1050
#> gene_3  909.2224  973.8122  851.0622  984.5271  978.6644  920.2616  984.4385
#> gene_4  939.2274 1015.6191  970.8708 1049.2400  943.0289  994.3796 1000.1460
#> gene_5  985.5078 1011.7561 1016.2642 1021.3856 1104.9176 1106.9404 1150.7218
#> gene_6 1002.4917 1037.0390 1016.5300  928.4179 1008.1346 1033.4585  931.2809
#>                                     
#> gene_1  984.2403 1015.4648  959.1034
#> gene_2  967.8603 1078.3930  948.5409
#> gene_3  980.1063 1038.9000  929.2623
#> gene_4 1031.0791  977.0814  991.7372
#> gene_5  994.9728  970.9883 1038.2173
#> gene_6 1006.9099  981.1736 1007.2674

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.

sessionInfo()
#> R version 4.3.0 RC (2023-04-18 r84287)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.2 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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.3.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] SummarizedExperiment_1.31.1 gtable_0.3.3               
#>  [3] xfun_0.39                   bslib_0.4.2                
#>  [5] ggplot2_3.4.2               Biobase_2.61.0             
#>  [7] lattice_0.21-8              vctrs_0.6.2                
#>  [9] tools_4.3.0                 bitops_1.0-7               
#> [11] generics_0.1.3              stats4_4.3.0               
#> [13] parallel_4.3.0              tibble_3.2.1               
#> [15] fansi_1.0.4                 highr_0.10                 
#> [17] pkgconfig_2.0.3             Matrix_1.5-4.1             
#> [19] data.table_1.14.8           RColorBrewer_1.1-3         
#> [21] S4Vectors_0.39.1            sparseMatrixStats_1.13.0   
#> [23] RcppParallel_5.1.7          lifecycle_1.0.3            
#> [25] GenomeInfoDbData_1.2.10     farver_2.1.1               
#> [27] compiler_4.3.0              munsell_0.5.0              
#> [29] codetools_0.2-19            GenomeInfoDb_1.37.1        
#> [31] htmltools_0.5.5             sass_0.4.6                 
#> [33] RCurl_1.98-1.12             yaml_2.3.7                 
#> [35] pillar_1.9.0                crayon_1.5.2               
#> [37] jquerylib_0.1.4             tidyr_1.3.0                
#> [39] BiocParallel_1.35.1         DelayedArray_0.27.3        
#> [41] cachem_1.0.8                tidyselect_1.2.0           
#> [43] digest_0.6.31               dplyr_1.1.2                
#> [45] purrr_1.0.1                 labeling_0.4.2             
#> [47] fastmap_1.1.1               grid_4.3.0                 
#> [49] colorspace_2.1-0            cli_3.6.1                  
#> [51] SparseArray_1.1.5           magrittr_2.0.3             
#> [53] S4Arrays_1.1.4              utf8_1.2.3                 
#> [55] withr_2.5.0                 scales_1.2.1               
#> [57] rmarkdown_2.21              XVector_0.41.1             
#> [59] matrixStats_0.63.0          proxyC_0.3.3               
#> [61] evaluate_0.21               knitr_1.42                 
#> [63] GenomicRanges_1.53.1        IRanges_2.35.1             
#> [65] rlang_1.1.1                 Rcpp_1.0.10                
#> [67] glue_1.6.2                  BiocGenerics_0.47.0        
#> [69] jsonlite_1.8.4              R6_2.5.1                   
#> [71] MatrixGenerics_1.13.0       zlibbioc_1.47.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.