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.254944 -0.68624064 0.68003271 0.27480272 0.005883444
## ENSMUSG00000000003 1.600655 1.70113214 3.03387115 -2.18713018 -2.892133999
## ENSMUSG00000000028 1.326387 -0.02107965 0.09199239 0.05194276 0.015911692
## ENSMUSG00000000037 1.020234 -5.30171036 14.15643078 -6.12772977 -2.737958417
## ENSMUSG00000000049 1.025713 -0.06609732 0.08317797 0.05857072 0.059142265
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.292717 14.636110 3.769187 1.965028
## ENSMUSG00000000003 24.102798 3.608860 5.927295 8.687267
## ENSMUSG00000000028 9.005452 7.272410 3.435735 2.538691
## ENSMUSG00000000037 8.682910 11.476334 7.513607 2.153677
## ENSMUSG00000000049 6.211663 7.958283 3.012329 1.237730
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.069941369 0.033425201 0.014848534 0.006174708
## ENSMUSG00000000028
## 0.005496859
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.265158 -0.54170526 0.4409279 0.26480738 6.396932e-02
## ENSMUSG00000000003 1.558915 1.28673718 3.1616723 -1.54105545 -3.192248e+00
## ENSMUSG00000000028 1.321757 -0.07391479 0.1450761 0.08784784 5.691108e-05
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.935776 14.351450 3.502479 1.716379
## ENSMUSG00000000003 24.254295 6.143700 5.548887 8.561998
## ENSMUSG00000000028 8.600133 6.722864 3.735727 2.597490
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9018146 -3.6837802 19.014784 -23.1873171 7.7165870
## ENSMUSG00000000003 -0.8323118 -0.8014328 2.422573 -0.6748569 -0.8747437
## ENSMUSG00000000028 2.3267649 -21.9483226 106.516574 -153.3382950 69.0295274
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.358080 6.059127 4.407252 1.513229
## ENSMUSG00000000003 6.375058 9.335479 4.701191 2.998822
## ENSMUSG00000000028 9.335452 6.870896 4.669407 3.286500
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)
## ENSMUSG00000000028 ENSMUSG00000000037 ENSMUSG00000000001 ENSMUSG00000000003
## 0.062564514 0.051749523 0.031259081 0.027272725
## ENSMUSG00000000049
## 0.009347043
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 Under development (unstable) (2025-03-01 r87860 ucrt)
## Platform: x86_64-w64-mingw32/x64
## Running under: Windows Server 2022 x64 (build 20348)
##
## Matrix products: default
## LAPACK version 3.12.0
##
## 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.1 SingleCellExperiment_1.29.2
## [3] SummarizedExperiment_1.37.0 Biobase_2.67.0
## [5] GenomicRanges_1.59.1 GenomeInfoDb_1.43.4
## [7] IRanges_2.41.3 S4Vectors_0.45.4
## [9] BiocGenerics_0.53.6 generics_0.1.3
## [11] MatrixGenerics_1.19.1 matrixStats_1.5.0
## [13] mist_0.99.18 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 farver_2.1.2 dplyr_1.1.4
## [4] Biostrings_2.75.4 bitops_1.0-9 fastmap_1.2.0
## [7] RCurl_1.98-1.17 GenomicAlignments_1.43.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.5
## [16] sass_0.4.9 tools_4.5.0 yaml_2.3.10
## [19] rtracklayer_1.67.1 knitr_1.50 labeling_0.4.3
## [22] S4Arrays_1.7.3 curl_6.2.2 DelayedArray_0.33.6
## [25] abind_1.4-8 BiocParallel_1.41.2 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.56
## [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.47.2 restfulr_0.0.15
## [46] vctrs_0.6.5 Matrix_1.7-3 jsonlite_1.9.1
## [49] SparseM_1.84-2 carData_3.0-5 bookdown_0.42
## [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.17.1
## [61] UCSC.utils_1.3.1 munsell_0.5.1 tibble_3.2.1
## [64] pillar_1.10.1 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-6 Rsamtools_2.23.1 bslib_0.9.0
## [73] MatrixModels_0.5-3 Rcpp_1.0.14 coda_0.19-4.1
## [76] SparseArray_1.7.7 xfun_0.51 pkgconfig_2.0.3