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.258209 -0.46543563 0.43272467 0.26916929 -0.01760002
## ENSMUSG00000000003 1.629619 -0.76478593 13.55127934 -15.90449821 2.81945939
## ENSMUSG00000000028 1.280700 -0.03056072 0.13948963 0.04315599 -0.01156884
## ENSMUSG00000000037 1.017325 -6.39317372 17.11613387 -8.16697273 -2.55774297
## ENSMUSG00000000049 1.022287 -0.08960261 0.07587963 0.09711572 0.08252785
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.540654 13.272853 3.412429 1.639169
## ENSMUSG00000000003 24.461810 6.221647 5.775931 9.284335
## ENSMUSG00000000028 7.323966 8.577230 3.218084 2.247009
## ENSMUSG00000000037 8.036366 14.776163 6.874131 2.274733
## ENSMUSG00000000049 5.913013 8.692354 3.016373 1.189801
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.081878212 0.038116510 0.011145959 0.007516655
## ENSMUSG00000000028
## 0.005974192
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.236942 -0.47761818 0.3446414 0.30917035 0.07552722
## ENSMUSG00000000003 1.630755 0.63933514 8.6628248 -9.83438931 0.14926298
## ENSMUSG00000000028 1.283961 -0.02715222 0.1445697 0.03822086 -0.02297181
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.228473 14.481626 3.108797 1.777443
## ENSMUSG00000000003 24.919380 5.496385 6.038831 8.671463
## ENSMUSG00000000028 7.446209 7.884764 3.187679 2.250759
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9129718 -0.7589824 5.081785 -3.6806882 -0.8044568
## ENSMUSG00000000003 -0.8484761 -2.7323846 8.353467 -5.1310115 -0.4862201
## ENSMUSG00000000028 2.3197197 -0.1517510 1.018528 -0.3452947 -0.3871734
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.538857 7.295441 3.137310 1.342026
## ENSMUSG00000000003 6.094796 9.270331 5.192968 3.022225
## ENSMUSG00000000028 11.514815 5.091852 3.127074 3.118985
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.056242170 0.038742724 0.022764549 0.009721450
## ENSMUSG00000000028
## 0.002546106
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 Under development (unstable) (2025-12-22 r89219)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.23-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.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.6.0
## [16] rlang_1.1.7 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.0 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.58 dichromat_2.0-0.1 cli_3.6.5
## [37] mvtnorm_1.3-3 rmarkdown_2.30 crayon_1.5.3
## [40] otel_0.2.0 httr_1.4.7 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.6.5 Matrix_1.7-4 jsonlite_2.0.0
## [52] SparseM_1.84-2 carData_3.0-5 bookdown_0.46
## [55] car_3.1-3 MCMCpack_1.7-1 Formula_1.2-5
## [58] magick_2.9.0 jquerylib_0.1.4 glue_1.8.0
## [61] codetools_0.2-20 gtable_0.3.6 BiocIO_1.21.0
## [64] tibble_3.3.0 pillar_1.11.1 htmltools_0.5.9
## [67] quantreg_6.1 R6_2.6.1 evaluate_1.0.5
## [70] lattice_0.22-7 Rsamtools_2.27.0 cigarillo_1.1.0
## [73] bslib_0.9.0 MatrixModels_0.5-4 Rcpp_1.1.0.8.2
## [76] coda_0.19-4.1 SparseArray_1.11.10 xfun_0.55
## [79] pkgconfig_2.0.3