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.246064 -1.04673226 1.02367671 0.52778596 -0.20097516
## ENSMUSG00000000003 1.545179 1.36227426 2.63514090 -1.29932837 -2.92951531
## ENSMUSG00000000028 1.281824 -0.04903440 0.17108851 0.03308171 -0.01279306
## ENSMUSG00000000037 1.028838 -3.67777676 10.49801115 -4.38186817 -2.46583517
## ENSMUSG00000000049 1.031958 -0.06233893 0.07670723 0.08360788 0.06956525
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.667150 14.694537 2.928585 1.942503
## ENSMUSG00000000003 24.491419 5.423209 6.061727 9.040520
## ENSMUSG00000000028 7.650019 6.533704 3.010469 2.395342
## ENSMUSG00000000037 9.195549 12.529765 6.837504 2.232754
## ENSMUSG00000000049 5.958322 9.494829 3.068395 1.332752
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.053249542 0.031890268 0.019735305 0.007512754
## ENSMUSG00000000028
## 0.006144577
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.252758 -0.72847914 0.7271970 0.32691032 -0.09992272
## ENSMUSG00000000003 1.590748 0.21346570 7.9809142 -7.81404502 -0.71830997
## ENSMUSG00000000028 1.288447 -0.04392353 0.1659242 0.02179529 -0.01827020
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.159453 15.024230 3.672005 1.732141
## ENSMUSG00000000003 25.298969 8.133583 5.025178 8.895314
## ENSMUSG00000000028 7.728085 6.944934 3.382666 2.181663
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9012285 -2.197820 11.363127 -11.708073 2.382945
## ENSMUSG00000000003 -0.8233923 -2.443719 7.096286 -3.450223 -1.214735
## ENSMUSG00000000028 2.3045281 -3.257554 14.124731 -17.047119 6.300813
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.310023 6.810922 3.703171 1.426427
## ENSMUSG00000000003 7.208322 10.916552 5.142645 2.770846
## ENSMUSG00000000028 10.681545 5.679282 3.242529 3.131001
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 ENSMUSG00000000028
## 0.040622368 0.034898783 0.027621137 0.010200423
## ENSMUSG00000000049
## 0.009947558
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 (2025-04-11 ucrt)
## Platform: x86_64-w64-mingw32/x64
## Running under: Windows Server 2022 x64 (build 20348)
##
## Matrix products: default
## LAPACK version 3.12.1
##
## locale:
## [1] LC_COLLATE=C
## [2] LC_CTYPE=English_United States.utf8
## [3] LC_MONETARY=English_United States.utf8
## [4] LC_NUMERIC=C
## [5] LC_TIME=English_United States.utf8
##
## 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.4
## [7] IRanges_2.43.0 S4Vectors_0.47.0
## [9] BiocGenerics_0.55.0 generics_0.1.4
## [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.1 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.1 curl_6.2.3 DelayedArray_0.35.1
## [25] RColorBrewer_1.1-3 abind_1.4-8 BiocParallel_1.43.3
## [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] R6_2.6.1 evaluate_1.0.3 lattice_0.22-7
## [70] Rsamtools_2.25.0 bslib_0.9.0 MatrixModels_0.5-4
## [73] Rcpp_1.0.14 coda_0.19-4.1 SparseArray_1.9.0
## [76] xfun_0.52 pkgconfig_2.0.3