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.268073 -0.641346066 0.62872696 0.26944551 0.01082976
## ENSMUSG00000000003 1.632751 1.445787891 1.65496070 -1.30422414 -2.01886258
## ENSMUSG00000000028 1.300790 0.001772445 0.07561411 0.02997227 0.01351331
## ENSMUSG00000000037 1.056125 -4.662113132 12.03859243 -3.99812648 -3.41446472
## ENSMUSG00000000049 1.034762 -0.058868386 0.07324100 0.07763724 0.06958267
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.422380 12.434563 3.647117 1.810457
## ENSMUSG00000000003 25.958193 6.388699 5.647371 8.519822
## ENSMUSG00000000028 8.414973 7.287698 3.465817 2.536087
## ENSMUSG00000000037 9.358945 15.943890 6.630018 2.191420
## ENSMUSG00000000049 6.327892 8.485897 3.075552 1.373523
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.062641849 0.026284848 0.014112123 0.007227834
## ENSMUSG00000000028
## 0.004864650
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.261034 -0.76943826 0.87074603 0.25452286 -0.07990514
## ENSMUSG00000000003 1.662519 1.40815265 1.74971303 -1.18990716 -2.21195650
## ENSMUSG00000000028 1.303897 -0.01451171 0.08711917 0.04415376 0.03647076
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.346253 12.371657 3.512322 1.859134
## ENSMUSG00000000003 25.547602 3.110734 5.810756 9.133135
## ENSMUSG00000000028 7.972939 7.609205 3.748072 2.631517
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9174541 -0.9308021 5.996186 -4.3859365 -0.87186378
## ENSMUSG00000000003 -0.8373239 -0.6674600 2.066530 -0.7120709 -0.61386741
## ENSMUSG00000000028 2.2879051 -1.1854237 4.320329 -2.9340149 -0.08564002
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.533410 5.997445 3.502361 1.381580
## ENSMUSG00000000003 6.967639 9.993848 4.737226 3.194577
## ENSMUSG00000000028 10.889958 5.549569 3.887735 3.306782
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 ENSMUSG00000000001 ENSMUSG00000000003 ENSMUSG00000000049
## 0.067221524 0.025516772 0.023276128 0.011661651
## ENSMUSG00000000028
## 0.006141644
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.2 (2025-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-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.2 SingleCellExperiment_1.32.0
## [3] SummarizedExperiment_1.40.0 Biobase_2.70.0
## [5] GenomicRanges_1.62.1 Seqinfo_1.0.0
## [7] IRanges_2.44.0 S4Vectors_0.48.0
## [9] BiocGenerics_0.56.0 generics_0.1.4
## [11] MatrixGenerics_1.22.0 matrixStats_1.5.0
## [13] mist_1.2.2 BiocStyle_2.38.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 dplyr_1.2.0 farver_2.1.2
## [4] Biostrings_2.78.0 S7_0.2.1 bitops_1.0-9
## [7] fastmap_1.2.0 RCurl_1.98-1.17 GenomicAlignments_1.46.0
## [10] XML_3.99-0.20 digest_0.6.39 lifecycle_1.0.5
## [13] survival_3.8-6 magrittr_2.0.4 compiler_4.5.2
## [16] rlang_1.1.7 sass_0.4.10 tools_4.5.2
## [19] yaml_2.3.12 rtracklayer_1.70.1 knitr_1.51
## [22] labeling_0.4.3 S4Arrays_1.10.1 curl_7.0.0
## [25] DelayedArray_0.36.0 RColorBrewer_1.1-3 abind_1.4-8
## [28] BiocParallel_1.44.0 withr_3.0.2 grid_4.5.2
## [31] scales_1.4.0 MASS_7.3-65 mcmc_0.9-8
## [34] tinytex_0.58 dichromat_2.0-0.1 cli_3.6.5
## [37] mvtnorm_1.3-3 rmarkdown_2.30 crayon_1.5.3
## [40] otel_0.2.0 httr_1.4.7 rjson_0.2.23
## [43] cachem_1.1.0 splines_4.5.2 parallel_4.5.2
## [46] BiocManager_1.30.27 XVector_0.50.0 restfulr_0.0.16
## [49] vctrs_0.7.1 Matrix_1.7-4 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.0 jquerylib_0.1.4 glue_1.8.0
## [61] codetools_0.2-20 gtable_0.3.6 BiocIO_1.20.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-7 Rsamtools_2.26.0 cigarillo_1.0.0
## [73] bslib_0.10.0 MatrixModels_0.5-4 Rcpp_1.1.1
## [76] coda_0.19-4.1 SparseArray_1.10.8 xfun_0.56
## [79] pkgconfig_2.0.3