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.2352078 -0.59382846 0.54097541 0.40644417 -0.06557406
## ENSMUSG00000000003 1.5098697 1.36654560 2.55242898 -1.16926780 -2.93567734
## ENSMUSG00000000028 1.2746883 -0.02308856 0.06186967 0.07037421 0.03212724
## ENSMUSG00000000037 0.9965868 -3.71335752 9.99483896 -2.59030542 -3.69039196
## ENSMUSG00000000049 1.0183113 -0.08971393 0.10528193 0.08348824 0.05579343
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.079774 14.067479 3.406856 2.081700
## ENSMUSG00000000003 24.093041 5.015443 7.556076 8.732858
## ENSMUSG00000000028 8.071954 8.360294 2.741420 2.310786
## ENSMUSG00000000037 9.001149 14.016126 7.301419 2.956061
## ENSMUSG00000000049 6.070573 7.941924 2.889617 1.141733
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.058575139 0.032309156 0.014733458 0.007207325
## ENSMUSG00000000028
## 0.005618742
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.242212 -0.73314285 0.69080375 0.43198701 -0.09798260
## ENSMUSG00000000003 1.446763 1.17919143 2.85267257 -1.29330612 -2.87003423
## ENSMUSG00000000028 1.277998 -0.01593982 0.06350528 0.06394311 0.04090137
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.220353 13.721410 3.501995 1.949820
## ENSMUSG00000000003 23.381762 9.227266 7.531978 8.890925
## ENSMUSG00000000028 7.938390 6.543336 2.965063 2.542555
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9033600 -0.91754081 4.895633 -2.78366005 -1.3183120
## ENSMUSG00000000003 -0.8433011 -0.73877692 2.453010 -1.04709872 -0.5968822
## ENSMUSG00000000028 2.3488960 0.03154822 0.458176 -0.03472785 -0.3497799
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.536664 5.320669 3.321906 1.627234
## ENSMUSG00000000003 6.516494 9.589168 4.961824 2.845080
## ENSMUSG00000000028 12.270334 5.645936 3.831583 2.894099
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.055433832 0.025234654 0.023444335 0.012106817
## ENSMUSG00000000028
## 0.003561124
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.2 SingleCellExperiment_1.32.0
## [3] SummarizedExperiment_1.40.0 Biobase_2.70.0
## [5] GenomicRanges_1.62.1 Seqinfo_1.0.0
## [7] IRanges_2.44.0 S4Vectors_0.48.0
## [9] BiocGenerics_0.56.0 generics_0.1.4
## [11] MatrixGenerics_1.22.0 matrixStats_1.5.0
## [13] mist_1.2.3 BiocStyle_2.38.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.0 farver_2.1.2
## [4] Biostrings_2.78.0 S7_0.2.1 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.17 GenomicAlignments_1.46.0
## [10] XML_3.99-0.20 digest_0.6.39 lifecycle_1.0.5
## [13] survival_3.8-6 magrittr_2.0.4 compiler_4.5.2
## [16] rlang_1.1.7 sass_0.4.10 tools_4.5.2
## [19] yaml_2.3.12 rtracklayer_1.70.1 knitr_1.51
## [22] S4Arrays_1.10.1 labeling_0.4.3 curl_7.0.0
## [25] DelayedArray_0.36.0 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.44.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.50.0
## [46] restfulr_0.0.16 vctrs_0.7.1 Matrix_1.7-4
## [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.0
## [58] codetools_0.2-20 gtable_0.3.6 BiocIO_1.20.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-7 Rsamtools_2.26.0 cigarillo_1.0.0
## [70] bslib_0.10.0 MatrixModels_0.5-4 coda_0.19-4.1
## [73] SparseArray_1.10.8 xfun_0.56 buildtools_1.0.0
## [76] pkgconfig_2.0.3
estiParamdmSingleplotGene
estiParamdmTwoGroups