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.242258 -0.6224564 0.67967881 0.32324240 -0.14276644
## ENSMUSG00000000003 1.502365 1.4918870 4.06645447 -2.85632808 -3.01950917
## ENSMUSG00000000028 1.267887 -0.0158547 0.09417752 0.02728282 0.03127849
## ENSMUSG00000000037 1.024013 -4.2421383 11.96386127 -5.74654225 -1.95860429
## ENSMUSG00000000049 1.017909 -0.1136836 0.11282333 0.10780710 0.07018410
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.631366 13.575626 3.537522 1.841812
## ENSMUSG00000000003 26.399930 7.216271 7.164505 9.651759
## ENSMUSG00000000028 7.849256 7.437372 3.053947 2.334055
## ENSMUSG00000000037 8.332989 13.097018 7.063658 2.269228
## ENSMUSG00000000049 5.650790 9.470500 3.129593 1.234002
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.058116750 0.037414456 0.013370305 0.008219157
## ENSMUSG00000000028
## 0.005558315
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.235635 -0.764457840 0.6779759 0.41321849 -0.03811891
## ENSMUSG00000000003 1.610273 2.002030574 3.2453532 -2.88297966 -2.80459836
## ENSMUSG00000000028 1.274080 -0.008190074 0.1025516 0.02024307 0.02125887
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.454641 13.907133 3.503710 1.961562
## ENSMUSG00000000003 28.612966 3.611619 6.653995 8.759865
## ENSMUSG00000000028 7.980330 7.284778 3.326566 2.380037
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.8954555 -1.33522987 7.5904504 -7.15178471 0.7415187
## ENSMUSG00000000003 -0.8131198 -0.31602501 0.9449495 0.02563909 -0.5633839
## ENSMUSG00000000028 2.3312309 -0.09544047 1.0415252 -0.39738174 -0.4125042
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.420610 5.852837 3.062897 1.280664
## ENSMUSG00000000003 8.511375 9.909064 5.159043 2.798039
## ENSMUSG00000000028 11.665818 6.630344 3.174623 3.350426
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.049203254 0.032734698 0.026651313 0.010588164
## ENSMUSG00000000028
## 0.003146898
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-04 r88126)
## Platform: aarch64-apple-darwin20
## Running under: macOS Ventura 13.7.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
##
## locale:
## [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
##
## 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.30.0
## [3] SummarizedExperiment_1.38.0 Biobase_2.68.0
## [5] GenomicRanges_1.60.0 GenomeInfoDb_1.44.0
## [7] IRanges_2.42.0 S4Vectors_0.46.0
## [9] BiocGenerics_0.54.0 generics_0.1.3
## [11] MatrixGenerics_1.20.0 matrixStats_1.5.0
## [13] mist_1.0.0 BiocStyle_2.36.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.2.1 farver_2.1.2 dplyr_1.1.4
## [4] Biostrings_2.76.0 bitops_1.0-9 fastmap_1.2.0
## [7] RCurl_1.98-1.17 GenomicAlignments_1.44.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.68.0 knitr_1.50 labeling_0.4.3
## [22] S4Arrays_1.8.0 curl_6.2.2 DelayedArray_0.34.0
## [25] abind_1.4-8 BiocParallel_1.42.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.48.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.18.0
## [61] UCSC.utils_1.4.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.24.0 bslib_0.9.0
## [73] MatrixModels_0.5-4 Rcpp_1.0.14 coda_0.19-4.1
## [76] SparseArray_1.8.0 xfun_0.52 pkgconfig_2.0.3