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 480 538 499 500 489 506 473 471 498 498
#> gene_2 520 515 528 483 554 511 526 551 507 508
#> gene_3 565 515 536 466 496 529 498 511 586 494
#> gene_4 504 464 514 501 498 498 508 548 459 535
#> gene_5 483 447 531 488 463 472 498 454 462 445
#> gene_6 522 481 483 489 551 511 482 462 487 479
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 917.3090 847.6256 929.4734 1050.3268 954.6312 994.4291 978.7278
#> gene_2 1007.0218 1003.8761 931.8601 963.2853 954.3600 925.9787 987.2465
#> gene_3 998.9605 1002.6937 991.5607 993.3463 960.9643 1023.5216 978.8270
#> gene_4 1032.9990 987.0389 960.1971 1033.8829 1040.3039 977.9487 996.1787
#> gene_5 954.1748 891.3941 999.4426 962.6647 910.8827 1006.8530 1002.0658
#> gene_6 961.0549 944.8124 981.4050 856.4421 927.2970 1069.7504 927.3532
#>
#> gene_1 930.8161 955.1814 972.8222
#> gene_2 990.2235 1018.3494 915.1991
#> gene_3 893.4930 1074.5208 922.6471
#> gene_4 967.7040 943.6689 1044.1654
#> gene_5 925.8495 951.5725 918.9086
#> gene_6 978.2388 1007.4616 967.2921
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.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.2 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.9.0 rmarkdown_2.29
#>
#> loaded via a namespace (and not attached):
#> [1] sass_0.4.9 generics_0.1.3
#> [3] tidyr_1.3.1 SparseArray_1.7.6
#> [5] lattice_0.22-6 magrittr_2.0.3
#> [7] digest_0.6.37 RColorBrewer_1.1-3
#> [9] evaluate_1.0.3 sparseMatrixStats_1.19.0
#> [11] grid_4.4.2 fastmap_1.2.0
#> [13] jsonlite_1.9.0 Matrix_1.7-2
#> [15] GenomeInfoDb_1.43.4 proxyC_0.4.1
#> [17] httr_1.4.7 purrr_1.0.4
#> [19] scales_1.3.0 UCSC.utils_1.3.1
#> [21] codetools_0.2-20 jquerylib_0.1.4
#> [23] abind_1.4-8 cli_3.6.4
#> [25] rlang_1.1.5 crayon_1.5.3
#> [27] XVector_0.47.2 Biobase_2.67.0
#> [29] munsell_0.5.1 withr_3.0.2
#> [31] cachem_1.1.0 DelayedArray_0.33.6
#> [33] yaml_2.3.10 S4Arrays_1.7.3
#> [35] tools_4.4.2 parallel_4.4.2
#> [37] BiocParallel_1.41.2 dplyr_1.1.4
#> [39] colorspace_2.1-1 ggplot2_3.5.1
#> [41] GenomeInfoDbData_1.2.13 SummarizedExperiment_1.37.0
#> [43] BiocGenerics_0.53.6 buildtools_1.0.0
#> [45] vctrs_0.6.5 R6_2.6.1
#> [47] matrixStats_1.5.0 stats4_4.4.2
#> [49] lifecycle_1.0.4 S4Vectors_0.45.4
#> [51] IRanges_2.41.3 pkgconfig_2.0.3
#> [53] gtable_0.3.6 pillar_1.10.1
#> [55] bslib_0.9.0 data.table_1.17.0
#> [57] glue_1.8.0 Rcpp_1.0.14
#> [59] tidyselect_1.2.1 tibble_3.2.1
#> [61] xfun_0.51 GenomicRanges_1.59.1
#> [63] sys_3.4.3 MatrixGenerics_1.19.1
#> [65] knitr_1.49 farver_2.1.2
#> [67] htmltools_0.5.8.1 labeling_0.4.3
#> [69] maketools_1.3.2 compiler_4.4.2