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.251235 -0.52615638 0.39023560 0.42033359 -0.02599344
## ENSMUSG00000000003 1.574206 1.37512009 2.98043018 -1.90416699 -2.75332804
## ENSMUSG00000000028 1.296041 -0.11933866 0.32250371 0.01783067 -0.06873215
## ENSMUSG00000000037 1.019019 -5.34840936 15.08173981 -7.55292794 -2.19672655
## ENSMUSG00000000049 1.030708 -0.04910964 0.08233504 0.08448160 0.04374166
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.302492 14.105541 3.168315 1.768467
## ENSMUSG00000000003 23.126506 3.678297 5.065776 8.225515
## ENSMUSG00000000028 8.011839 7.335632 3.133549 2.474678
## ENSMUSG00000000037 7.911133 12.757499 6.926176 2.288500
## ENSMUSG00000000049 5.957227 9.363629 3.675403 1.251859
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.070602135 0.032255830 0.013094748 0.007337344
## ENSMUSG00000000028
## 0.007173050
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.231054 -0.3937353 0.2563001 0.31254017 0.08997116
## ENSMUSG00000000003 1.589262 1.7666945 2.0598294 -1.87113319 -2.28233559
## ENSMUSG00000000028 1.279408 -0.1063847 0.3053657 0.01241552 -0.05933072
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.157675 13.200676 3.481040 1.983486
## ENSMUSG00000000003 23.404223 3.164957 4.816241 8.137398
## ENSMUSG00000000028 7.697486 8.123531 3.225649 2.504152
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9062137 -1.281048 7.456582 -7.233937 0.9351495
## ENSMUSG00000000003 -0.7985744 -1.395990 4.233355 -2.179590 -0.5941231
## ENSMUSG00000000028 2.3322469 -2.168141 9.505404 -11.870873 4.7186224
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.318967 6.845212 3.175845 1.434011
## ENSMUSG00000000003 6.212816 9.234135 4.477713 2.792110
## ENSMUSG00000000028 11.495604 6.533785 3.152168 3.372169
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.063200015 0.031270699 0.022519373 0.010453134
## ENSMUSG00000000028
## 0.005222698
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.6.0 (2026-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 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=en_US.UTF-8
## [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.3 SingleCellExperiment_1.34.0
## [3] SummarizedExperiment_1.42.0 Biobase_2.72.0
## [5] GenomicRanges_1.64.0 Seqinfo_1.2.0
## [7] IRanges_2.46.0 S4Vectors_0.50.0
## [9] BiocGenerics_0.58.0 generics_0.1.4
## [11] MatrixGenerics_1.24.0 matrixStats_1.5.0
## [13] mist_1.4.0 BiocStyle_2.40.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.1 farver_2.1.2
## [4] Biostrings_2.80.0 S7_0.2.2 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.18 GenomicAlignments_1.48.0
## [10] XML_3.99-0.23 digest_0.6.39 lifecycle_1.0.5
## [13] survival_3.8-6 magrittr_2.0.5 compiler_4.6.0
## [16] rlang_1.2.0 sass_0.4.10 tools_4.6.0
## [19] yaml_2.3.12 rtracklayer_1.72.0 knitr_1.51
## [22] S4Arrays_1.12.0 labeling_0.4.3 curl_7.1.0
## [25] DelayedArray_0.38.0 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.46.0 withr_3.0.2 sys_3.4.3
## [31] grid_4.6.0 scales_1.4.0 MASS_7.3-65
## [34] mcmc_0.9-8 cli_3.6.6 mvtnorm_1.3-7
## [37] rmarkdown_2.31 crayon_1.5.3 httr_1.4.8
## [40] rjson_0.2.23 cachem_1.1.0 splines_4.6.0
## [43] parallel_4.6.0 BiocManager_1.30.27 XVector_0.52.0
## [46] restfulr_0.0.16 vctrs_0.7.3 Matrix_1.7-5
## [49] jsonlite_2.0.0 SparseM_1.84-2 carData_3.0-6
## [52] car_3.1-5 MCMCpack_1.7-1 Formula_1.2-5
## [55] maketools_1.3.2 jquerylib_0.1.4 glue_1.8.1
## [58] codetools_0.2-20 gtable_0.3.6 BiocIO_1.22.0
## [61] tibble_3.3.1 pillar_1.11.1 htmltools_0.5.9
## [64] quantreg_6.1 R6_2.6.1 evaluate_1.0.5
## [67] lattice_0.22-9 Rsamtools_2.28.0 cigarillo_1.2.0
## [70] bslib_0.10.0 MatrixModels_0.5-4 coda_0.19-4.1
## [73] SparseArray_1.12.0 xfun_0.57 buildtools_1.0.0
## [76] pkgconfig_2.0.3
estiParamdmSingleplotGene
estiParamdmTwoGroups