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.229565 -0.906256230 0.9449425 0.35219988 -0.10106566
## ENSMUSG00000000003 1.577763 1.731825389 1.7484500 -0.93807626 -2.85938265
## ENSMUSG00000000028 1.300256 0.001863788 0.1011794 0.02698799 -0.00596132
## ENSMUSG00000000037 1.034239 -3.035795349 8.4141149 -3.02967490 -2.31930084
## ENSMUSG00000000049 1.019115 -0.094677922 0.1213672 0.09206358 0.06155379
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.446638 13.691372 2.975291 1.843944
## ENSMUSG00000000003 25.142947 3.850606 4.932253 9.010189
## ENSMUSG00000000028 8.051416 8.230083 3.207845 2.425138
## ENSMUSG00000000037 8.414027 12.808300 6.652033 2.324699
## ENSMUSG00000000049 5.719849 8.771707 3.386213 1.337142
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.046127342 0.031131709 0.017782415 0.008339259
## ENSMUSG00000000028
## 0.005156739
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.234728 -0.755156479 0.83968604 0.30669009 -0.14096961
## ENSMUSG00000000003 1.614334 2.016753535 1.14998350 -1.09656827 -2.38203882
## ENSMUSG00000000028 1.295590 0.008139507 0.09355377 0.02067561 -0.01258234
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.416661 13.691300 3.321143 1.725017
## ENSMUSG00000000003 25.582497 2.993087 5.228878 8.759563
## ENSMUSG00000000028 7.976525 8.361564 3.614812 2.267689
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.909795 -0.9023612 5.114487 -3.265287 -1.0988102
## ENSMUSG00000000003 -0.852604 -1.1100448 3.286372 -1.426551 -0.7038767
## ENSMUSG00000000028 2.319481 -2.6649750 11.983062 -13.691708 4.4514022
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.730759 6.271297 3.069641 1.437109
## ENSMUSG00000000003 6.545853 10.123953 4.588248 3.216118
## ENSMUSG00000000028 10.944550 5.346765 3.831982 3.309672
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.047265158 0.028838316 0.022959994 0.010442367
## ENSMUSG00000000049
## 0.009436176
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-01-21 r87610 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.1
## [3] SummarizedExperiment_1.37.0 Biobase_2.67.0
## [5] GenomicRanges_1.59.1 GenomeInfoDb_1.43.3
## [7] IRanges_2.41.2 S4Vectors_0.45.2
## [9] BiocGenerics_0.53.4 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.3 bitops_1.0-9 fastmap_1.2.0
## [7] RCurl_1.98-1.16 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.0 knitr_1.49 labeling_0.4.3
## [22] S4Arrays_1.7.1 curl_6.1.0 DelayedArray_0.33.4
## [25] abind_1.4-8 BiocParallel_1.41.0 withr_3.0.2
## [28] grid_4.5.0 colorspace_2.1-1 scales_1.3.0
## [31] MASS_7.3-64 mcmc_0.9-8 tinytex_0.54
## [34] cli_3.6.3 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-1 jsonlite_1.8.9
## [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.5 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_5.99.1
## [67] GenomeInfoDbData_1.2.13 R6_2.5.1 evaluate_1.0.3
## [70] lattice_0.22-6 Rsamtools_2.23.1 bslib_0.8.0
## [73] MatrixModels_0.5-3 Rcpp_1.0.14 coda_0.19-4.1
## [76] SparseArray_1.7.4 xfun_0.50 pkgconfig_2.0.3