Title: | Differentially Methylated Cytosines via a Bayesian Functional Approach |
---|---|
Description: | DMCFB is a pipeline for identifying differentially methylated cytosines using a Bayesian functional regression model in bisulfite sequencing data. By using a functional regression data model, it tries to capture position-specific, group-specific and other covariates-specific methylation patterns as well as spatial correlation patterns and unknown underlying models of methylation data. It is robust and flexible with respect to the true underlying models and inclusion of any covariates, and the missing values are imputed using spatial correlation between positions and samples. A Bayesian approach is adopted for estimation and inference in the proposed method. |
Authors: | Farhad Shokoohi [aut, cre] |
Maintainer: | Farhad Shokoohi <[email protected]> |
License: | GPL-3 |
Version: | 1.21.0 |
Built: | 2024-11-29 05:26:09 UTC |
Source: | https://github.com/bioc/DMCFB |
DMCFB is a profiling tool for identifying differentially methylated cytosines using Functional Bayesian Model in bisulfite sequencing data.
DMCFB
methodsfindDMCFB
,
plotDMCFB
,
cBSDMC
,
readBismark
.
BSDMC
objectsThe BSDMC
object is an S4 class that represents
differentially methylated CpG sites (DMCs) in BS-Seq Data.
methReads |
The matrix |
totalReads |
The matrix |
methLevels |
The matrix |
A BSDMC-class
object
methReads
An integer matrix
totalReads
An integer matrix
methLevels
A numeric matrix
Farhad Shokoohi <[email protected]>
RangedSummarizedExperiment-class
GRanges-class
nr <- 500 nc <- 16 metht <- matrix(as.integer(runif(nr * nc, 0, nr)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) meths <- matrix(as.integer(runif(nr * nc, 0, 10)), nr) methl <- methc / metht methv <- matrix((runif(nr * nc, 0.1, 0.5)), nr) r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame(Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc]) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, methStates = meths, methVars = methv, colData = cd1 ) OBJ2
nr <- 500 nc <- 16 metht <- matrix(as.integer(runif(nr * nc, 0, nr)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) meths <- matrix(as.integer(runif(nr * nc, 0, 10)), nr) methl <- methc / metht methv <- matrix((runif(nr * nc, 0.1, 0.5)), nr) r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame(Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc]) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, methStates = meths, methVars = methv, colData = cd1 ) OBJ2
Creates a BSDMC-class
object
cBSDMC( methReads, totalReads, methLevels, rowRanges, colData = DataFrame(row.names = colnames(methReads)), metadata = list(), ... ) ## S4 method for signature 'matrix,matrix,matrix,GRanges' cBSDMC( methReads, totalReads, methLevels, rowRanges, colData = DataFrame(row.names = colnames(methReads)), metadata = list(), ... )
cBSDMC( methReads, totalReads, methLevels, rowRanges, colData = DataFrame(row.names = colnames(methReads)), metadata = list(), ... ) ## S4 method for signature 'matrix,matrix,matrix,GRanges' cBSDMC( methReads, totalReads, methLevels, rowRanges, colData = DataFrame(row.names = colnames(methReads)), metadata = list(), ... )
methReads |
The matrix |
totalReads |
The matrix |
methLevels |
The matrix |
rowRanges |
A |
colData |
Object of class |
metadata |
A |
... |
other possible parameters |
The rows of a BSDMC
object represent ranges (in genomic
coordinates) of interest. The ranges of interest are described by a
GRanges
or a GRangesList
object,
accessible using the rowRanges
function.
The GRanges
and GRangesList
classes
contains sequence (e.g., chromosome) name, genomic coordinates, and strand
information. Each range can be annotated with additional data; this data
might be used to describe the range or to
summarize results (e.g., statistics of differential abundance) relevant to
the range. Rows may or may not have row names; they often will not.
Farhad Shokoohi <[email protected]>
set.seed(1980) nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) meths <- matrix(as.integer(runif(nr * nc, 0, 10)), nr) methl <- methc / metht methv <- matrix((runif(nr * nc, 0.1, 0.5)), nr) r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, methStates = meths, methVars = methv, colData = cd1 ) OBJ2
set.seed(1980) nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) meths <- matrix(as.integer(runif(nr * nc, 0, 10)), nr) methl <- methc / metht methv <- matrix((runif(nr * nc, 0.1, 0.5)), nr) r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, methStates = meths, methVars = methv, colData = cd1 ) OBJ2
combine two BSDMC-class
or
two BSDMC-class
combine(obj1, obj2) ## S4 method for signature 'BSDMC,BSDMC' combine(obj1, obj2)
combine(obj1, obj2) ## S4 method for signature 'BSDMC,BSDMC' combine(obj1, obj2)
obj1 |
|
obj2 |
A BSDMC-class
or BSDMC-class
Farhad Shokoohi <[email protected]>
set.seed(1980) nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc * 2, 0, nr)), nr) methc <- matrix( rbinom(n = nr * nc, c(metht), prob = runif(nr * nc * 2)), nr, nc * 2 ) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame(Group = rep("G1", each = nc), row.names = LETTERS[1:nc]) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc[, 1:nc], totalReads = metht[, 1:nc], methLevels = methl[, 1:nc], colData = cd1 ) cd2 <- DataFrame( Group = rep("G2", each = nc), row.names = LETTERS[nc + 1:nc] ) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc[, nc + 1:nc], totalReads = metht[, nc + 1:nc], methLevels = methl[, nc + 1:nc], colData = cd2 ) OBJ3 <- combine(OBJ1, OBJ2) OBJ3
set.seed(1980) nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc * 2, 0, nr)), nr) methc <- matrix( rbinom(n = nr * nc, c(metht), prob = runif(nr * nc * 2)), nr, nc * 2 ) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame(Group = rep("G1", each = nc), row.names = LETTERS[1:nc]) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc[, 1:nc], totalReads = metht[, 1:nc], methLevels = methl[, 1:nc], colData = cd1 ) cd2 <- DataFrame( Group = rep("G2", each = nc), row.names = LETTERS[nc + 1:nc] ) OBJ2 <- cBSDMC( rowRanges = r1, methReads = methc[, nc + 1:nc], totalReads = metht[, nc + 1:nc], methLevels = methl[, nc + 1:nc], colData = cd2 ) OBJ3 <- combine(OBJ1, OBJ2) OBJ3
DMC identification via Bayesian functional regression models
findDMCFB( object, bwa, bwb, nBurn, nMC, nThin, alpha, sdv, nCores, pSize, sfiles ) ## S4 method for signature 'BSDMC' findDMCFB( object, bwa, bwb, nBurn, nMC, nThin, alpha, sdv, nCores, pSize, sfiles )
findDMCFB( object, bwa, bwb, nBurn, nMC, nThin, alpha, sdv, nCores, pSize, sfiles ) ## S4 method for signature 'BSDMC' findDMCFB( object, bwa, bwb, nBurn, nMC, nThin, alpha, sdv, nCores, pSize, sfiles )
object |
A |
bwa |
An integer value specifying the band-width size of B-spline basis matrix for a natural cubic spline for the group-specific effects of the Bayesian functional regression model |
bwb |
An integer value specifying the band-width size of B-spline basis matrix for a natural cubic spline for the individual-specific effects of the Bayesian functional regression model |
nBurn |
An integer value specifying the number of burn-in samples |
nMC |
An integer value specifying the number of MCMC samples after burn-in |
nThin |
An integer value specifying the thining number in MCMC |
alpha |
A numeric value specifying the level of |
sdv |
An double value specifying the standard deviation of priors |
nCores |
An integer value specifying the number of machine cores for parallel computing |
pSize |
An integer value specifying the number of cytosines in a regrion to be used in a Bayesian functiona regression model for DMC detection |
sfiles |
A logical value indicating whether files to be saved or not. |
BSDMC-class
object
Farhad Shokoohi <[email protected]>
set.seed(1980) nr <- 1000 nc <- 4 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) OBJ2 <- findDMCFB(OBJ1, bwa = 10, bwb = 10, nBurn = 50, nMC = 50, nThin = 1, alpha = 0.05, nCores = 2, pSize = 500, sfiles = FALSE ) OBJ2
set.seed(1980) nr <- 1000 nc <- 4 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) OBJ2 <- findDMCFB(OBJ1, bwa = 10, bwb = 10, nBurn = 50, nMC = 50, nThin = 1, alpha = 0.05, nCores = 2, pSize = 500, sfiles = FALSE ) OBJ2
Returns methLevels
stored in BSDMC-class
Assigns methLevels
to BSDMC-class
methLevels(object) methLevels(object) <- value ## S4 method for signature 'BSDMC' methLevels(object) ## S4 replacement method for signature 'BSDMC,matrix' methLevels(object) <- value
methLevels(object) methLevels(object) <- value ## S4 method for signature 'BSDMC' methLevels(object) ## S4 replacement method for signature 'BSDMC,matrix' methLevels(object) <- value
object |
A |
value |
An integer matrix |
A matrix
A BSDMC-class
object
Farhad Shokoohi <[email protected]>
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) methLevels(OBJ1) methLevels(OBJ1) <- methl
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) methLevels(OBJ1) methLevels(OBJ1) <- methl
Returns methReads
stored in BSDMC-class
Assigns methReads
to BSDMC-class
methReads(object) methReads(object) <- value ## S4 method for signature 'BSDMC' methReads(object) ## S4 replacement method for signature 'BSDMC,matrix' methReads(object) <- value
methReads(object) methReads(object) <- value ## S4 method for signature 'BSDMC' methReads(object) ## S4 replacement method for signature 'BSDMC,matrix' methReads(object) <- value
object |
A |
value |
An integer matrix |
A matrix
A BSDMC-class
object
Farhad Shokoohi <[email protected]>
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) methReads(OBJ1) methReads(OBJ1) <- methc
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) methReads(OBJ1) methReads(OBJ1) <- methc
parameters name and their descriptions
methReads |
The matrix |
totalReads |
The matrix |
methLevels |
The matrix |
rowRanges |
A |
colData |
Object of class |
metadata |
A |
object |
A |
value |
An integer matrix |
name |
A character list |
obj1 |
|
obj2 |
|
files |
A character list |
file |
A character |
nCores |
An integer value specifying the number of machine cores for parallel computing |
mc.cores |
An integer greater than 0 |
pSize |
An integer value specifying the number of cytosines in a regrion to be used in a Bayesian functiona regression model for DMC detection |
bwa |
An integer value specifying the band-width size of B-spline basis matrix for a natural cubic spline for the group-specific effects of the Bayesian functional regression model |
bwb |
An integer value specifying the band-width size of B-spline basis matrix for a natural cubic spline for the individual-specific effects of the Bayesian functional regression model |
nBurn |
An integer value specifying the number of burn-in samples |
nThin |
An integer value specifying the thining number in MCMC |
nMC |
An integer value specifying the number of MCMC samples after burn-in |
sdv |
An double value specifying the standard deviation of priors |
alpha |
A numeric value specifying the level of |
col |
A character vector indicating which colors to alternate. |
sfiles |
A logical value indicating whether files to be saved or not. |
region |
An integer vector of length two specifying which subset of the object to be plotted |
nSplit |
A integer value specifying the number of subsets must be done for plotting the results of DMC identification |
parList |
A list specifying plots parameters, see |
... |
other possible parameters |
Farhad Shokoohi <[email protected]>
Plotting the results of DMC identifation stored in
a BSDMC-class
object
plotDMCFB(object, region, nSplit, parList) ## S4 method for signature 'BSDMC' plotDMCFB(object, region, nSplit, parList)
plotDMCFB(object, region, nSplit, parList) ## S4 method for signature 'BSDMC' plotDMCFB(object, region, nSplit, parList)
object |
A |
region |
An integer vector of length two specifying which subset of the object to be plotted |
nSplit |
A integer value specifying the number of subsets must be done for plotting the results of DMC identification |
parList |
A list specifying plots parameters, see |
Plot
Farhad Shokoohi <[email protected]>
set.seed(1980) nr <- 1000 nc <- 4 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) OBJ2 <- findDMCFB(OBJ1, bwa = 10, bwb = 10, nBurn = 50, nMC = 50, nThin = 1, alpha = 0.05, nCores = 2, pSize = 500, sfiles = FALSE ) plotDMCFB(OBJ2)
set.seed(1980) nr <- 1000 nc <- 4 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) OBJ2 <- findDMCFB(OBJ1, bwa = 10, bwb = 10, nBurn = 50, nMC = 50, nThin = 1, alpha = 0.05, nCores = 2, pSize = 500, sfiles = FALSE ) plotDMCFB(OBJ2)
reads BS-Seq data
readBismark(files, colData, mc.cores) ## S4 method for signature 'character,DataFrame,numeric' readBismark(files, colData, mc.cores) ## S4 method for signature 'character,data.frame,numeric' readBismark(files, colData, mc.cores) ## S4 method for signature 'character,character,numeric' readBismark(files, colData, mc.cores)
readBismark(files, colData, mc.cores) ## S4 method for signature 'character,DataFrame,numeric' readBismark(files, colData, mc.cores) ## S4 method for signature 'character,data.frame,numeric' readBismark(files, colData, mc.cores) ## S4 method for signature 'character,character,numeric' readBismark(files, colData, mc.cores)
files |
A character list |
colData |
Object of class |
mc.cores |
An integer greater than 0 |
A BSDMC-class
object
Farhad Shokoohi <[email protected]>
fn <- list.files(system.file("extdata", package = "DMCFB")) fn.f <- list.files(system.file("extdata", package = "DMCFB"), full.names = TRUE ) OBJ <- readBismark(fn.f, fn, mc.cores=1) cdOBJ <- DataFrame(Cell = factor(c("BC", "TC", "Mono"), labels = c("BC", "TC", "Mono") ), row.names = c("BCU1568", "BCU173", "BCU551")) colData(OBJ) <- cdOBJ OBJ
fn <- list.files(system.file("extdata", package = "DMCFB")) fn.f <- list.files(system.file("extdata", package = "DMCFB"), full.names = TRUE ) OBJ <- readBismark(fn.f, fn, mc.cores=1) cdOBJ <- DataFrame(Cell = factor(c("BC", "TC", "Mono"), labels = c("BC", "TC", "Mono") ), row.names = c("BCU1568", "BCU173", "BCU551")) colData(OBJ) <- cdOBJ OBJ
Returns totalReads
stored in BSDMC-class
Assigns totalReads
to BSDMC-class
totalReads(object) totalReads(object) <- value ## S4 method for signature 'BSDMC' totalReads(object) ## S4 replacement method for signature 'BSDMC,matrix' totalReads(object) <- value
totalReads(object) totalReads(object) <- value ## S4 method for signature 'BSDMC' totalReads(object) ## S4 replacement method for signature 'BSDMC,matrix' totalReads(object) <- value
object |
A |
value |
An integer matrix |
A matrix
A BSDMC-class
object
Farhad Shokoohi <[email protected]>
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) totalReads(OBJ1) totalReads(OBJ1) <- metht
nr <- 150 nc <- 8 metht <- matrix(as.integer(runif(nr * nc, 0, 100)), nr) methc <- matrix(rbinom(n = nr * nc, c(metht), prob = runif(nr * nc)), nr, nc) methl <- methc / metht r1 <- GRanges(rep("chr1", nr), IRanges(1:nr, width = 1), strand = "*") names(r1) <- 1:nr cd1 <- DataFrame( Group = rep(c("G1", "G2"), each = nc / 2), row.names = LETTERS[1:nc] ) OBJ1 <- cBSDMC( rowRanges = r1, methReads = methc, totalReads = metht, methLevels = methl, colData = cd1 ) totalReads(OBJ1) totalReads(OBJ1) <- metht