cellCellSimulate
functionscTensor 1.2.1
Here, we explain the way to generate CCI simulation data.
scTensor has a function cellCellSimulate
to generate the simulation data.
The simplest way to generate such data is cellCellSimulate
with default parameters.
suppressPackageStartupMessages(library("scTensor"))
sim <- cellCellSimulate()
## Getting the values of params...
## Setting random seed...
## Generating simulation data...
## Done!
This function internally generate the parameter sets by newCCSParams
,
and the values of the parameter can be changed, and specified as the input of cellCellSimulate
by users as follows.
# Default parameters
params <- newCCSParams()
str(params)
## Formal class 'CCSParams' [package "scTensor"] with 5 slots
## ..@ nGene : num 1000
## ..@ nCell : num [1:3] 50 50 50
## ..@ cciInfo:List of 4
## .. ..$ nPair: num 500
## .. ..$ CCI1 :List of 4
## .. .. ..$ LPattern: num [1:3] 1 0 0
## .. .. ..$ RPattern: num [1:3] 0 1 0
## .. .. ..$ nGene : num 50
## .. .. ..$ fc : chr "E10"
## .. ..$ CCI2 :List of 4
## .. .. ..$ LPattern: num [1:3] 0 1 0
## .. .. ..$ RPattern: num [1:3] 0 0 1
## .. .. ..$ nGene : num 50
## .. .. ..$ fc : chr "E10"
## .. ..$ CCI3 :List of 4
## .. .. ..$ LPattern: num [1:3] 0 0 1
## .. .. ..$ RPattern: num [1:3] 1 0 0
## .. .. ..$ nGene : num 50
## .. .. ..$ fc : chr "E10"
## ..@ lambda : num 1
## ..@ seed : num 1234
# Setting different parameters
# No. of genes : 1000
setParam(params, "nGene") <- 1000
# 3 cell types, 20 cells in each cell type
setParam(params, "nCell") <- c(20, 20, 20)
# Setting for Ligand-Receptor pair list
setParam(params, "cciInfo") <- list(
nPair=500, # Total number of L-R pairs
# 1st CCI
CCI1=list(
LPattern=c(1,0,0), # Only 1st cell type has this pattern
RPattern=c(0,1,0), # Only 2nd cell type has this pattern
nGene=50, # 50 pairs are generated as CCI1
fc="E10"), # Degree of differential expression (Fold Change)
# 2nd CCI
CCI2=list(
LPattern=c(0,1,0),
RPattern=c(0,0,1),
nGene=30,
fc="E100")
)
# Degree of Dropout
setParam(params, "lambda") <- 10
# Random number seed
setParam(params, "seed") <- 123
# Simulation data
sim <- cellCellSimulate(params)
## Getting the values of params...
## Setting random seed...
## Generating simulation data...
## Done!
The output object sim has some attributes as follows.
Firstly, sim$input contains a synthetic gene expression matrix. The size can be changed by nGene and nCell parameters described above.
dim(sim$input)
## [1] 1000 60
sim$input[1:2,1:3]
## Cell1 Cell2 Cell3
## Gene1 9105 2 0
## Gene2 4 37 850
Next, sim$LR contains a ligand-receptor (L-R) pair list. The size can be changed by nPair parameter of cciInfo, and the differentially expressed (DE) L-R pairs are saved in the upper side of this matrix. Here, two DE L-R patterns are specified as cciInfo, and each number of pairs is 50 and 30, respectively.
dim(sim$LR)
## [1] 500 2
sim$LR[1:10,]
## GENEID_L GENEID_R
## 1 Gene1 Gene81
## 2 Gene2 Gene82
## 3 Gene3 Gene83
## 4 Gene4 Gene84
## 5 Gene5 Gene85
## 6 Gene6 Gene86
## 7 Gene7 Gene87
## 8 Gene8 Gene88
## 9 Gene9 Gene89
## 10 Gene10 Gene90
sim$LR[46:55,]
## GENEID_L GENEID_R
## 46 Gene46 Gene126
## 47 Gene47 Gene127
## 48 Gene48 Gene128
## 49 Gene49 Gene129
## 50 Gene50 Gene130
## 51 Gene51 Gene131
## 52 Gene52 Gene132
## 53 Gene53 Gene133
## 54 Gene54 Gene134
## 55 Gene55 Gene135
sim$LR[491:500,]
## GENEID_L GENEID_R
## 491 Gene571 Gene991
## 492 Gene572 Gene992
## 493 Gene573 Gene993
## 494 Gene574 Gene994
## 495 Gene575 Gene995
## 496 Gene576 Gene996
## 497 Gene577 Gene997
## 498 Gene578 Gene998
## 499 Gene579 Gene999
## 500 Gene580 Gene1000
Finally, sim$celltypes contains a cell type vector. Since nCell is specified as “c(20, 20, 20)” described above, three cell types are generated.
length(sim$celltypes)
## [1] 60
head(sim$celltypes)
## Celltype1 Celltype1 Celltype1 Celltype1 Celltype1 Celltype1
## "Cell1" "Cell2" "Cell3" "Cell4" "Cell5" "Cell6"
table(names(sim$celltypes))
##
## Celltype1 Celltype2 Celltype3
## 20 20 20
## R version 3.6.1 (2019-07-05)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.10-bioc/R/lib/libRlapack.so
##
## 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
##
## attached base packages:
## [1] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] AnnotationHub_2.18.0
## [2] BiocFileCache_1.10.2
## [3] dbplyr_1.4.2
## [4] Homo.sapiens_1.3.1
## [5] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [6] org.Hs.eg.db_3.10.0
## [7] GO.db_3.10.0
## [8] OrganismDbi_1.28.0
## [9] GenomicFeatures_1.38.0
## [10] AnnotationDbi_1.48.0
## [11] MeSH.Mmu.eg.db_1.13.0
## [12] LRBase.Mmu.eg.db_1.2.0
## [13] MeSH.Hsa.eg.db_1.13.0
## [14] MeSHDbi_1.22.0
## [15] SingleCellExperiment_1.8.0
## [16] SummarizedExperiment_1.16.0
## [17] DelayedArray_0.12.0
## [18] BiocParallel_1.20.0
## [19] matrixStats_0.55.0
## [20] Biobase_2.46.0
## [21] GenomicRanges_1.38.0
## [22] GenomeInfoDb_1.22.0
## [23] IRanges_2.20.0
## [24] S4Vectors_0.24.0
## [25] BiocGenerics_0.32.0
## [26] scTensor_1.2.1
## [27] RSQLite_2.1.2
## [28] LRBase.Hsa.eg.db_1.2.0
## [29] LRBaseDbi_1.4.0
## [30] BiocStyle_2.14.0
##
## loaded via a namespace (and not attached):
## [1] rappdirs_0.3.1 rtracklayer_1.46.0
## [3] AnnotationForge_1.28.0 tidyr_1.0.0
## [5] ggplot2_3.2.1 acepack_1.4.1
## [7] bit64_0.9-7 knitr_1.25
## [9] data.table_1.12.6 rpart_4.1-15
## [11] RCurl_1.95-4.12 AnnotationFilter_1.10.0
## [13] cowplot_1.0.0 europepmc_0.3
## [15] bit_1.1-14 enrichplot_1.6.0
## [17] webshot_0.5.1 xml2_1.2.2
## [19] httpuv_1.5.2 assertthat_0.2.1
## [21] viridis_0.5.1 xfun_0.10
## [23] hms_0.5.2 evaluate_0.14
## [25] promises_1.1.0 TSP_1.1-7
## [27] progress_1.2.2 caTools_1.17.1.2
## [29] dendextend_1.12.0 Rgraphviz_2.30.0
## [31] igraph_1.2.4.1 DBI_1.0.0
## [33] htmlwidgets_1.5.1 MeSH.db_1.13.0
## [35] purrr_0.3.3 dplyr_0.8.3
## [37] backports_1.1.5 bookdown_0.14
## [39] annotate_1.64.0 biomaRt_2.42.0
## [41] vctrs_0.2.0 ensembldb_2.10.0
## [43] abind_1.4-5 ggforce_0.3.1
## [45] Gviz_1.30.0 triebeard_0.3.0
## [47] BSgenome_1.54.0 checkmate_1.9.4
## [49] GenomicAlignments_1.22.0 gclus_1.3.2
## [51] fdrtool_1.2.15 prettyunits_1.0.2
## [53] cluster_2.1.0 DOSE_3.12.0
## [55] dotCall64_1.0-0 lazyeval_0.2.2
## [57] crayon_1.3.4 genefilter_1.68.0
## [59] pkgconfig_2.0.3 tweenr_1.0.1
## [61] ProtGenerics_1.18.0 seriation_1.2-8
## [63] nnet_7.3-12 rlang_0.4.1
## [65] lifecycle_0.1.0 meshr_1.22.0
## [67] registry_0.5-1 MeSH.PCR.db_1.13.0
## [69] rTensor_1.4 GOstats_2.52.0
## [71] dichromat_2.0-0 polyclip_1.10-0
## [73] graph_1.64.0 Matrix_1.2-17
## [75] urltools_1.7.3 base64enc_0.1-3
## [77] ggridges_0.5.1 viridisLite_0.3.0
## [79] MeSH.AOR.db_1.13.0 bitops_1.0-6
## [81] visNetwork_2.0.8 KernSmooth_2.23-16
## [83] spam_2.4-0 MeSH.Bsu.168.eg.db_1.13.0
## [85] Biostrings_2.54.0 blob_1.2.0
## [87] stringr_1.4.0 qvalue_2.18.0
## [89] nnTensor_1.0.2 gridGraphics_0.4-1
## [91] reactome.db_1.70.0 scales_1.0.0
## [93] graphite_1.32.0 memoise_1.1.0
## [95] GSEABase_1.48.0 magrittr_1.5
## [97] plyr_1.8.4 gplots_3.0.1.1
## [99] gdata_2.18.0 zlibbioc_1.32.0
## [101] compiler_3.6.1 RColorBrewer_1.1-2
## [103] plotrix_3.7-6 Rsamtools_2.2.1
## [105] XVector_0.26.0 Category_2.52.1
## [107] MeSH.Aca.eg.db_1.13.0 htmlTable_1.13.2
## [109] Formula_1.2-3 MASS_7.3-51.4
## [111] tidyselect_0.2.5 stringi_1.4.3
## [113] highr_0.8 MeSH.Syn.eg.db_1.13.0
## [115] yaml_2.2.0 GOSemSim_2.12.0
## [117] askpass_1.1 latticeExtra_0.6-28
## [119] ggrepel_0.8.1 grid_3.6.1
## [121] VariantAnnotation_1.32.0 fastmatch_1.1-0
## [123] tools_3.6.1 rstudioapi_0.10
## [125] foreach_1.4.7 foreign_0.8-72
## [127] tagcloud_0.6 outliers_0.14
## [129] gridExtra_2.3 farver_1.1.0
## [131] ggraph_2.0.0 rvcheck_0.1.6
## [133] digest_0.6.22 BiocManager_1.30.9
## [135] shiny_1.4.0 Rcpp_1.0.3
## [137] BiocVersion_3.10.1 later_1.0.0
## [139] httr_1.4.1 cummeRbund_2.28.0
## [141] biovizBase_1.34.0 colorspace_1.4-1
## [143] XML_3.98-1.20 splines_3.6.1
## [145] fields_9.9 RBGL_1.62.1
## [147] graphlayouts_0.5.0 ggplotify_0.0.4
## [149] plotly_4.9.1 xtable_1.8-4
## [151] jsonlite_1.6 heatmaply_0.16.0
## [153] tidygraph_1.1.2 zeallot_0.1.0
## [155] R6_2.4.0 Hmisc_4.3-0
## [157] pillar_1.4.2 htmltools_0.4.0
## [159] mime_0.7 glue_1.3.1
## [161] fastmap_1.0.1 interactiveDisplayBase_1.24.0
## [163] codetools_0.2-16 maps_3.3.0
## [165] fgsea_1.12.0 lattice_0.20-38
## [167] tibble_2.1.3 curl_4.2
## [169] gtools_3.8.1 ReactomePA_1.30.0
## [171] misc3d_0.8-4 openssl_1.4.1
## [173] survival_3.1-7 rmarkdown_1.16
## [175] munsell_0.5.0 DO.db_2.9
## [177] GenomeInfoDbData_1.2.2 plot3D_1.1.1
## [179] iterators_1.0.12 reshape2_1.4.3
## [181] gtable_0.3.0