estiParamdmSingleplotGeneestiParamdmTwoGroupsmist (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.2676800 -0.598878847 0.5100346 0.34358245 -0.005613921
## ENSMUSG00000000003 1.6181894 1.495399767 3.0556513 -3.18590289 -1.554493079
## ENSMUSG00000000028 1.2936013 -0.008897203 0.1012861 0.04795074 -0.017341694
## ENSMUSG00000000037 0.9915283 -5.956975464 17.1808289 -9.39725965 -1.876998969
## ENSMUSG00000000049 1.0190144 -0.105142793 0.1280139 0.09604899 0.042100028
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.432892 14.480615 3.268779 1.945273
## ENSMUSG00000000003 25.695878 5.750383 6.302808 9.662156
## ENSMUSG00000000028 8.407409 8.195303 2.857582 2.391422
## ENSMUSG00000000037 7.882963 14.108995 7.045180 2.233891
## ENSMUSG00000000049 5.968150 9.445924 3.136424 1.192521
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.077244885 0.029381493 0.013297119 0.007622820
## ENSMUSG00000000028
## 0.005178364
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.262652 -0.452434383 0.3431416 0.27430053 0.08524158
## ENSMUSG00000000003 1.674666 1.393614284 3.4069355 -3.21173312 -1.83042343
## ENSMUSG00000000028 1.298124 -0.009471835 0.1142839 0.04447859 -0.03703512
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.671606 16.180691 3.619820 1.992995
## ENSMUSG00000000003 25.697486 3.932150 5.527550 10.250138
## ENSMUSG00000000028 8.008246 8.066523 3.113934 2.268636
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9205187 -1.557707 8.648568 -7.6387291 0.3656197
## ENSMUSG00000000003 -0.8275165 -1.040965 2.865798 -0.7613642 -0.9971975
## ENSMUSG00000000028 2.3610669 -2.429621 11.186509 -14.0509483 5.4355475
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.479661 6.481981 3.967230 1.557449
## ENSMUSG00000000003 6.900045 13.135785 5.203956 2.840029
## ENSMUSG00000000028 11.553943 4.913955 3.770720 3.202984
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.066868571 0.027146463 0.025592710 0.008996873
## ENSMUSG00000000028
## 0.007147953
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.6.0 RC (2026-04-17 r89917)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.24-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 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_4.0.3 SingleCellExperiment_1.35.0
## [3] SummarizedExperiment_1.43.0 Biobase_2.73.0
## [5] GenomicRanges_1.65.0 Seqinfo_1.3.0
## [7] IRanges_2.47.0 S4Vectors_0.51.0
## [9] BiocGenerics_0.59.0 generics_0.1.4
## [11] MatrixGenerics_1.25.0 matrixStats_1.5.0
## [13] mist_1.5.0 BiocStyle_2.41.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.1 farver_2.1.2
## [4] Biostrings_2.81.0 S7_0.2.2 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.18 GenomicAlignments_1.49.0
## [10] XML_3.99-0.23 digest_0.6.39 lifecycle_1.0.5
## [13] survival_3.8-6 magrittr_2.0.5 compiler_4.6.0
## [16] rlang_1.2.0 sass_0.4.10 tools_4.6.0
## [19] yaml_2.3.12 rtracklayer_1.73.0 knitr_1.51
## [22] labeling_0.4.3 S4Arrays_1.13.0 curl_7.1.0
## [25] DelayedArray_0.39.0 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.47.0 withr_3.0.2 grid_4.6.0
## [31] scales_1.4.0 MASS_7.3-65 mcmc_0.9-8
## [34] tinytex_0.59 dichromat_2.0-0.1 cli_3.6.6
## [37] mvtnorm_1.3-7 rmarkdown_2.31 crayon_1.5.3
## [40] otel_0.2.0 httr_1.4.8 rjson_0.2.23
## [43] cachem_1.1.0 splines_4.6.0 parallel_4.6.0
## [46] BiocManager_1.30.27 XVector_0.53.0 restfulr_0.0.16
## [49] vctrs_0.7.3 Matrix_1.7-5 jsonlite_2.0.0
## [52] SparseM_1.84-2 carData_3.0-6 bookdown_0.46
## [55] car_3.1-5 MCMCpack_1.7-1 Formula_1.2-5
## [58] magick_2.9.1 jquerylib_0.1.4 glue_1.8.1
## [61] codetools_0.2-20 gtable_0.3.6 BiocIO_1.23.0
## [64] tibble_3.3.1 pillar_1.11.1 htmltools_0.5.9
## [67] quantreg_6.1 R6_2.6.1 evaluate_1.0.5
## [70] lattice_0.22-9 Rsamtools_2.29.0 cigarillo_1.3.0
## [73] bslib_0.10.0 MatrixModels_0.5-4 Rcpp_1.1.1-1.1
## [76] coda_0.19-4.1 SparseArray_1.13.0 xfun_0.57
## [79] pkgconfig_2.0.3