Contents

0.1 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.

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

0.3 Example Workflow for Single-Group Analysis

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

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

2 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.254570 -0.66428126  0.5890856  0.32220381  0.02693085
## ENSMUSG00000000003 1.587354  1.81821615  1.7484529 -1.42884779 -2.40340675
## ENSMUSG00000000028 1.298808 -0.01870284  0.1195747  0.03277182 -0.01212524
## ENSMUSG00000000037 1.016930 -4.89567929 13.2765761 -5.84533783 -2.53502124
## ENSMUSG00000000049 1.018608 -0.08200303  0.1084319  0.07221070  0.05115999
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.534338 13.642131 3.734521 1.963797
## ENSMUSG00000000003 25.883281  3.656273 6.115603 8.737515
## ENSMUSG00000000028  7.466891  7.260971 3.797514 2.277551
## ENSMUSG00000000037  8.788247 13.565056 6.795162 2.303394
## ENSMUSG00000000049  5.943850  8.263436 3.067052 1.184817

3 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.065803405        0.030241906        0.014715770        0.006987605 
## ENSMUSG00000000028 
##        0.005141305

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

4.1 Example Workflow for Two-Group Analysis

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

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

6 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.242179 -0.678591685 0.76905422  0.30193115 -0.104913967
## ENSMUSG00000000003 1.565390  1.753905814 2.61003512 -1.89259634 -2.751706140
## ENSMUSG00000000028 1.297740 -0.008799848 0.08440254  0.02836873 -0.001974097
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.484694 13.926588 3.432080 2.237804
## ENSMUSG00000000003 24.393702  4.177242 6.405818 9.494624
## ENSMUSG00000000028  7.957134  7.774183 2.993894 2.207018
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
##                        Beta_0    Beta_1    Beta_2     Beta_3    Beta_4
## ENSMUSG00000000001  1.9017308 -3.696403 20.312050 -26.381107 9.6227984
## ENSMUSG00000000003 -0.8381918 -2.396281  7.951101  -6.373335 0.8607604
## ENSMUSG00000000028  2.3491996 -2.456967 10.653108 -11.927032 3.8272402
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.539401  6.406628 4.075735 1.336894
## ENSMUSG00000000003  6.772438 10.604057 4.962369 3.162399
## ENSMUSG00000000028 11.075112  4.658120 3.647120 3.547162

7 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.062070012        0.033195085        0.031254865        0.010393782 
## ENSMUSG00000000028 
##        0.009343981

7.1 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 Under development (unstable) (2025-10-20 r88955)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 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.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.6.0          
## [16] rlang_1.1.6              sass_0.4.10              tools_4.6.0             
## [19] yaml_2.3.10              rtracklayer_1.71.0       knitr_1.50              
## [22] labeling_0.4.3           S4Arrays_1.11.0          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              grid_4.6.0              
## [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.6.0            parallel_4.6.0           BiocManager_1.30.26     
## [46] XVector_0.51.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.21.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.27.0         cigarillo_1.1.0          bslib_0.9.0             
## [73] MatrixModels_0.5-4       Rcpp_1.1.0               coda_0.19-4.1           
## [76] SparseArray_1.11.1       xfun_0.54                pkgconfig_2.0.3