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(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)
)
)
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: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 504 518 503 549 550 538 501 514 482 474
#> gene_2 460 455 460 492 500 490 480 501 490 525
#> gene_3 515 505 534 499 537 528 499 484 491 524
#> gene_4 432 486 466 476 489 511 441 508 442 458
#> gene_5 501 480 456 484 475 505 502 533 551 499
#> gene_6 462 489 462 496 500 470 437 516 491 477
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 1081.0202 1041.4781 989.2049 946.5158 952.9356 1010.4092 1165.0720
#> gene_2 1117.7805 966.4502 1067.3581 1095.3402 1039.8625 953.3147 1070.9810
#> gene_3 1098.9304 1064.2898 1015.6385 1040.7121 1029.9330 976.8305 1046.8770
#> gene_4 1065.4897 1010.4520 999.5653 1037.7010 927.7590 1071.7857 1041.7859
#> gene_5 1021.2263 986.0555 1001.3601 933.7968 986.1226 952.7440 1003.3185
#> gene_6 956.7748 972.7914 1028.7646 1006.2795 946.8972 1015.6459 993.3265
#>
#> gene_1 973.3584 894.2724 939.9585
#> gene_2 1010.6160 1003.8041 980.0991
#> gene_3 952.1582 1022.9981 1138.2899
#> gene_4 976.7598 981.0899 1069.3504
#> gene_5 974.0344 965.5192 940.6701
#> gene_6 1021.6327 982.1775 982.0193
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.
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: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] SimBu_1.7.0 rmarkdown_2.28
#>
#> loaded via a namespace (and not attached):
#> [1] SummarizedExperiment_1.35.2 gtable_0.3.5
#> [3] xfun_0.47 bslib_0.8.0
#> [5] ggplot2_3.5.1 Biobase_2.65.1
#> [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-0 data.table_1.16.0
#> [19] RColorBrewer_1.1-3 S4Vectors_0.43.2
#> [21] sparseMatrixStats_1.17.2 lifecycle_1.0.4
#> [23] GenomeInfoDbData_1.2.12 farver_2.1.2
#> [25] compiler_4.4.1 munsell_0.5.1
#> [27] codetools_0.2-20 GenomeInfoDb_1.41.1
#> [29] htmltools_0.5.8.1 sys_3.4.2
#> [31] buildtools_1.0.0 sass_0.4.9
#> [33] yaml_2.3.10 pillar_1.9.0
#> [35] crayon_1.5.3 jquerylib_0.1.4
#> [37] tidyr_1.3.1 BiocParallel_1.39.0
#> [39] DelayedArray_0.31.12 cachem_1.1.0
#> [41] abind_1.4-8 tidyselect_1.2.1
#> [43] digest_0.6.37 dplyr_1.1.4
#> [45] purrr_1.0.2 labeling_0.4.3
#> [47] maketools_1.3.0 fastmap_1.2.0
#> [49] grid_4.4.1 colorspace_2.1-1
#> [51] cli_3.6.3 SparseArray_1.5.41
#> [53] magrittr_2.0.3 S4Arrays_1.5.10
#> [55] utf8_1.2.4 withr_3.0.1
#> [57] UCSC.utils_1.1.0 scales_1.3.0
#> [59] XVector_0.45.0 httr_1.4.7
#> [61] matrixStats_1.4.1 proxyC_0.4.1
#> [63] evaluate_1.0.0 knitr_1.48
#> [65] GenomicRanges_1.57.1 IRanges_2.39.2
#> [67] rlang_1.1.4 Rcpp_1.0.13
#> [69] glue_1.7.0 BiocGenerics_0.51.2
#> [71] jsonlite_1.8.9 R6_2.5.1
#> [73] MatrixGenerics_1.17.0 zlibbioc_1.51.1