mist:methylation inference for single-cell along trajectory

Introduction

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.

Installation

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")

Example Workflow for Single-Group Analysis

In this section, we will estimate parameters and perform differential methylation analysis using single-group data.

Step 1: Load Example 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"))

Step 2: Estimate Parameters Using 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.270574 -0.485643061  0.36447120  0.27602409  0.07925594
## ENSMUSG00000000003 1.593984  1.939530623  2.40750973 -2.34737563 -2.32462268
## ENSMUSG00000000028 1.258473 -0.003330174  0.05312787  0.04032000  0.01704900
## ENSMUSG00000000037 1.027970 -4.052447595 10.67916017 -3.74948496 -2.90765206
## ENSMUSG00000000049 1.023655 -0.050343965  0.06407348  0.06382444  0.05492184
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.901728 15.237998 4.000153 1.786390
## ENSMUSG00000000003 25.609861  3.278955 6.290538 8.483097
## ENSMUSG00000000028  8.242223  6.456884 2.788160 2.314968
## ENSMUSG00000000037  8.826890 13.901417 6.467071 2.009865
## ENSMUSG00000000049  6.176580  8.325357 3.089389 1.112546

Step 3: Perform Differential Methylation Analysis Using 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.055823432        0.032403815        0.011713180        0.006007808 
## ENSMUSG00000000028 
##        0.004367364

Step 4: Perform Differential Methylation Analysis Using 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")

Example Workflow for Two-Group Analysis

In this section, we will estimate parameters and perform DM analysis using data from two phenotypic groups.

Step 1: Load Two-Group Data

# 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"))

Step 2: Estimate Parameters Using 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.265787 -0.42085936 0.39426888  0.22568896  0.01314925
## ENSMUSG00000000003 1.500251  2.08369483 1.80585647 -1.88777307 -2.23973495
## ENSMUSG00000000028 1.251440 -0.01322191 0.06060946  0.05428661  0.02969226
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.922087 14.133139 4.328806 1.792165
## ENSMUSG00000000003 24.692213  4.721189 7.236784 8.766995
## ENSMUSG00000000028  7.841481  7.374258 3.148056 2.402697
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
##                        Beta_0     Beta_1    Beta_2     Beta_3     Beta_4
## ENSMUSG00000000001  1.9019991 -4.4327459 22.720611 -28.693149 10.2740699
## ENSMUSG00000000003 -0.8380432 -2.1429040  6.844693  -5.030565  0.3803865
## ENSMUSG00000000028  2.3244878 -0.5472775  2.914222  -2.735423  0.5029152
##                     Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.168388 7.461159 3.944114 1.589582
## ENSMUSG00000000003  6.724073 9.878083 4.243961 3.175602
## ENSMUSG00000000028 11.679581 4.777432 4.931178 3.333007

Step 3: Perform Differential Methylation Analysis for Two-Group Comparison Using 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.046873628        0.034437051        0.031160088        0.009233123 
## ENSMUSG00000000028 
##        0.004434071

Conclusion

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.

Session info

## R version 4.6.0 (2026-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [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: Etc/UTC
## 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.35.1
##  [3] SummarizedExperiment_1.43.0 Biobase_2.73.1             
##  [5] GenomicRanges_1.65.0        Seqinfo_1.3.0              
##  [7] IRanges_2.47.2              S4Vectors_0.51.3           
##  [9] BiocGenerics_0.59.6         generics_0.1.4             
## [11] MatrixGenerics_1.25.0       matrixStats_1.5.0          
## [13] mist_1.5.0                  BiocStyle_2.41.0           
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1         dplyr_1.2.1              farver_2.1.2            
##  [4] Biostrings_2.81.2        S7_0.2.2                 bitops_1.0-9            
##  [7] fastmap_1.2.0            RCurl_1.98-1.18          GenomicAlignments_1.49.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.73.0       knitr_1.51              
## [22] labeling_0.4.3           S4Arrays_1.13.0          curl_7.1.0              
## [25] DelayedArray_0.39.3      RColorBrewer_1.1-3       abind_1.4-8             
## [28] BiocParallel_1.47.0      withr_3.0.2              sys_3.4.3               
## [31] grid_4.6.0               scales_1.4.0             MASS_7.3-65             
## [34] mcmc_0.9-8               cli_3.6.6                mvtnorm_1.3-7           
## [37] rmarkdown_2.31           crayon_1.5.3             httr_1.4.8              
## [40] rjson_0.2.23             BiocBaseUtils_1.15.1     cachem_1.1.0            
## [43] splines_4.6.0            parallel_4.6.0           BiocManager_1.30.27     
## [46] XVector_0.53.0           restfulr_0.0.16          vctrs_0.7.3             
## [49] Matrix_1.7-5             jsonlite_2.0.0           SparseM_1.84-2          
## [52] carData_3.0-6            car_3.1-5                MCMCpack_1.7-1          
## [55] Formula_1.2-5            maketools_1.3.2          jquerylib_0.1.4         
## [58] glue_1.8.1               codetools_0.2-20         gtable_0.3.6            
## [61] BiocIO_1.23.3            tibble_3.3.1             pillar_1.11.1           
## [64] htmltools_0.5.9          quantreg_6.1             R6_2.6.1                
## [67] evaluate_1.0.5           lattice_0.22-9           Rsamtools_2.29.0        
## [70] cigarillo_1.3.0          bslib_0.11.0             MatrixModels_0.5-4      
## [73] coda_0.19-4.1            SparseArray_1.13.2       xfun_0.57               
## [76] buildtools_1.0.0         pkgconfig_2.0.3