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.268368 -0.493016964  0.5097054  0.31225427 -0.11255822
## ENSMUSG00000000003 1.544941  1.542133694  2.8549180 -1.75754011 -2.86007066
## ENSMUSG00000000028 1.294364  0.009540156  0.1068241  0.01793735 -0.02237140
## ENSMUSG00000000037 1.053545 -4.628421848 12.9348016 -5.61607029 -2.73367975
## ENSMUSG00000000049 1.010878 -0.124899476  0.1274451  0.11493762  0.06949836
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.505557 14.999090 3.594427 1.770791
## ENSMUSG00000000003 24.719346  3.384022 7.509212 9.111115
## ENSMUSG00000000028  7.464531  7.644023 3.635891 2.303272
## ENSMUSG00000000037  8.511641 12.921891 6.980450 2.221980
## ENSMUSG00000000049  5.652139  9.128325 2.901747 1.317275dmSingle# 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.062606403        0.032962403        0.011378271        0.008744585 
## ENSMUSG00000000028 
##        0.004827250plotGene# 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.264373 -0.50393635 0.4464748  0.31479556 -0.03650538
## ENSMUSG00000000003 1.507949  1.62616534 3.1017785 -1.83044990 -3.17190169
## ENSMUSG00000000028 1.288418  0.01005007 0.1025316  0.01597944 -0.02290974
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.837911 13.247077 3.809745 1.742368
## ENSMUSG00000000003 24.455615  4.481056 7.310687 9.337942
## ENSMUSG00000000028  7.438908  7.608484 3.761568 2.261506head(rowData(Dat_sce_g2)$mist_pars, n = 3)##                        Beta_0      Beta_1    Beta_2      Beta_3     Beta_4
## ENSMUSG00000000001  1.9203339  -0.6839444  4.537461   -2.727702 -1.2855591
## ENSMUSG00000000003 -0.8179274  -1.3738261  3.783131   -1.392957 -0.9510838
## ENSMUSG00000000028  2.2684264 -17.2559992 82.462057 -117.792832 52.8611262
##                    Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.679195 6.387089 3.278315 1.441077
## ENSMUSG00000000003 7.597635 9.543038 4.561997 2.940527
## ENSMUSG00000000028 8.620777 6.288073 3.675560 3.304914dmTwoGroups# 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)## ENSMUSG00000000028 ENSMUSG00000000037 ENSMUSG00000000003 ENSMUSG00000000001 
##         0.04972029         0.03736648         0.03265349         0.02132975 
## ENSMUSG00000000049 
##         0.01246003mist 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.32.0
##  [3] SummarizedExperiment_1.40.0 Biobase_2.70.0             
##  [5] GenomicRanges_1.62.0        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.0                  BiocStyle_2.38.0           
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1         dplyr_1.1.4              farver_2.1.2            
##  [4] Biostrings_2.78.0        S7_0.2.0                 bitops_1.0-9            
##  [7] fastmap_1.2.0            RCurl_1.98-1.17          GenomicAlignments_1.46.0
## [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.70.0       knitr_1.50              
## [22] labeling_0.4.3           S4Arrays_1.10.0          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              grid_4.5.1              
## [31] scales_1.4.0             MASS_7.3-65              mcmc_0.9-8              
## [34] tinytex_0.57             dichromat_2.0-0.1        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.50.0           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.20.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.26.0         cigarillo_1.0.0          bslib_0.9.0             
## [73] MatrixModels_0.5-4       Rcpp_1.1.0               coda_0.19-4.1           
## [76] SparseArray_1.10.0       xfun_0.53                pkgconfig_2.0.3