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 482 534 491 499 495 507 542 458 546 564
#> gene_2 459 519 493 509 494 529 479 464 506 461
#> gene_3 536 527 506 441 548 518 477 515 516 554
#> gene_4 487 529 470 501 481 547 506 524 497 551
#> gene_5 577 526 568 538 561 522 532 557 591 526
#> gene_6 505 500 488 477 517 527 483 517 509 516
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 1033.8223 1027.2117 1076.3390 967.9239 966.5618 1057.1133 1027.6049
#> gene_2 1087.6159 1030.5861 975.0244 1028.7160 899.5300 1044.8301 1064.0081
#> gene_3 1011.8519 994.6910 1028.7709 976.7467 880.6966 947.5179 951.2934
#> gene_4 892.3603 1004.6317 1018.6388 975.8361 1004.4825 1051.9136 956.5659
#> gene_5 1119.4544 1018.7277 1002.1062 1044.6554 1024.8747 923.9191 1112.8411
#> gene_6 887.5217 938.2594 955.7645 1059.9833 952.8989 960.0112 916.8834
#>
#> gene_1 1123.6891 1023.6998 1058.0969
#> gene_2 1058.7990 1009.8753 1024.6929
#> gene_3 1046.5420 993.3665 936.7591
#> gene_4 956.3624 991.8566 1025.3053
#> gene_5 1125.5781 949.4014 1004.2567
#> gene_6 946.3113 920.0235 993.7817
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.3.0 RC (2023-04-13 r84269)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.2 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.17-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.2.0
#>
#> loaded via a namespace (and not attached):
#> [1] sass_0.4.5 utf8_1.2.3
#> [3] generics_0.1.3 tidyr_1.3.0
#> [5] bitops_1.0-7 lattice_0.21-8
#> [7] digest_0.6.31 magrittr_2.0.3
#> [9] RColorBrewer_1.1-3 evaluate_0.20
#> [11] sparseMatrixStats_1.12.0 grid_4.3.0
#> [13] fastmap_1.1.1 jsonlite_1.8.4
#> [15] Matrix_1.5-4 GenomeInfoDb_1.36.0
#> [17] proxyC_0.3.3 purrr_1.0.1
#> [19] fansi_1.0.4 scales_1.2.1
#> [21] codetools_0.2-19 jquerylib_0.1.4
#> [23] cli_3.6.1 rlang_1.1.0
#> [25] XVector_0.40.0 Biobase_2.60.0
#> [27] munsell_0.5.0 withr_2.5.0
#> [29] cachem_1.0.7 DelayedArray_0.26.0
#> [31] yaml_2.3.7 tools_4.3.0
#> [33] parallel_4.3.0 BiocParallel_1.34.0
#> [35] dplyr_1.1.2 colorspace_2.1-0
#> [37] ggplot2_3.4.2 GenomeInfoDbData_1.2.10
#> [39] SummarizedExperiment_1.30.0 BiocGenerics_0.46.0
#> [41] vctrs_0.6.2 R6_2.5.1
#> [43] matrixStats_0.63.0 stats4_4.3.0
#> [45] lifecycle_1.0.3 zlibbioc_1.46.0
#> [47] S4Vectors_0.38.0 IRanges_2.34.0
#> [49] pkgconfig_2.0.3 gtable_0.3.3
#> [51] RcppParallel_5.1.7 bslib_0.4.2
#> [53] pillar_1.9.0 data.table_1.14.8
#> [55] glue_1.6.2 Rcpp_1.0.10
#> [57] highr_0.10 tidyselect_1.2.0
#> [59] tibble_3.2.1 xfun_0.39
#> [61] GenomicRanges_1.52.0 MatrixGenerics_1.12.0
#> [63] knitr_1.42 farver_2.1.1
#> [65] htmltools_0.5.5 labeling_0.4.2
#> [67] rmarkdown_2.21 compiler_4.3.0
#> [69] RCurl_1.98-1.12
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.