mist (Methylation Inference for Single-cell along
Trajectory) is an R package for differential methylation (DM) analysis
of single-cell DNA methylation (scDNAm) data. The package employs a
Bayesian approach to model methylation changes along pseudotime and
estimates developmental-stage-specific biological variations. It
supports both single-group and two-group analyses, enabling users to
identify genomic features exhibiting temporal changes in methylation
levels or different methylation patterns between groups.
This vignette demonstrates how to use mist for: 1.
Single-group analysis. 2. Two-group analysis.
To install the latest version of mist, run the following
commands:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
# Install mist from GitHub
BiocManager::install("https://github.com/dxd429/mist")From Bioconductor:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("mist")To view the package vignette in HTML format, run the following lines in R:
In this section, we will estimate parameters and perform differential methylation analysis using single-group data.
Here we load the example data from GSE121708.
estiParam# Estimate parameters for single-group
Dat_sce <- estiParam(
Dat_sce = Dat_sce,
Dat_name = "Methy_level_group1",
ptime_name = "pseudotime"
)
# Check the output
head(rowData(Dat_sce)$mist_pars)## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.217390 -0.46226325 0.38944755 0.27614882 0.07305266
## ENSMUSG00000000003 1.650459 1.93019825 2.26120361 -1.78927609 -2.79609101
## ENSMUSG00000000028 1.281953 -0.03635691 0.13788683 0.04751644 -0.01347517
## ENSMUSG00000000037 1.029595 -3.78230882 10.48311799 -4.31565453 -2.39659535
## ENSMUSG00000000049 1.021413 -0.07113480 0.07937013 0.08449837 0.06680611
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.306936 13.123691 3.245110 2.237024
## ENSMUSG00000000003 25.511253 2.273001 5.832437 8.652112
## ENSMUSG00000000028 7.768132 7.829904 2.937132 2.292439
## ENSMUSG00000000037 8.779129 11.899942 7.693535 1.997782
## ENSMUSG00000000049 6.010624 7.935351 2.772729 1.242864
dmSingle# Perform differential methylation analysis for the single-group
Dat_sce <- dmSingle(Dat_sce)
# View the top genomic features with drastic methylation changes
head(rowData(Dat_sce)$mist_int)## ENSMUSG00000000037 ENSMUSG00000000003 ENSMUSG00000000001 ENSMUSG00000000049
## 0.053262720 0.031822853 0.013143407 0.007246521
## ENSMUSG00000000028
## 0.005785958
plotGene# Produce scatterplot with fitted curve of a specific gene
library(ggplot2)
plotGene(Dat_sce = Dat_sce,
Dat_name = "Methy_level_group1",
ptime_name = "pseudotime",
gene_name = "ENSMUSG00000000037")In this section, we will estimate parameters and perform DM analysis using data from two phenotypic groups.
estiParam# Estimate parameters for both groups
Dat_sce_g1 <- estiParam(
Dat_sce = Dat_sce_g1,
Dat_name = "Methy_level_group1",
ptime_name = "pseudotime"
)
Dat_sce_g2 <- estiParam(
Dat_sce = Dat_sce_g2,
Dat_name = "Methy_level_group2",
ptime_name = "pseudotime"
)
# Check the output
head(rowData(Dat_sce_g1)$mist_pars, n = 3)## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.253698 -0.76574991 0.6644385 0.37726119 2.379373e-05
## ENSMUSG00000000003 1.653018 1.90330028 2.3744679 -1.78950434 -2.883644e+00
## ENSMUSG00000000028 1.280757 -0.03122829 0.1429273 0.04818946 -1.070677e-02
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.521048 13.646851 3.378493 1.831297
## ENSMUSG00000000003 25.664334 2.936989 5.531436 8.720741
## ENSMUSG00000000028 7.864438 6.555444 3.088594 2.562817
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9049074 -0.1432088 3.1243501 -1.48762646 -1.6788618
## ENSMUSG00000000003 -0.8348991 -2.6803680 7.3362687 -3.53179360 -1.0982579
## ENSMUSG00000000028 2.3427540 0.0409317 0.5633593 -0.06590714 -0.4400482
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.682541 5.269209 3.740168 1.193827
## ENSMUSG00000000003 6.197854 10.796634 4.331385 2.794119
## ENSMUSG00000000028 12.186660 5.305368 3.739195 2.973789
dmTwoGroups# Perform DM analysis to compare the two groups
dm_results <- dmTwoGroups(
Dat_sce_g1 = Dat_sce_g1,
Dat_sce_g2 = Dat_sce_g2
)
# View the top genomic features with different temporal patterns between groups
head(dm_results)## ENSMUSG00000000037 ENSMUSG00000000003 ENSMUSG00000000001 ENSMUSG00000000049
## 0.048712299 0.042169292 0.026100611 0.010987740
## ENSMUSG00000000028
## 0.003689232
mist provides a comprehensive suite of tools for
analyzing scDNAm data along pseudotime, whether you are working with a
single group or comparing two phenotypic groups. With the combination of
Bayesian modeling and differential methylation analysis,
mist is a powerful tool for identifying significant genomic
features in scDNAm data.
## R version 4.5.2 (2025-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] ggplot2_4.0.1 SingleCellExperiment_1.33.0
## [3] SummarizedExperiment_1.41.0 Biobase_2.71.0
## [5] GenomicRanges_1.63.1 Seqinfo_1.1.0
## [7] IRanges_2.45.0 S4Vectors_0.49.0
## [9] BiocGenerics_0.57.0 generics_0.1.4
## [11] MatrixGenerics_1.23.0 matrixStats_1.5.0
## [13] mist_1.3.1 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.1.4 farver_2.1.2
## [4] Biostrings_2.79.4 S7_0.2.1 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.17 GenomicAlignments_1.47.0
## [10] XML_3.99-0.20 digest_0.6.39 lifecycle_1.0.5
## [13] survival_3.8-3 magrittr_2.0.4 compiler_4.5.2
## [16] rlang_1.1.6 sass_0.4.10 tools_4.5.2
## [19] yaml_2.3.12 rtracklayer_1.71.3 knitr_1.51
## [22] S4Arrays_1.11.1 labeling_0.4.3 curl_7.0.0
## [25] DelayedArray_0.37.0 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.45.0 withr_3.0.2 sys_3.4.3
## [31] grid_4.5.2 scales_1.4.0 MASS_7.3-65
## [34] mcmc_0.9-8 cli_3.6.5 mvtnorm_1.3-3
## [37] rmarkdown_2.30 crayon_1.5.3 httr_1.4.7
## [40] rjson_0.2.23 cachem_1.1.0 splines_4.5.2
## [43] parallel_4.5.2 BiocManager_1.30.27 XVector_0.51.0
## [46] restfulr_0.0.16 vctrs_0.6.5 Matrix_1.7-4
## [49] jsonlite_2.0.0 SparseM_1.84-2 carData_3.0-5
## [52] car_3.1-3 MCMCpack_1.7-1 Formula_1.2-5
## [55] maketools_1.3.2 jquerylib_0.1.4 glue_1.8.0
## [58] codetools_0.2-20 gtable_0.3.6 BiocIO_1.21.0
## [61] tibble_3.3.0 pillar_1.11.1 htmltools_0.5.9
## [64] quantreg_6.1 R6_2.6.1 evaluate_1.0.5
## [67] lattice_0.22-7 Rsamtools_2.27.0 cigarillo_1.1.0
## [70] bslib_0.9.0 MatrixModels_0.5-4 coda_0.19-4.1
## [73] SparseArray_1.11.10 xfun_0.55 buildtools_1.0.0
## [76] pkgconfig_2.0.3
estiParamdmSingleplotGene
estiParamdmTwoGroups