To install the developmental version of the package, run:
To install from Bioconductor:
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.
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)
)
)
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:
filter_genes: if TRUE, genes which have expression values of 0 in all cells will be removed.
variance_cutoff: remove all genes with a expression variance below the chosen cutoff.
type_abundance_cutoff: remove all cells, which belong to a cell type that appears less the the given amount.
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:
even: this creates samples, where all existing cell-types in the dataset appear in the same proportions. So using a dataset with 3 cell-types, this will simulate samples, where all cell-type fractions are 1/3. In order to still have a slight variation between cell type fractions, you can increase the balance_uniform_mirror_scenario
parameter (default to 0.01). Setting it to 0 will generate simulations with exactly the same cell type fractions.
random: this scenario will create random cell type fractions using all present types for each sample. The random sampling is based on the uniform distribution.
mirror_db: this scenario will mirror the exact fractions of cell types which are present in the provided dataset. If it consists of 20% T cells, 30% B cells and 50% NK cells, all simulated samples will mirror these fractions. Similar to the uniform scenario, you can add a small variation to these fractions with the balance_uniform_mirror_scenario
parameter.
weighted: here you need to set two additional parameters for the simulate_bulk()
function: weighted_cell_type
sets the cell-type you want to be over-representing and weighted_amount
sets the fraction of this cell-type. You could for example use B-cell
and 0.5
to create samples, where 50% are B-cells and the rest is filled randomly with other cell-types.
pure: this creates simulations of only one single cell-type. You have to provide the name of this cell-type with the pure_cell_type
parameter.
custom: here you are able to create your own set of cell-type fractions. When using this scenario, you additionally need to provide a dataframe in the custom_scenario_data
parameter, where each row represents one sample (therefore the number of rows need to match the nsamples
parameter). Each column has to represent one cell-type, which also occurs in the dataset and describes the fraction of this cell-type in a sample. The fractions per sample need to sum up to 1. An example can be seen here:
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
The simulation
object contains three named entries:
bulk
: a SummarizedExperiment object with the pseudo-bulk dataset(s) stored in the assays
. They can be accessed like this: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 478 497 453 491 500 481 448 517 496 477
#> gene_2 518 470 492 498 494 492 493 522 499 518
#> gene_3 537 466 479 453 477 446 486 492 517 461
#> gene_4 460 519 481 481 498 491 449 508 484 491
#> gene_5 514 505 494 502 486 525 483 510 508 555
#> gene_6 535 473 479 506 472 461 476 497 491 485
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 1060.6120 944.4662 1072.4908 925.7692 1005.7658 962.9583 920.2136
#> gene_2 1094.6510 1036.7257 1075.6130 1097.0251 1010.2219 960.6225 1070.5615
#> gene_3 940.6247 1110.4760 918.4488 1027.7519 911.2649 1043.9975 1047.5467
#> gene_4 1025.1214 920.1798 1000.6698 1006.5665 1028.0672 916.1504 1006.8949
#> gene_5 984.1692 900.6989 1006.5150 1051.2152 1022.9485 972.7817 1127.2808
#> gene_6 890.4835 987.7674 1013.4043 894.7700 936.3616 982.0254 947.8275
#>
#> gene_1 951.4859 922.7575 950.0023
#> gene_2 1009.9925 1017.3430 1024.6977
#> gene_3 979.4026 985.0034 975.6505
#> gene_4 1003.8724 948.9098 1062.2918
#> gene_5 1024.9971 949.6046 1085.0718
#> gene_6 1040.2142 964.7090 878.3131
If only a single matrix was given to the dataset initially, only one assay is filled.
cell_fractions
: a table where rows represent the simulated samples and columns represent the different simulated cell-types. The entries in the table store the specific cell-type fraction per sample.
scaling_vector
: a named list, with the used scaling value for each cell from the single cell dataset.
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:
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.2.2 (2022-10-31)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.5 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
#>
#> 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
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] SimBu_1.0.2
#>
#> loaded via a namespace (and not attached):
#> [1] SummarizedExperiment_1.28.0 tidyselect_1.2.0
#> [3] xfun_0.37 bslib_0.4.2
#> [5] purrr_1.0.1 lattice_0.20-45
#> [7] colorspace_2.1-0 vctrs_0.5.2
#> [9] generics_0.1.3 htmltools_0.5.4
#> [11] stats4_4.2.2 yaml_2.3.7
#> [13] utf8_1.2.3 rlang_1.0.6
#> [15] jquerylib_0.1.4 pillar_1.8.1
#> [17] glue_1.6.2 withr_2.5.0
#> [19] BiocParallel_1.32.5 RColorBrewer_1.1-3
#> [21] BiocGenerics_0.44.0 matrixStats_0.63.0
#> [23] GenomeInfoDbData_1.2.9 lifecycle_1.0.3
#> [25] zlibbioc_1.44.0 MatrixGenerics_1.10.0
#> [27] munsell_0.5.0 gtable_0.3.1
#> [29] proxyC_0.3.3 codetools_0.2-19
#> [31] evaluate_0.20 labeling_0.4.2
#> [33] Biobase_2.58.0 knitr_1.42
#> [35] IRanges_2.32.0 fastmap_1.1.1
#> [37] GenomeInfoDb_1.34.9 parallel_4.2.2
#> [39] fansi_1.0.4 highr_0.10
#> [41] Rcpp_1.0.10 scales_1.2.1
#> [43] cachem_1.0.7 DelayedArray_0.24.0
#> [45] S4Vectors_0.36.2 RcppParallel_5.1.7
#> [47] jsonlite_1.8.4 XVector_0.38.0
#> [49] farver_2.1.1 ggplot2_3.4.1
#> [51] digest_0.6.31 dplyr_1.1.0
#> [53] GenomicRanges_1.50.2 grid_4.2.2
#> [55] cli_3.6.0 tools_4.2.2
#> [57] bitops_1.0-7 magrittr_2.0.3
#> [59] sass_0.4.5 RCurl_1.98-1.10
#> [61] tibble_3.2.0 tidyr_1.3.0
#> [63] pkgconfig_2.0.3 Matrix_1.5-3
#> [65] data.table_1.14.8 sparseMatrixStats_1.10.0
#> [67] rmarkdown_2.20 R6_2.5.1
#> [69] compiler_4.2.2
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.