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.229358 -0.79756013 0.71079382  0.48500686 -0.05384609
## ENSMUSG00000000003 1.619208  1.39649835 3.21481717 -2.92635838 -1.91747149
## ENSMUSG00000000028 1.289347 -0.05432853 0.11266552  0.04860176  0.03355944
## ENSMUSG00000000037 1.053122 -2.80965904 7.47428896 -1.85139688 -2.79623749
## ENSMUSG00000000049 1.014522 -0.10578958 0.09803603  0.10502108  0.07339598
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.363294 12.393336 3.259193 2.280623
## ENSMUSG00000000003 24.911560  5.027838 5.874613 8.546568
## ENSMUSG00000000028  8.150715  6.979663 2.941894 2.301342
## ENSMUSG00000000037  9.244015 14.408558 6.787336 2.322025
## ENSMUSG00000000049  5.862191  9.491265 2.878722 1.189985

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.043843088        0.030109795        0.018454756        0.007884151 
## ENSMUSG00000000028 
##        0.005706508

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.242371 -0.48803127 0.3108058  0.36450019  0.077777241
## ENSMUSG00000000003 1.561865  1.53239077 2.3139768 -1.31965981 -2.750816467
## ENSMUSG00000000028 1.277313 -0.02187773 0.1104212  0.02668885 -0.003651618
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.417616 13.607535 3.679199 1.900187
## ENSMUSG00000000003 24.096764  5.591455 5.942960 9.131139
## ENSMUSG00000000028  7.823941  6.474997 2.969883 2.133866
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
##                        Beta_0     Beta_1    Beta_2     Beta_3     Beta_4
## ENSMUSG00000000001  1.9088445  -1.205858  6.792645  -5.441605 -0.3016342
## ENSMUSG00000000003 -0.8338466  -3.597110 10.450255  -6.204598 -0.6458703
## ENSMUSG00000000028  2.3006964 -11.520292 54.271785 -76.637126 34.1651335
##                    Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.407086  6.18872 3.421180 1.519962
## ENSMUSG00000000003 6.015730 10.58480 4.623419 2.889453
## ENSMUSG00000000028 9.407655  6.44427 3.619148 3.388481

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.046689237        0.045513589        0.030662015        0.024507811 
## ENSMUSG00000000049 
##        0.009830224

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-apple-darwin20
## Running under: macOS Monterey 12.7.6
## 
## Matrix products: default
## BLAS:   /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRblas.0.dylib 
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRlapack.dylib;  LAPACK version 3.12.1
## 
## locale:
## [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
## 
## time zone: America/New_York
## tzcode source: internal
## 
## 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