estiParamdmSingleplotGeneestiParamdmTwoGroupsmist (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:
library(mist)
vignette("mist")
In this section, we will estimate parameters and perform differential methylation analysis using single-group data.
Here we load the example data from GSE121708.
library(mist)
library(SingleCellExperiment)
# Load sample scDNAm data
Dat_sce <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist"))
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.241749 -0.76800331 0.7346722 0.55171213 -0.28001415
## ENSMUSG00000000003 1.633451 1.60283171 3.1555495 -2.57154929 -2.48211186
## ENSMUSG00000000028 1.288420 -0.02085989 0.1092266 0.02767581 -0.01446369
## ENSMUSG00000000037 1.043325 -2.86032345 8.4659260 -4.04928792 -1.49510886
## ENSMUSG00000000049 1.013552 -0.15644863 0.1489367 0.13059243 0.07849163
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.683649 12.376846 3.041475 1.685727
## ENSMUSG00000000003 24.588812 3.150670 6.452057 9.132812
## ENSMUSG00000000028 7.544277 8.062820 3.061076 2.122944
## ENSMUSG00000000037 8.502920 12.973362 7.752955 2.409681
## ENSMUSG00000000049 5.637839 9.709056 2.723889 1.344920
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.043483439 0.031599806 0.015057152 0.009486068
## ENSMUSG00000000028
## 0.004390350
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.
# Load two-group scDNAm data
Dat_sce_g1 <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist"))
Dat_sce_g2 <- readRDS(system.file("extdata", "group2_sampleData_sce.rds", package = "mist"))
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.240868 -0.66883273 0.6046788 0.44378928 -0.155298009
## ENSMUSG00000000003 1.584353 1.96038547 3.4332142 -3.12988211 -2.675816290
## ENSMUSG00000000028 1.287033 -0.02661283 0.1222689 0.03102503 0.005224782
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.821696 12.615866 3.086289 1.591238
## ENSMUSG00000000003 25.020793 4.600922 6.569482 8.729799
## ENSMUSG00000000028 7.734164 6.787281 3.036936 2.490455
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.8996838 -1.37778290 8.3190057 -8.39606794 1.2827175
## ENSMUSG00000000003 -0.8120776 -0.24898106 0.8628573 0.04759776 -0.5586017
## ENSMUSG00000000028 2.3225449 -0.08148821 0.8862270 -0.14389036 -0.5379907
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.608289 5.447918 3.424949 1.283065
## ENSMUSG00000000003 6.867410 11.097983 4.233418 2.829036
## ENSMUSG00000000028 11.227907 6.323533 3.259783 3.084103
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.052810531 0.032454854 0.025570388 0.009248057
## ENSMUSG00000000028
## 0.003154980
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 RC (2026-04-17 r89917)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.24-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## 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
##
## time zone: America/New_York
## 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.33.2
## [3] SummarizedExperiment_1.41.1 Biobase_2.71.0
## [5] GenomicRanges_1.63.2 Seqinfo_1.1.0
## [7] IRanges_2.45.0 S4Vectors_0.49.2
## [9] BiocGenerics_0.57.1 generics_0.1.4
## [11] MatrixGenerics_1.23.0 matrixStats_1.5.0
## [13] mist_1.3.3 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.1 farver_2.1.2
## [4] Biostrings_2.79.5 S7_0.2.1-1 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.18 GenomicAlignments_1.47.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.71.3 knitr_1.51
## [22] labeling_0.4.3 S4Arrays_1.11.1 curl_7.0.0
## [25] DelayedArray_0.37.1 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.45.0 withr_3.0.2 grid_4.6.0
## [31] scales_1.4.0 MASS_7.3-65 mcmc_0.9-8
## [34] tinytex_0.59 dichromat_2.0-0.1 cli_3.6.6
## [37] mvtnorm_1.3-7 rmarkdown_2.31 crayon_1.5.3
## [40] otel_0.2.0 httr_1.4.8 rjson_0.2.23
## [43] cachem_1.1.0 splines_4.6.0 parallel_4.6.0
## [46] BiocManager_1.30.27 XVector_0.51.0 restfulr_0.0.16
## [49] vctrs_0.7.3 Matrix_1.7-5 jsonlite_2.0.0
## [52] SparseM_1.84-2 carData_3.0-6 bookdown_0.46
## [55] car_3.1-5 MCMCpack_1.7-1 Formula_1.2-5
## [58] magick_2.9.1 jquerylib_0.1.4 glue_1.8.1
## [61] codetools_0.2-20 gtable_0.3.6 BiocIO_1.21.0
## [64] tibble_3.3.1 pillar_1.11.1 htmltools_0.5.9
## [67] quantreg_6.1 R6_2.6.1 evaluate_1.0.5
## [70] lattice_0.22-9 Rsamtools_2.27.2 cigarillo_1.1.0
## [73] bslib_0.10.0 MatrixModels_0.5-4 Rcpp_1.1.1-1
## [76] coda_0.19-4.1 SparseArray_1.11.13 xfun_0.57
## [79] pkgconfig_2.0.3