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", "small_sampleData_sce.rds", package = "mist"))

2 Step 2: Estimate Parameters Using estiParamSingle

# Estimate parameters for single-group
Dat_sce <- estiParamSingle(
    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.262857 -0.82570332 0.89383567  0.38510549 -0.19217520
## ENSMUSG00000000003 1.622052  1.41981016 3.75355224 -2.49211873 -2.99933414
## ENSMUSG00000000028 1.299436 -0.02152701 0.08056626  0.03691574  0.01987437
## ENSMUSG00000000037 1.051256 -2.95245765 8.18448693 -3.03664622 -2.15036790
## ENSMUSG00000000049 1.031172 -0.06726572 0.08582297  0.07263554  0.07443740
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.986163 12.640717 3.517811 1.909103
## ENSMUSG00000000003 25.525006  5.033066 6.250292 9.664823
## ENSMUSG00000000028  8.220065  7.519973 2.856204 2.323349
## ENSMUSG00000000037  8.722130 14.250900 7.438585 2.224646
## ENSMUSG00000000049  6.166812  9.835495 3.535402 1.262893

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.044579925        0.033251354        0.015868666        0.007449015 
## ENSMUSG00000000028 
##        0.004671754

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 <- readRDS(system.file("extdata", "small_sampleData_sce.rds", package = "mist"))

6 Step 2: Estimate Parameters Using estiParamTwoGroups

# Estimate parameters for both groups
Dat_sce <- estiParamTwo(
    Dat_sce = Dat_sce,
    Dat_name_g1 = "Methy_level_group1",
    Dat_name_g2 = "Methy_level_group2",
    ptime_name_g1 = "pseudotime",
    ptime_name_g2 = "pseudotime_g2"
)

# Check the output
head(rowData(Dat_sce)$mist_pars_group1, n = 3)
##                      Beta_0      Beta_1     Beta_2      Beta_3      Beta_4
## ENSMUSG00000000001 1.259949 -0.68095332 0.75093242  0.30069231 -0.14511948
## ENSMUSG00000000003 1.608136  1.74889594 3.42319369 -2.45517532 -3.08373867
## ENSMUSG00000000028 1.300544 -0.01443848 0.07967227  0.03416459  0.01082773
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.634544 14.498278 3.532108 1.705306
## ENSMUSG00000000003 25.004417  3.700018 6.667944 9.987166
## ENSMUSG00000000028  8.005642  6.835942 2.872354 2.264894
head(rowData(Dat_sce)$mist_pars_group2, n = 3)
##                        Beta_0     Beta_1   Beta_2      Beta_3     Beta_4
## ENSMUSG00000000001  1.8895202 -2.9927836 16.46334 -20.5952085  6.9936294
## ENSMUSG00000000003 -0.8169146 -0.7723547  2.05888  -0.3357485 -0.8569916
## ENSMUSG00000000028  2.3281356 -5.2038293 24.72615 -32.4925459 13.0296047
##                     Sigma2_1  Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001  5.289114  6.574531 3.880439 1.339768
## ENSMUSG00000000003  7.116410 11.432390 4.337887 3.086491
## ENSMUSG00000000028 10.827256  5.851526 4.218941 3.052746

7 Step 3: Perform Differential Methylation Analysis for Two-Group Comparison Using dmTwoGroups

# Perform DM analysis to compare the two groups
Dat_sce <- dmTwoGroups(Dat_sce)

# View the top genomic features with different temporal patterns between groups
head(rowData(Dat_sce)$mist_int_2group)
## ENSMUSG00000000037 ENSMUSG00000000003 ENSMUSG00000000001 ENSMUSG00000000028 
##        0.054991089        0.032006806        0.028626331        0.016534433 
## ENSMUSG00000000049 
##        0.008731698

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) (2024-10-21 r87258)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 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
## 
## 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.1               SingleCellExperiment_1.29.1
##  [3] SummarizedExperiment_1.37.0 Biobase_2.67.0             
##  [5] GenomicRanges_1.59.1        GenomeInfoDb_1.43.2        
##  [7] IRanges_2.41.2              S4Vectors_0.45.2           
##  [9] BiocGenerics_0.53.3         generics_0.1.3             
## [11] MatrixGenerics_1.19.1       matrixStats_1.5.0          
## [13] mist_0.99.17                BiocStyle_2.35.0           
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1         farver_2.1.2             dplyr_1.1.4             
##  [4] Biostrings_2.75.3        bitops_1.0-9             fastmap_1.2.0           
##  [7] RCurl_1.98-1.16          GenomicAlignments_1.43.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.4             
## [16] sass_0.4.9               tools_4.5.0              yaml_2.3.10             
## [19] rtracklayer_1.67.0       knitr_1.49               labeling_0.4.3          
## [22] S4Arrays_1.7.1           curl_6.1.0               DelayedArray_0.33.3     
## [25] abind_1.4-8              BiocParallel_1.41.0      withr_3.0.2             
## [28] grid_4.5.0               colorspace_2.1-1         scales_1.3.0            
## [31] MASS_7.3-64              mcmc_0.9-8               tinytex_0.54            
## [34] cli_3.6.3                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.47.2           restfulr_0.0.15         
## [46] vctrs_0.6.5              Matrix_1.7-1             jsonlite_1.8.9          
## [49] SparseM_1.84-2           carData_3.0-5            bookdown_0.42           
## [52] car_3.1-3                MCMCpack_1.7-1           Formula_1.2-5           
## [55] magick_2.8.5             jquerylib_0.1.4          snow_0.4-4              
## [58] glue_1.8.0               codetools_0.2-20         gtable_0.3.6            
## [61] BiocIO_1.17.1            UCSC.utils_1.3.1         munsell_0.5.1           
## [64] tibble_3.2.1             pillar_1.10.1            htmltools_0.5.8.1       
## [67] quantreg_5.99.1          GenomeInfoDbData_1.2.13  R6_2.5.1                
## [70] evaluate_1.0.3           lattice_0.22-6           Rsamtools_2.23.1        
## [73] bslib_0.8.0              MatrixModels_0.5-3       Rcpp_1.0.14             
## [76] coda_0.19-4.1            SparseArray_1.7.3        xfun_0.50               
## [79] pkgconfig_2.0.3