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.247744 -0.55621185 0.4408920 0.42459559 -0.06974521
## ENSMUSG00000000003 1.497222 0.97550437 4.5068271 -2.86596635 -2.87579933
## ENSMUSG00000000028 1.294696 -0.01566068 0.1184678 0.03410627 -0.01652968
## ENSMUSG00000000037 1.015887 -5.17374634 15.3327487 -8.05074296 -2.17326840
## ENSMUSG00000000049 1.022185 -0.14994702 0.1258588 0.10717922 0.09677429
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.558307 14.186320 3.163354 1.714499
## ENSMUSG00000000003 24.172354 8.325426 5.890633 9.361806
## ENSMUSG00000000028 7.742754 7.597192 3.477485 2.310449
## ENSMUSG00000000037 8.399742 13.263543 8.821118 2.381865
## ENSMUSG00000000049 6.051868 8.289543 3.177404 1.287998
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.069857734 0.036176230 0.012853305 0.008392260
## ENSMUSG00000000028
## 0.005127856
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.238881 -0.812006739 0.58507949 0.56058502 -0.03838835
## ENSMUSG00000000003 1.580430 1.471180653 2.62491775 -1.97977662 -2.39063563
## ENSMUSG00000000028 1.294058 -0.002085751 0.09029263 0.02822138 -0.01091758
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.559921 12.970724 3.020887 1.767316
## ENSMUSG00000000003 25.384368 4.568840 5.218424 8.873856
## ENSMUSG00000000028 8.130446 7.802975 3.136585 2.333930
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4 Sigma2_1
## ENSMUSG00000000001 1.9072704 -0.559547 3.955544 -2.15724 -1.384450 5.604358
## ENSMUSG00000000003 -0.8314969 -1.753746 4.948880 -1.62998 -1.563977 7.304354
## ENSMUSG00000000028 2.3172556 -10.503731 49.781084 -69.34902 30.284748 9.928704
## Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.167275 3.004716 1.374840
## ENSMUSG00000000003 8.948249 4.122950 3.341158
## ENSMUSG00000000028 5.202400 4.136011 3.751942
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 ENSMUSG00000000028 ENSMUSG00000000001
## 0.054309643 0.031626641 0.028384590 0.026200395
## ENSMUSG00000000049
## 0.009942994
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.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.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] labeling_0.4.3 S4Arrays_1.12.0 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 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.52.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.22.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.28.0 cigarillo_1.2.0
## [73] bslib_0.10.0 MatrixModels_0.5-4 Rcpp_1.1.1-1.1
## [76] coda_0.19-4.1 SparseArray_1.12.0 xfun_0.57
## [79] pkgconfig_2.0.3