This vignette is made for users that are already familiar with the basic condiments workflow described in the first vignette. Here, we will show how to modify the default parameters for the first two steps of the workflow
# For analysis
library(condiments)
library(slingshot)
set.seed(21)
We rely on the same toy dataset as the first vignette
data("toy_dataset", package = "condiments")
df <- toy_dataset$sd
rd <- as.matrix(df[, c("Dim1", "Dim2")])
sds <- slingshot(rd, df$cl)
By default, the topologyTest function requires only two inputs, the sds
object and the condition
labels. To limit run time for the vignette, we also change the default number of permutations used to generate trajectories under the null by setting the rep
argument to \(10\) instead of the default \(100\). As such, the test statistics might be more variable.
top_res <- topologyTest(sds = sds, conditions = df$conditions, rep = 10)
## Generating permuted trajectories
## Running KS-mean test
knitr::kable(top_res)
method | thresh | statistic | p.value |
---|---|---|---|
KS_mean | 0.01 | 0 | 1 |
The topologyTest function can be relatively slow on large datasets. Moreover, when changing the method used to test the null hypothesis that a common trajectory should be fitted, the first permutation part of generating rep
trajectories under the null is identical. Therefore, we allow users to specify more than one method and one value of the threshold. Here, we will use both the Kolmogorov-Smirnov test test(Smirnov 1939) and the classifier-test(Lopez-Paz and Oquab 2016).
top_res <- topologyTest(sds = sds, conditions = df$conditions, rep = 10,
methods = c("KS_mean", "Classifier"),
threshs = c(0, .01, .05, .1))
## Generating permuted trajectories
## Running KS-mean test
## Running Classifier test
knitr::kable(top_res)
method | thresh | statistic | p.value |
---|---|---|---|
KS_mean | 0 | 0.0070000 | 1.0000000 |
KS_mean | 0.01 | 0.0000000 | 1.0000000 |
KS_mean | 0.05 | 0.0000000 | 1.0000000 |
KS_mean | 0.1 | 0.0000000 | 1.0000000 |
Classifier | 0 | 0.4150000 | 0.9999821 |
Classifier | 0.01 | 0.3800000 | 1.0000000 |
Classifier | 0.05 | 0.3333333 | 1.0000000 |
Classifier | 0.1 | 0.2833333 | 1.0000000 |
To see all methods avaible, use /tmp/RtmpFJHOA3/Rinst163f0a2edcb04/condiments/help/topologyTest and look at the methods
argument.
For all methods but the KS test, additional paramters can be specified, using a custom argument: args_classifier
, args_wass
or args_mmd
. See the help file for given test more information on those parameters. For example, since the default test based on the wasserstein distance and permutation test is quite slow, we can pass a fast
argument.
top_res <- topologyTest(sds = sds, conditions = df$conditions, rep = 10,
methods = "wasserstein_permutation",
args_wass = list(fast = TRUE, S = 100, iterations = 10^2))
## Generating permuted trajectories
## Running wassertsein permutation test
knitr::kable(top_res)
method | thresh | statistic | p.value |
---|---|---|---|
wasserstein_permutation | NA | 1.356887 | 0.85 |
For now, the first part of the topologyTest has been designed for parallelisation using the BiocParallel package. For example, to run with 4 cores, you can run the following command
library(BiocParallel)
BPPARAM <- bpparam()
BPPARAM$progressbar <- TRUE
BPPARAM$workers <- 4
top_res <- topologyTest(sds = sds, conditions = df$conditions, rep = 100,
parallel = TRUE, BPPARAM = BPPARAM)
knitr::kable(top_res)
The tests for the second test are much less compute-intensive, therefore there is no parallelisation. However, the other changes introduce in the previous section are still possible
prog_res <- progressionTest(sds, conditions = df$conditions)
knitr::kable(prog_res)
lineage | statistic | p.value |
---|---|---|
All | 5.506366 | 0 |
dif_res <- fateSelectionTest(sds, conditions = df$conditions)
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
##
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
knitr::kable(dif_res)
pair | statistic | p.value |
---|---|---|
1vs2 | 0.6611712 | 8e-07 |
prog_res <- progressionTest(sds, conditions = df$conditions, method = "Classifier")
knitr::kable(prog_res)
lineage | statistic | p.value |
---|---|---|
All | 0.6026126 | 0.0012246 |
dif_res <- fateSelectionTest(sds, conditions = df$conditions, thresh = .05)
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
##
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
knitr::kable(dif_res)
pair | statistic | p.value |
---|---|---|
1vs2 | 0.6301802 | 6.01e-05 |
prog_res <- progressionTest(sds, conditions = df$conditions, method = "Classifier",
args_classifier = list(method = "rf"))
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## Warning in randomForest.default(x, y, mtry = param$mtry, ...): invalid mtry:
## reset to within valid range
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
knitr::kable(prog_res)
lineage | statistic | p.value |
---|---|---|
All | 0.5890991 | 0.0043539 |
dif_res <- fateSelectionTest(sds, conditions = df$conditions)
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
##
## note: only 1 unique complexity parameters in default grid. Truncating the grid to 1 .
knitr::kable(dif_res)
pair | statistic | p.value |
---|---|---|
1vs2 | 0.6341441 | 3.41e-05 |
For all of the above procedures, it is important to note that we are making multiple comparisons. The p-values we obtain from these tests should be corrected for multiple testing, especially for trajectories with a large number of lineages.
That said, trajectory inference is often one of the last computational methods in a very long analysis pipeline (generally including gene-level quantification, gene filtering / feature selection, and dimensionality reduction). Hence, we strongly discourage the reader from putting too much faith in any p-value that comes out of this analysis. Such values may be useful suggestions, indicating particular features or cells for follow-up study, but should generally not be treated as meaningful statistical quantities.
If some commands and parameters are still unclear after going through this vignette, do not hesitate to open an issue on the condiments Github repository.
sessionInfo()
## R version 4.2.1 (2022-06-23)
## 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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] caret_6.0-93 lattice_0.20-45
## [3] viridis_0.6.2 viridisLite_0.4.1
## [5] RColorBrewer_1.1-3 ggplot2_3.3.6
## [7] tidyr_1.2.1 dplyr_1.0.10
## [9] slingshot_2.6.0 TrajectoryUtils_1.6.0
## [11] SingleCellExperiment_1.20.0 SummarizedExperiment_1.28.0
## [13] Biobase_2.58.0 GenomicRanges_1.50.0
## [15] GenomeInfoDb_1.34.0 IRanges_2.32.0
## [17] S4Vectors_0.36.0 BiocGenerics_0.44.0
## [19] MatrixGenerics_1.10.0 matrixStats_0.62.0
## [21] princurve_2.1.6 condiments_1.6.0
## [23] knitr_1.40
##
## loaded via a namespace (and not attached):
## [1] plyr_1.8.7 igraph_1.3.5
## [3] splines_4.2.1 BiocParallel_1.32.0
## [5] listenv_0.8.0 scater_1.26.0
## [7] digest_0.6.30 foreach_1.5.2
## [9] htmltools_0.5.3 fansi_1.0.3
## [11] magrittr_2.0.3 ScaledMatrix_1.6.0
## [13] tensor_1.5 doParallel_1.0.17
## [15] limma_3.54.0 recipes_1.0.2
## [17] globals_0.16.1 gower_1.0.0
## [19] hardhat_1.2.0 spatstat.sparse_3.0-0
## [21] colorspace_2.0-3 ggrepel_0.9.1
## [23] xfun_0.34 RCurl_1.98-1.9
## [25] jsonlite_1.8.3 distinct_1.10.0
## [27] spatstat_2.3-4 spatstat.data_3.0-0
## [29] survival_3.4-0 iterators_1.0.14
## [31] glue_1.6.2 polyclip_1.10-4
## [33] gtable_0.3.1 ipred_0.9-13
## [35] zlibbioc_1.44.0 XVector_0.38.0
## [37] DelayedArray_0.24.0 Ecume_0.9.1
## [39] BiocSingular_1.14.0 kernlab_0.9-31
## [41] future.apply_1.9.1 abind_1.4-5
## [43] scales_1.2.1 DBI_1.1.3
## [45] rngtools_1.5.2 spatstat.random_3.0-0
## [47] Rcpp_1.0.9 spatstat.core_2.4-4
## [49] rsvd_1.0.5 proxy_0.4-27
## [51] lava_1.7.0 prodlim_2019.11.13
## [53] ellipsis_0.3.2 farver_2.1.1
## [55] pkgconfig_2.0.3 scuttle_1.8.0
## [57] nnet_7.3-18 sass_0.4.2
## [59] deldir_1.0-6 utf8_1.2.2
## [61] labeling_0.4.2 tidyselect_1.2.0
## [63] rlang_1.0.6 reshape2_1.4.4
## [65] munsell_0.5.0 tools_4.2.1
## [67] cachem_1.0.6 cli_3.4.1
## [69] generics_0.1.3 evaluate_0.17
## [71] stringr_1.4.1 fastmap_1.1.0
## [73] yaml_2.3.6 goftest_1.2-3
## [75] ModelMetrics_1.2.2.2 transport_0.13-0
## [77] randomForest_4.7-1.1 purrr_0.3.5
## [79] RANN_2.6.1 pbapply_1.5-0
## [81] future_1.28.0 nlme_3.1-160
## [83] doRNG_1.8.2 sparseMatrixStats_1.10.0
## [85] compiler_4.2.1 beeswarm_0.4.0
## [87] e1071_1.7-12 spatstat.utils_3.0-1
## [89] spatstat.linnet_2.3-2 tibble_3.1.8
## [91] bslib_0.4.0 stringi_1.7.8
## [93] highr_0.9 Matrix_1.5-1
## [95] vctrs_0.5.0 pillar_1.8.1
## [97] lifecycle_1.0.3 spatstat.geom_3.0-3
## [99] jquerylib_0.1.4 BiocNeighbors_1.16.0
## [101] irlba_2.3.5.1 data.table_1.14.4
## [103] bitops_1.0-7 R6_2.5.1
## [105] gridExtra_2.3 vipor_0.4.5
## [107] parallelly_1.32.1 codetools_0.2-18
## [109] MASS_7.3-58.1 assertthat_0.2.1
## [111] withr_2.5.0 GenomeInfoDbData_1.2.9
## [113] mgcv_1.8-41 parallel_4.2.1
## [115] grid_4.2.1 rpart_4.1.19
## [117] beachmat_2.14.0 timeDate_4021.106
## [119] class_7.3-20 rmarkdown_2.17
## [121] DelayedMatrixStats_1.20.0 pROC_1.18.0
## [123] lubridate_1.8.0 ggbeeswarm_0.6.0
Lopez-Paz, David, and Maxime Oquab. 2016. “Revisiting Classifier Two-Sample Tests.” Arxiv, October, 1–15. http://arxiv.org/abs/1610.06545.
Smirnov, Nikolai V. 1939. “On the Estimation of the Discrepancy Between Empirical Curves of Distribution for Two Independent Samples.” Bull. Math. Univ. Moscou 2 (2): 3–14.