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.251150 -0.61626932 0.54613938 0.39381622 -0.0558824692
## ENSMUSG00000000003 1.575272 1.69745688 2.27707392 -1.19814446 -3.1021146891
## ENSMUSG00000000028 1.263647 0.00387133 0.08405289 0.03523659 -0.0008817753
## ENSMUSG00000000037 1.037173 -3.54260648 10.04499516 -3.92642023 -2.5784949204
## ENSMUSG00000000049 1.025288 -0.09082411 0.10021899 0.10097271 0.0776860050
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.877484 14.604223 3.341745 2.098561
## ENSMUSG00000000003 24.698764 5.224976 5.502739 8.782775
## ENSMUSG00000000028 7.281718 7.123456 3.281830 1.982603
## ENSMUSG00000000037 8.511323 13.618053 7.033326 2.317321
## ENSMUSG00000000049 5.679876 8.327787 4.105190 1.358472
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.052539138 0.032628419 0.014177131 0.008529279
## ENSMUSG00000000028
## 0.005117429
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.269095 -0.665525983 0.48006403 0.41808292 0.014796448
## ENSMUSG00000000003 1.627367 1.977047017 1.65723731 -1.58832545 -2.364856449
## ENSMUSG00000000028 1.257413 0.009031809 0.08605852 0.03394364 -0.006492632
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.143304 13.719028 3.458299 1.736353
## ENSMUSG00000000003 25.369275 2.571538 6.165069 8.421154
## ENSMUSG00000000028 7.314130 7.255099 2.775601 2.128494
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9033706 -0.866442571 5.1239610 -2.98822880 -1.4352154
## ENSMUSG00000000003 -0.8215392 -2.197668240 6.3166754 -3.08689121 -1.0277497
## ENSMUSG00000000028 2.3102473 0.001104077 0.7137061 -0.04944793 -0.5475033
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.638006 6.071890 3.914513 1.429074
## ENSMUSG00000000003 6.449451 10.622588 4.344451 2.857273
## ENSMUSG00000000028 11.240047 4.849653 3.647838 2.913902
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.048973535 0.036600534 0.025400760 0.009067063
## ENSMUSG00000000028
## 0.002845009
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 RC (2025-04-03 r88103 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.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 farver_2.1.2 dplyr_1.1.4
## [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.0
## [25] abind_1.4-8 BiocParallel_1.43.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.49.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.19.0
## [61] UCSC.utils_1.5.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.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