estiParam
dmSingle
plotGene
estiParam
dmTwoGroups
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.
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")
In this section, we will estimate parameters and perform differential methylation analysis using single-group 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"))
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.268944 -3.818520e-01 0.2914901 0.28139231 0.04596664
## ENSMUSG00000000003 1.536817 9.611772e-01 3.1903305 -1.68775946 -2.67142591
## ENSMUSG00000000028 1.291467 -3.080192e-05 0.1071733 0.04867304 -0.02096125
## ENSMUSG00000000037 1.048652 -2.930170e+00 7.5497965 -1.13978426 -3.47385148
## ENSMUSG00000000049 1.015904 -2.043675e-01 0.2003745 0.09771344 0.09414582
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.010588 13.693273 3.561577 1.895999
## ENSMUSG00000000003 23.944624 8.232259 5.111561 8.478882
## ENSMUSG00000000028 7.855017 6.991929 3.295653 2.606817
## ENSMUSG00000000037 8.980761 11.935988 8.357927 2.295164
## ENSMUSG00000000049 6.007453 8.490531 3.048808 1.282860
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.046967872 0.031149846 0.010903139 0.009118700
## ENSMUSG00000000028
## 0.005671512
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")
In this section, we will estimate parameters and perform DM analysis using data from two phenotypic groups.
# 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"))
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.258047 -0.488922810 0.3055507 0.36697657 0.07563823
## ENSMUSG00000000003 1.509964 1.613839833 1.9517978 -1.64480205 -2.14040747
## ENSMUSG00000000028 1.298791 -0.002039645 0.1396054 0.04658121 -0.05638440
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.693721 15.830394 2.987427 1.794977
## ENSMUSG00000000003 24.121958 4.366282 6.092656 8.346126
## ENSMUSG00000000028 7.686776 7.074796 3.538004 2.451185
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9081492 0.05127459 2.4760443 -1.030953035 -1.655542
## ENSMUSG00000000003 -0.8297118 -0.92028330 2.7930387 -1.025659000 -0.757990
## ENSMUSG00000000028 2.3467502 0.02682188 0.5189964 0.002194115 -0.388745
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.551244 8.211444 3.466165 1.272790
## ENSMUSG00000000003 6.306315 10.356262 4.548397 3.171298
## ENSMUSG00000000028 11.242020 5.927492 3.463530 3.335002
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.045390762 0.028042676 0.022744730 0.010874014
## ENSMUSG00000000028
## 0.001933617
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.
## R version 4.5.0 Patched (2025-04-21 r88169)
## Platform: aarch64-apple-darwin20
## Running under: macOS Ventura 13.7.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/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.31.0
## [3] SummarizedExperiment_1.39.0 Biobase_2.69.0
## [5] GenomicRanges_1.61.0 GenomeInfoDb_1.45.0
## [7] IRanges_2.43.0 S4Vectors_0.47.0
## [9] BiocGenerics_0.55.0 generics_0.1.3
## [11] MatrixGenerics_1.21.0 matrixStats_1.5.0
## [13] mist_1.1.0 BiocStyle_2.37.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.1.4 farver_2.1.2
## [4] Biostrings_2.77.0 bitops_1.0-9 fastmap_1.2.0
## [7] RCurl_1.98-1.17 GenomicAlignments_1.45.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.69.0 knitr_1.50 labeling_0.4.3
## [22] S4Arrays_1.9.0 curl_6.2.2 DelayedArray_0.35.1
## [25] RColorBrewer_1.1-3 abind_1.4-8 BiocParallel_1.43.0
## [28] withr_3.0.2 grid_4.5.0 scales_1.4.0
## [31] MASS_7.3-65 mcmc_0.9-8 tinytex_0.57
## [34] dichromat_2.0-0.1 cli_3.6.5 mvtnorm_1.3-3
## [37] rmarkdown_2.29 crayon_1.5.3 httr_1.4.7
## [40] rjson_0.2.23 cachem_1.1.0 splines_4.5.0
## [43] parallel_4.5.0 BiocManager_1.30.25 XVector_0.49.0
## [46] restfulr_0.0.15 vctrs_0.6.5 Matrix_1.7-3
## [49] jsonlite_2.0.0 SparseM_1.84-2 carData_3.0-5
## [52] bookdown_0.43 car_3.1-3 MCMCpack_1.7-1
## [55] Formula_1.2-5 magick_2.8.6 jquerylib_0.1.4
## [58] glue_1.8.0 codetools_0.2-20 gtable_0.3.6
## [61] BiocIO_1.19.0 UCSC.utils_1.5.0 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.25.0 bslib_0.9.0
## [73] MatrixModels_0.5-4 Rcpp_1.0.14 coda_0.19-4.1
## [76] SparseArray_1.9.0 xfun_0.52 pkgconfig_2.0.3