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.257692 -0.6952493382 0.61567251 0.53515602 -0.2344914767
## ENSMUSG00000000003 1.581106 0.8747367422 3.81659066 -1.70037645 -3.2294237166
## ENSMUSG00000000028 1.284364 -0.0001793846 0.08159529 0.03056432 0.0001473274
## ENSMUSG00000000037 1.056139 -1.7023323624 4.66504214 -0.32599238 -2.6042201663
## ENSMUSG00000000049 1.031498 -0.0829472601 0.11482518 0.07411521 0.0443292554
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.511900 15.427292 3.299756 1.589265
## ENSMUSG00000000003 24.565991 6.307187 5.536823 9.050213
## ENSMUSG00000000028 8.324742 7.357275 2.935553 2.299307
## ENSMUSG00000000037 9.376057 11.740765 6.979171 2.227382
## ENSMUSG00000000049 6.073160 8.638733 3.099675 1.195544
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.033068314 0.032421935 0.013717945 0.007009606
## ENSMUSG00000000028
## 0.004656495
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.258457 -0.94912891 0.9374746 0.69322657 -0.421447748
## ENSMUSG00000000003 1.626355 1.52270690 3.3993763 -2.32428415 -3.001825689
## ENSMUSG00000000028 1.290358 -0.02222768 0.1070045 0.04051844 -0.004007954
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.691391 13.641890 3.284727 1.802250
## ENSMUSG00000000003 26.513402 2.946628 6.139663 8.403902
## ENSMUSG00000000028 7.834738 8.306904 3.154024 2.149072
head(rowData(Dat_sce_g2)$mist_pars, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4 Sigma2_1
## ENSMUSG00000000001 1.932540 -0.2894189 3.366104 -2.050875 -1.1790194 6.258029
## ENSMUSG00000000003 -0.828484 -1.1237174 3.507944 -1.332851 -1.0222934 7.184402
## ENSMUSG00000000028 2.391034 -0.6979583 3.626524 -3.433485 0.6108754 12.248834
## Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.184055 3.333934 1.403501
## ENSMUSG00000000003 9.709888 4.291719 3.281514
## ENSMUSG00000000028 4.413611 4.053735 3.273897
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.030358735 0.028399948 0.024677912 0.011094711
## ENSMUSG00000000028
## 0.004754208
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-22 r87618)
## Platform: x86_64-apple-darwin20
## Running under: macOS Monterey 12.7.6
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.0
##
## 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.1 SingleCellExperiment_1.29.1
## [3] SummarizedExperiment_1.37.0 Biobase_2.67.0
## [5] GenomicRanges_1.59.1 GenomeInfoDb_1.43.4
## [7] IRanges_2.41.2 S4Vectors_0.45.2
## [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.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.2 curl_6.2.0 DelayedArray_0.33.5
## [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-2 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_6.00
## [67] GenomeInfoDbData_1.2.13 R6_2.5.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.5 xfun_0.50 pkgconfig_2.0.3