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.266768 -0.54051250  0.4633536  0.34662497 -0.024151074
## ENSMUSG00000000003 1.603219  1.84966880  2.6015515 -1.98519867 -2.744601583
## ENSMUSG00000000028 1.305001 -0.02077755  0.1132451  0.07266097 -0.003556298
## ENSMUSG00000000037 1.035405 -4.58219283 11.8943629 -3.75262040 -3.573892974
## ENSMUSG00000000049 1.023193 -0.09989114  0.1023212  0.10308438  0.061061125
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.650841 15.371655 3.760407 1.878646
## ENSMUSG00000000003 25.745003  2.723271 8.249524 9.024874
## ENSMUSG00000000028  7.690074  7.232867 3.163674 2.619686
## ENSMUSG00000000037  8.695636 12.634229 6.901861 2.271397
## ENSMUSG00000000049  5.841854  9.098041 2.980254 1.181453

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.064027300        0.032310814        0.012637190        0.007706739 
## ENSMUSG00000000028 
##        0.006554595

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.270070 -0.433327056 0.36464672  0.32180744 -0.03780113
## ENSMUSG00000000003 1.568221  2.175535303 2.92579092 -2.58183609 -2.89633591
## ENSMUSG00000000028 1.312524  0.001362241 0.09942835  0.04455538 -0.02964313
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  6.043138 14.410455 3.814026 1.835891
## ENSMUSG00000000003 26.192101  4.216544 8.074808 8.750057
## ENSMUSG00000000028  7.861299  6.729265 3.780578 2.204434
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
##                       Beta_0     Beta_1    Beta_2      Beta_3     Beta_4
## ENSMUSG00000000001  1.917650 -0.7529896  4.997538  -3.1101008 -1.3125902
## ENSMUSG00000000003 -0.851246 -0.8168955  2.508121  -0.9595917 -0.6373476
## ENSMUSG00000000028  2.343069 -3.2832118 14.665955 -17.4207489  6.1187961
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  6.043115  5.779985 4.274549 1.360212
## ENSMUSG00000000003  6.891536 10.126726 4.484020 3.437034
## ENSMUSG00000000028 10.660622  6.318389 3.605409 3.193554

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 ENSMUSG00000000028 
##         0.04351797         0.03267508         0.02168978         0.01123635 
## ENSMUSG00000000049 
##         0.01081597

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.5.1 (2025-06-13)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 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=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: 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.0               SingleCellExperiment_1.33.0
##  [3] SummarizedExperiment_1.41.0 Biobase_2.71.0             
##  [5] GenomicRanges_1.63.0        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.0                  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.1        S7_0.2.0                 bitops_1.0-9            
##  [7] fastmap_1.2.0            RCurl_1.98-1.17          GenomicAlignments_1.47.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.69.1       knitr_1.50              
## [22] S4Arrays_1.11.0          labeling_0.4.3           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              sys_3.4.3               
## [31] grid_4.5.1               scales_1.4.0             MASS_7.3-65             
## [34] mcmc_0.9-8               cli_3.6.5                mvtnorm_1.3-3           
## [37] rmarkdown_2.30           crayon_1.5.3             httr_1.4.7              
## [40] rjson_0.2.23             cachem_1.1.0             splines_4.5.1           
## [43] parallel_4.5.1           BiocManager_1.30.26      XVector_0.51.0          
## [46] restfulr_0.0.16          vctrs_0.6.5              Matrix_1.7-4            
## [49] jsonlite_2.0.0           SparseM_1.84-2           carData_3.0-5           
## [52] car_3.1-3                MCMCpack_1.7-1           Formula_1.2-5           
## [55] maketools_1.3.2          jquerylib_0.1.4          glue_1.8.0              
## [58] codetools_0.2-20         gtable_0.3.6             BiocIO_1.21.0           
## [61] tibble_3.3.0             pillar_1.11.1            htmltools_0.5.8.1       
## [64] quantreg_6.1             R6_2.6.1                 evaluate_1.0.5          
## [67] lattice_0.22-7           Rsamtools_2.27.0         cigarillo_1.1.0         
## [70] bslib_0.9.0              MatrixModels_0.5-4       coda_0.19-4.1           
## [73] SparseArray_1.11.1       xfun_0.54                buildtools_1.0.0        
## [76] pkgconfig_2.0.3