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.246536 -0.60276933 0.43725142 0.36695993 0.10442036
## ENSMUSG00000000003 1.616021 1.64864662 2.04802630 -1.30756888 -2.66930977
## ENSMUSG00000000028 1.303839 0.01640244 0.07552145 0.03472030 -0.00919962
## ENSMUSG00000000037 1.030833 -4.28095498 12.03480766 -5.33645684 -2.41328212
## ENSMUSG00000000049 1.033205 -0.06294789 0.05955992 0.09157738 0.07914562
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.265453 14.642869 3.484166 2.385946
## ENSMUSG00000000003 26.420517 4.954308 5.828470 9.537878
## ENSMUSG00000000028 8.052628 8.623274 3.761634 2.349188
## ENSMUSG00000000037 8.510447 14.051668 7.761888 2.413510
## ENSMUSG00000000049 6.446328 9.044015 3.076302 1.391624
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.059857691 0.029988308 0.015191511 0.007427256
## ENSMUSG00000000028
## 0.004836925
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.267955 -0.542186546 0.42490014 0.32521796 0.04777201
## ENSMUSG00000000003 1.572868 1.599509728 2.14581318 -1.21342912 -2.81759409
## ENSMUSG00000000028 1.308759 0.002820363 0.06428554 0.04945468 0.01406688
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 6.572918 15.294667 3.261281 2.101758
## ENSMUSG00000000003 25.153584 6.300158 5.289404 10.068562
## ENSMUSG00000000028 8.470304 8.217529 3.482690 2.352629
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.9168166 -0.8273217 4.969761 -2.601536 -1.6993554
## ENSMUSG00000000003 -0.8079073 -1.4453384 4.209237 -2.196539 -0.5054973
## ENSMUSG00000000028 2.3457638 -3.6881028 16.872630 -21.640109 8.5908510
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.802272 7.981131 4.100534 1.597460
## ENSMUSG00000000003 7.289530 11.373625 4.657467 3.033767
## ENSMUSG00000000028 11.482579 5.783470 3.333607 3.209432
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 ENSMUSG00000000028
## 0.034947938 0.030388788 0.023699210 0.010886676
## ENSMUSG00000000049
## 0.008980788
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 beta (2025-04-02 r88102)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 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_3.5.2 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.6
## [16] sass_0.4.10 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.5 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_2.0.0
## [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.2
## [61] UCSC.utils_1.3.1 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.23.1 bslib_0.9.0
## [73] MatrixModels_0.5-4 Rcpp_1.0.14 coda_0.19-4.1
## [76] SparseArray_1.7.7 xfun_0.52 pkgconfig_2.0.3