estiParam
dmSingle
plotGene
estiParam
dmTwoGroups
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:
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.256282 -0.711084488 0.73955969 0.32579364 -0.08136019
## ENSMUSG00000000003 1.602679 1.876446907 3.18206926 -2.40104786 -3.02899059
## ENSMUSG00000000028 1.308807 0.003071711 0.12859765 0.04399916 -0.04538110
## ENSMUSG00000000037 1.030676 -4.038589590 11.24514796 -5.27373857 -1.90219441
## ENSMUSG00000000049 1.026477 -0.075422250 0.08410717 0.06897879 0.07256543
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.707243 14.195256 3.135102 2.013360
## ENSMUSG00000000003 25.120886 6.268850 6.220360 9.360707
## ENSMUSG00000000028 7.561472 7.855669 3.558540 2.520030
## ENSMUSG00000000037 8.446140 14.810222 7.001826 2.231406
## ENSMUSG00000000049 6.474145 8.623758 3.203047 1.193284
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.055354394 0.034567822 0.015089743 0.006825270
## ENSMUSG00000000028
## 0.005636038
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.269016 -0.596814300 0.5568271 0.28612544 -0.006004751
## ENSMUSG00000000003 1.671126 1.762495221 3.2459344 -2.17861495 -3.229399714
## ENSMUSG00000000028 1.298523 0.002442056 0.1141796 0.04411878 -0.036442865
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.143100 13.032641 3.648959 1.852924
## ENSMUSG00000000003 25.841130 3.765555 5.982952 9.502159
## ENSMUSG00000000028 7.862449 8.652110 3.276898 2.623437
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9309537 -0.7319410 5.4615493 -4.3930030 -0.5124358
## ENSMUSG00000000003 -0.7968563 -0.6065772 1.9868299 -0.9432772 -0.3801065
## ENSMUSG00000000028 2.3464725 -0.1260123 0.8124242 -0.2217749 -0.2669067
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.923749 7.179324 3.382531 1.498424
## ENSMUSG00000000003 8.034342 10.038632 4.136010 2.782803
## ENSMUSG00000000028 11.453475 5.426425 3.557650 3.502671
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.053531553 0.028425198 0.023637467 0.009945321
## ENSMUSG00000000028
## 0.001298918
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.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-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.0 SingleCellExperiment_1.31.1
## [3] SummarizedExperiment_1.39.2 Biobase_2.69.1
## [5] GenomicRanges_1.61.5 Seqinfo_0.99.2
## [7] IRanges_2.43.5 S4Vectors_0.47.4
## [9] BiocGenerics_0.55.1 generics_0.1.4
## [11] MatrixGenerics_1.21.0 matrixStats_1.5.0
## [13] mist_1.1.0 BiocStyle_2.37.1
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.1.4 farver_2.1.2
## [4] Biostrings_2.77.2 S7_0.2.0 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.17 GenomicAlignments_1.45.4
## [10] XML_3.99-0.19 digest_0.6.37 lifecycle_1.0.4
## [13] survival_3.8-3 magrittr_2.0.4 compiler_4.5.1
## [16] rlang_1.1.6 sass_0.4.10 tools_4.5.1
## [19] yaml_2.3.10 rtracklayer_1.69.1 knitr_1.50
## [22] S4Arrays_1.9.1 labeling_0.4.3 curl_7.0.0
## [25] DelayedArray_0.35.3 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.43.4 withr_3.0.2 grid_4.5.1
## [31] scales_1.4.0 MASS_7.3-65 mcmc_0.9-8
## [34] dichromat_2.0-0.1 tinytex_0.57 cli_3.6.5
## [37] mvtnorm_1.3-3 rmarkdown_2.30 crayon_1.5.3
## [40] httr_1.4.7 rjson_0.2.23 cachem_1.1.0
## [43] splines_4.5.1 parallel_4.5.1 BiocManager_1.30.26
## [46] XVector_0.49.1 restfulr_0.0.16 vctrs_0.6.5
## [49] Matrix_1.7-4 jsonlite_2.0.0 SparseM_1.84-2
## [52] carData_3.0-5 bookdown_0.45 car_3.1-3
## [55] MCMCpack_1.7-1 Formula_1.2-5 magick_2.9.0
## [58] jquerylib_0.1.4 glue_1.8.0 codetools_0.2-20
## [61] gtable_0.3.6 BiocIO_1.19.0 tibble_3.3.0
## [64] pillar_1.11.1 htmltools_0.5.8.1 quantreg_6.1
## [67] R6_2.6.1 evaluate_1.0.5 lattice_0.22-7
## [70] Rsamtools_2.25.3 bslib_0.9.0 MatrixModels_0.5-4
## [73] Rcpp_1.1.0 coda_0.19-4.1 SparseArray_1.9.1
## [76] xfun_0.53 pkgconfig_2.0.3