--- title: "mist:methylation inference for single-cell along trajectory" output: BiocStyle::html_document: toc: true toc_depth: 2 vignette: > %\VignetteIndexEntry{mist_vignette} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} --- ## Introduction `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. ## Installation To install the latest version of `mist`, run the following commands: ```{r, eval = FALSE, warning=FALSE, message=FALSE} if (!requireNamespace("BiocManager", quietly = TRUE)) { install.packages("BiocManager") } # Install mist from GitHub BiocManager::install("https://github.com/dxd429/mist") ``` From Bioconductor: ```{r, eval = FALSE, warning=FALSE, message=FALSE} 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: ```{r eval=FALSE,warning=FALSE,message=FALSE} library(mist) vignette("mist") ``` ## Example Workflow for Single-Group Analysis In this section, we will estimate parameters and perform differential methylation analysis using single-group data. # Step 1: Load Example Data Here we load the example data from GSE121708. ```{r, eval=TRUE, message=FALSE, warning=FALSE} library(mist) library(SingleCellExperiment) # Load sample scDNAm data Dat_sce <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) ``` # Step 2: Estimate Parameters Using `estiParam` ```{r, eval=TRUE, message=FALSE, warning=FALSE} # 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) ``` # Step 3: Perform Differential Methylation Analysis Using `dmSingle` ```{r, eval=TRUE, message=FALSE, warning=FALSE} # 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) ``` # Step 4: Perform Differential Methylation Analysis Using `plotGene` ```{r eval=TRUE, fig.height=5, fig.width=8, message=FALSE, warning=FALSE, out.width='100%'} # 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") ``` ## Example Workflow for Two-Group Analysis In this section, we will estimate parameters and perform DM analysis using data from two phenotypic groups. # Step 1: Load Two-Group Data ```{r, eval=TRUE, message=FALSE, warning=FALSE} # 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")) ``` # Step 2: Estimate Parameters Using `estiParam` ```{r, eval=TRUE, message=FALSE, warning=FALSE} # 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) head(rowData(Dat_sce_g2)$mist_pars, n = 3) ``` # Step 3: Perform Differential Methylation Analysis for Two-Group Comparison Using `dmTwoGroups` ```{r, eval=TRUE, message=FALSE, warning=FALSE} # 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) ``` ## Conclusion `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. # Session info {.unnumbered} ```{r sessionInfo, echo=FALSE} sessionInfo() ```