estiParamSingle
dmSingle
plotGene
estiParamTwoGroups
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", "small_sampleData_sce.rds", package = "mist"))
estiParamSingle
# Estimate parameters for single-group
Dat_sce <- estiParamSingle(
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.239793 -0.812292554 0.7905570 0.36113628 -0.05642390
## ENSMUSG00000000003 1.636478 1.507091937 3.5987541 -2.59831691 -2.89088729
## ENSMUSG00000000028 1.308739 0.005595154 0.1052620 0.02305620 -0.02521009
## ENSMUSG00000000037 1.006877 -4.894850046 13.0707335 -4.61242636 -3.62477282
## ENSMUSG00000000049 1.005655 -0.158609861 0.3103474 -0.08733549 0.13228932
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.184032 15.183924 3.114117 1.826032
## ENSMUSG00000000003 25.416840 2.596218 6.030290 8.890647
## ENSMUSG00000000028 8.170899 6.888983 3.502161 2.283669
## ENSMUSG00000000037 8.567505 14.261554 6.471687 2.197505
## ENSMUSG00000000049 5.350651 8.705501 2.954788 1.314802
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.068134864 0.033241657 0.016583384 0.009446557
## ENSMUSG00000000028
## 0.004667178
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 <- readRDS(system.file("extdata", "small_sampleData_sce.rds", package = "mist"))
estiParamTwoGroups
# Estimate parameters for both groups
Dat_sce <- estiParamTwo(
Dat_sce = Dat_sce,
Dat_name_g1 = "Methy_level_group1",
Dat_name_g2 = "Methy_level_group2",
ptime_name_g1 = "pseudotime",
ptime_name_g2 = "pseudotime_g2"
)
# Check the output
head(rowData(Dat_sce)$mist_pars_group1, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.232429 -0.89294484 1.090386 0.38436298 -0.320267435
## ENSMUSG00000000003 1.616777 1.87283400 2.299257 -2.47445437 -1.995626860
## ENSMUSG00000000028 1.273690 -0.01300002 0.134100 0.04248286 -0.003700711
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.301442 13.851255 3.004876 1.865224
## ENSMUSG00000000003 25.768155 3.428453 6.035565 8.754384
## ENSMUSG00000000028 7.520042 8.172288 3.303499 2.762588
head(rowData(Dat_sce)$mist_pars_group2, n = 3)
## Beta_0 Beta_1 Beta_2 Beta_3 Beta_4
## ENSMUSG00000000001 1.8935580 -0.9014650 5.934080 -4.8783653 -0.3131641
## ENSMUSG00000000003 -0.8394897 -1.5851271 4.290100 -1.3204273 -1.3365637
## ENSMUSG00000000028 2.3327407 -0.6291666 2.476819 -0.6818228 -1.0461705
## Sigma2_1 Sigma2_2 Sigma2_3 Sigma2_4
## ENSMUSG00000000001 5.616122 6.646827 3.938307 1.303361
## ENSMUSG00000000003 6.661643 10.441203 4.594972 2.874717
## ENSMUSG00000000028 10.943867 4.978352 3.707316 3.413476
dmTwoGroups
# Perform DM analysis to compare the two groups
Dat_sce <- dmTwoGroups(Dat_sce)
# View the top genomic features with different temporal patterns between groups
head(rowData(Dat_sce)$mist_int_2group)
## ENSMUSG00000000037 ENSMUSG00000000003 ENSMUSG00000000001 ENSMUSG00000000049
## 0.04206112 0.03358577 0.02480387 0.01174847
## ENSMUSG00000000028
## 0.00482716
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) (2024-10-26 r87273 ucrt)
## Platform: x86_64-w64-mingw32/x64
## Running under: Windows Server 2022 x64 (build 20348)
##
## Matrix products: default
##
##
## 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.2
## [7] IRanges_2.41.2 S4Vectors_0.45.2
## [9] BiocGenerics_0.53.3 generics_0.1.3
## [11] MatrixGenerics_1.19.1 matrixStats_1.5.0
## [13] mist_0.99.17 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.4
## [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.3
## [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 snow_0.4-4
## [58] glue_1.8.0 codetools_0.2-20 gtable_0.3.6
## [61] BiocIO_1.17.1 UCSC.utils_1.3.1 munsell_0.5.1
## [64] tibble_3.2.1 pillar_1.10.1 htmltools_0.5.8.1
## [67] quantreg_5.99.1 GenomeInfoDbData_1.2.13 R6_2.5.1
## [70] evaluate_1.0.3 lattice_0.22-6 Rsamtools_2.23.1
## [73] bslib_0.8.0 MatrixModels_0.5-3 Rcpp_1.0.14
## [76] coda_0.19-4.1 SparseArray_1.7.3 xfun_0.50
## [79] pkgconfig_2.0.3