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.241292 -0.62490620  0.51857314  0.33853903  0.050143708
## ENSMUSG00000000003 1.583575  1.70737077  3.25832457 -2.57507155 -2.684669917
## ENSMUSG00000000028 1.286394 -0.01733794  0.10311323  0.03097087 -0.003909041
## ENSMUSG00000000037 1.033816 -3.73638198 10.30363194 -3.27337060 -3.336715718
## ENSMUSG00000000049 1.032356 -0.06913694  0.07363689  0.07877664  0.091803005
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.602400 15.452745 3.124282 1.948904
## ENSMUSG00000000003 25.776135  3.708126 7.970795 8.993499
## ENSMUSG00000000028  8.401626  6.767990 2.853196 2.156276
## ENSMUSG00000000037  8.955372 14.419091 7.011103 2.407441
## ENSMUSG00000000049  6.294127  9.362584 3.411021 1.373972

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.055857718        0.033659140        0.014710837        0.007765350 
## ENSMUSG00000000028 
##        0.004759708

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.262988 -0.56669940 0.5616504  0.26777072 -0.02316883
## ENSMUSG00000000003 1.550099  1.16247844 3.9640714 -2.31759563 -3.07198215
## ENSMUSG00000000028 1.276870 -0.04839283 0.1346460  0.04474141  0.01704771
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.811874 13.983816 3.735228 1.784201
## ENSMUSG00000000003 24.519159  7.826680 6.209487 9.607665
## ENSMUSG00000000028  7.814366  7.737493 2.905221 2.324463
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
##                        Beta_0    Beta_1    Beta_2      Beta_3    Beta_4
## ENSMUSG00000000001  1.9258027 -1.144325  6.715761  -6.5918728  0.942872
## ENSMUSG00000000003 -0.8507717 -1.008402  2.702022  -0.5298433 -1.114373
## ENSMUSG00000000028  2.3100747 -7.982103 35.501428 -46.4747577 19.084240
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.614912  6.423551 2.962210 1.908368
## ENSMUSG00000000003  6.639421 10.882583 4.628323 2.977315
## ENSMUSG00000000028 10.258500  6.289765 3.882762 3.080582

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 ENSMUSG00000000028 ENSMUSG00000000001 
##         0.04818409         0.02892537         0.02253438         0.02091419 
## ENSMUSG00000000049 
##         0.01051324

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 version 4.5.0 RC (2025-04-04 r88126)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.21-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_3.5.2               SingleCellExperiment_1.30.0
##  [3] SummarizedExperiment_1.38.0 Biobase_2.68.0             
##  [5] GenomicRanges_1.60.0        GenomeInfoDb_1.44.0        
##  [7] IRanges_2.42.0              S4Vectors_0.46.0           
##  [9] BiocGenerics_0.54.0         generics_0.1.3             
## [11] MatrixGenerics_1.20.0       matrixStats_1.5.0          
## [13] mist_1.0.0                  BiocStyle_2.36.0           
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1         farver_2.1.2             dplyr_1.1.4             
##  [4] Biostrings_2.76.0        bitops_1.0-9             fastmap_1.2.0           
##  [7] RCurl_1.98-1.17          GenomicAlignments_1.44.0 XML_3.99-0.18           
## [10] digest_0.6.37            lifecycle_1.0.4          survival_3.8-3          
## [13] magrittr_2.0.3           compiler_4.5.0           rlang_1.1.6             
## [16] sass_0.4.10              tools_4.5.0              yaml_2.3.10             
## [19] rtracklayer_1.68.0       knitr_1.50               labeling_0.4.3          
## [22] S4Arrays_1.8.0           curl_6.2.2               DelayedArray_0.34.0     
## [25] abind_1.4-8              BiocParallel_1.42.0      withr_3.0.2             
## [28] grid_4.5.0               colorspace_2.1-1         scales_1.3.0            
## [31] MASS_7.3-65              mcmc_0.9-8               tinytex_0.57            
## [34] cli_3.6.4                mvtnorm_1.3-3            rmarkdown_2.29          
## [37] crayon_1.5.3             httr_1.4.7               rjson_0.2.23            
## [40] cachem_1.1.0             splines_4.5.0            parallel_4.5.0          
## [43] BiocManager_1.30.25      XVector_0.48.0           restfulr_0.0.15         
## [46] vctrs_0.6.5              Matrix_1.7-3             jsonlite_2.0.0          
## [49] SparseM_1.84-2           carData_3.0-5            bookdown_0.43           
## [52] car_3.1-3                MCMCpack_1.7-1           Formula_1.2-5           
## [55] magick_2.8.6             jquerylib_0.1.4          glue_1.8.0              
## [58] codetools_0.2-20         gtable_0.3.6             BiocIO_1.18.0           
## [61] UCSC.utils_1.4.0         munsell_0.5.1            tibble_3.2.1            
## [64] pillar_1.10.2            htmltools_0.5.8.1        quantreg_6.1            
## [67] GenomeInfoDbData_1.2.14  R6_2.6.1                 evaluate_1.0.3          
## [70] lattice_0.22-7           Rsamtools_2.24.0         bslib_0.9.0             
## [73] MatrixModels_0.5-4       Rcpp_1.0.14              coda_0.19-4.1           
## [76] SparseArray_1.8.0        xfun_0.52                pkgconfig_2.0.3