Title: | Differential Methylation Analysis for scDNAm Data |
---|---|
Description: | mist (Methylation Inference for Single-cell along Trajectory) is a hierarchical Bayesian framework for modeling DNA methylation trajectories and performing differential methylation (DM) analysis in single-cell DNA methylation (scDNAm) data. It estimates developmental-stage-specific variations, identifies genomic features with drastic changes along pseudotime, and, for two phenotypic groups, detects features with distinct temporal methylation patterns. mist uses Gibbs sampling to estimate parameters for temporal changes and stage-specific variations. |
Authors: | Daoyu Duan [aut, cre] |
Maintainer: | Daoyu Duan <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.99.18 |
Built: | 2025-02-22 03:26:37 UTC |
Source: | https://github.com/bioc/mist |
This function performs DM analysis to identify genomic features showing drastic changes along pseudotime, by modeling the methylation changes along pseudotime and estimating the developmental-stage-specific biological variations.
dmSingle(Dat_sce, BPPARAM = MulticoreParam())
dmSingle(Dat_sce, BPPARAM = MulticoreParam())
Dat_sce |
The updated sce object with A numeric matrix of estimated parameters for all genomic features in the rowData, including:
|
BPPARAM |
A |
The updated sce object with A named numeric vector where each value corresponds to a genomic feature (e.g., a gene) in the rowData. The values represent the minimum area between the fitted curve of scDNA methylation levels along pseudotime and a constant horizontal line for each feature. This metric measures the magnitude of methylation level changes along pseudotime. The values are sorted in descending order, with larger values indicating more drastic changes.
library(mist) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) dm_sce <- dmSingle(Dat_sce_new)
library(mist) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) dm_sce <- dmSingle(Dat_sce_new)
This function performs differential methylation (DM) analysis to identify genomic features showing significant changes between two groups along pseudotime. The function models methylation changes and compares the fitted curves for each group, calculating the integral of the differences between the curves.
dmTwoGroups(Dat_sce_g1, Dat_sce_g2, BPPARAM = MulticoreParam())
dmTwoGroups(Dat_sce_g1, Dat_sce_g2, BPPARAM = MulticoreParam())
Dat_sce_g1 |
A
|
Dat_sce_g2 |
A
|
BPPARAM |
A |
A named numeric vector where each value corresponds to a genomic feature (e.g., a gene). The values represent the integral of the differences between the fitted curves of scDNA methylation levels for the two groups. The vector is sorted in descending order, with larger values indicating more drastic differences between the groups.
library(SingleCellExperiment) data_g1 <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) data_g2 <- readRDS(system.file("extdata", "group2_sampleData_sce.rds", package = "mist")) Dat_sce_g1 <- estiParam( Dat_sce = data_g1, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) Dat_sce_g2 <- estiParam( Dat_sce = data_g2, Dat_name = "Methy_level_group2", ptime_name = "pseudotime" ) # Run differential methylation analysis dm_results <- dmTwoGroups( Dat_sce_g1 = Dat_sce_g1, Dat_sce_g2 = Dat_sce_g2 )
library(SingleCellExperiment) data_g1 <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) data_g2 <- readRDS(system.file("extdata", "group2_sampleData_sce.rds", package = "mist")) Dat_sce_g1 <- estiParam( Dat_sce = data_g1, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) Dat_sce_g2 <- estiParam( Dat_sce = data_g2, Dat_name = "Methy_level_group2", ptime_name = "pseudotime" ) # Run differential methylation analysis dm_results <- dmTwoGroups( Dat_sce_g1 = Dat_sce_g1, Dat_sce_g2 = Dat_sce_g2 )
This function performs the Gibbs sampling procedure based on hierarchical Bayesian modeling to produce the parameters required for differential methylation analysis.
estiParam( Dat_sce, Dat_name, ptime_name, BPPARAM = MulticoreParam(), verbose = TRUE )
estiParam( Dat_sce, Dat_name, ptime_name, BPPARAM = MulticoreParam(), verbose = TRUE )
Dat_sce |
A |
Dat_name |
A character string specifying the name of the assay to extract the methylation level data. |
ptime_name |
A character string specifying the name of the column in |
BPPARAM |
A |
verbose |
A logical value indicating whether to print progress messages to the console.
Defaults to |
The updated sce object with A numeric matrix of estimated parameters for all genomic features in the rowData, including:
to
: Estimated coefficients for the polynomial of degree 4.
to
: Estimated variances for each stage along the pseudotime.
library(SingleCellExperiment) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" )
library(SingleCellExperiment) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" )
Generates a scatter plot of methylation levels vs. pseudotime for a specific gene, overlaid with the fitted curve based on the estimated coefficients.
plotGene(Dat_sce, Dat_name, ptime_name, gene_name)
plotGene(Dat_sce, Dat_name, ptime_name, gene_name)
Dat_sce |
The updated sce object with A numeric matrix of estimated parameters for all genomic features in the rowData, including:
|
Dat_name |
A character string specifying the name of the assay to extract the methylation data. |
ptime_name |
A character string specifying the name of the colData to extract the pseudotime vector. |
gene_name |
A character string specifying the gene name to plot. |
A ggplot2 scatter plot with an overlayed fitted curve.
library(SingleCellExperiment) library(ggplot2) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) plotGene(Dat_sce = Dat_sce_new, Dat_name = "Methy_level_group1", ptime_name = "pseudotime", gene_name = "ENSMUSG00000000037")
library(SingleCellExperiment) library(ggplot2) data <- readRDS(system.file("extdata", "group1_sampleData_sce.rds", package = "mist")) Dat_sce_new <- estiParam( Dat_sce = data, Dat_name = "Methy_level_group1", ptime_name = "pseudotime" ) plotGene(Dat_sce = Dat_sce_new, Dat_name = "Methy_level_group1", ptime_name = "pseudotime", gene_name = "ENSMUSG00000000037")
A small single-cell dataset for testing and illustration purposes, derived from the GEO dataset GSE121708. The data contains DNA methylation information from mouse gastrulation experiments at the single-cell level.
A SingleCellExperiment
object containing:
Metadata about cells and genes.
DNA methylation levels across a subset of genes and cells.
Data Processing Steps:
Annotation: Data was annotated using the mm10 reference genome.
Data Cleaning: Genes expressed in less than 10% of cells were removed, resulting in a scDNAm dataset containing 986 cells and 18,220 genes.
Pseudotime Inference: Pseudotime was inferred using the default settings in Monocle3.
Subsampling:
For illustration purposes (e.g., vignette generation), a random sample of 5 genes was selected to create group1_sampleData_sce.rds
.
The data originates from the GEO dataset GSE121708, as described in:
Argelaguet R, Clark SJ, Mohammed H, Stapel LC, et al. Multi-omics profiling of mouse gastrulation at single-cell resolution. Nature 2019 Dec;576(7787):487-491. PMID: 31827285.
Kapourani CA, Argelaguet R, Sanguinetti G, Vallejos CA. scMET: Bayesian modeling of DNA methylation heterogeneity at single-cell resolution. Genome Biol 2021 Apr 20;22(1):114. PMID: 33879195.
file_path <- system.file("extdata", "group1_sampleData_sce.rds", package = "mist") Dat_sce <- readRDS(file_path)
file_path <- system.file("extdata", "group1_sampleData_sce.rds", package = "mist") Dat_sce <- readRDS(file_path)
A small single-cell dataset for testing and illustration purposes, derived from the GEO dataset GSE121708. The data contains DNA methylation information from mouse gastrulation experiments at the single-cell level.
A SingleCellExperiment
object containing:
Metadata about cells and genes.
DNA methylation levels across a subset of genes and cells.
Data Processing Steps:
Annotation: Data was annotated using the mm10 reference genome.
Data Cleaning: Genes expressed in less than 10% of cells were removed, resulting in a scDNAm dataset containing 986 cells and 18,220 genes.
Pseudotime Inference: Pseudotime was inferred using the default settings in Monocle3.
Subsampling:
For illustration purposes (e.g., vignette generation), a random sample of 5 genes was selected to create group2_sampleData_sce.rds
.
The data originates from the GEO dataset GSE121708, as described in:
Argelaguet R, Clark SJ, Mohammed H, Stapel LC, et al. Multi-omics profiling of mouse gastrulation at single-cell resolution. Nature 2019 Dec;576(7787):487-491. PMID: 31827285.
Kapourani CA, Argelaguet R, Sanguinetti G, Vallejos CA. scMET: Bayesian modeling of DNA methylation heterogeneity at single-cell resolution. Genome Biol 2021 Apr 20;22(1):114. PMID: 33879195.
file_path <- system.file("extdata", "group2_sampleData_sce.rds", package = "mist") Dat_sce <- readRDS(file_path)
file_path <- system.file("extdata", "group2_sampleData_sce.rds", package = "mist") Dat_sce <- readRDS(file_path)