Package 'DMCHMM'

Title: Differentially Methylated CpG using Hidden Markov Model
Description: A pipeline for identifying differentially methylated CpG sites using Hidden Markov Model in bisulfite sequencing data. DNA methylation studies have enabled researchers to understand methylation patterns and their regulatory roles in biological processes and disease. However, only a limited number of statistical approaches have been developed to provide formal quantitative analysis. Specifically, a few available methods do identify differentially methylated CpG (DMC) sites or regions (DMR), but they suffer from limitations that arise mostly due to challenges inherent in bisulfite sequencing data. These challenges include: (1) that read-depths vary considerably among genomic positions and are often low; (2) both methylation and autocorrelation patterns change as regions change; and (3) CpG sites are distributed unevenly. Furthermore, there are several methodological limitations: almost none of these tools is capable of comparing multiple groups and/or working with missing values, and only a few allow continuous or multiple covariates. The last of these is of great interest among researchers, as the goal is often to find which regions of the genome are associated with several exposures and traits. To tackle these issues, we have developed an efficient DMC identification method based on Hidden Markov Models (HMMs) called “DMCHMM” which is a three-step approach (model selection, prediction, testing) aiming to address the aforementioned drawbacks.
Authors: Farhad Shokoohi
Maintainer: Farhad Shokoohi <[email protected]>
License: GPL-3
Version: 1.27.0
Built: 2024-06-30 05:28:51 UTC
Source: https://github.com/bioc/DMCHMM

Help Index


Differentially Methylated CpG using Hidden Markov Model

Description

DMCHMM is a novel profiling tool for identifying differentially methylated CpG sites using Hidden Markov Model in bisulfite sequencing data.

DMCHMM methods

cBSData, cBSDMCs, methHMEM, methHMMCMC, findDMCs, qqDMCs, manhattanDMCs, readBismark, writeBED.

DMCHMM objects

BSData-class, BSDMCs-class


BSData object

Description

The BSData object is an S4 class that represents BS-Seq Data.

Arguments

methReads

The matrix methReads contains the number of methylated reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

totalReads

The matrix totalReads contains the number of reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

Value

A BSData-class object

Slots

methReads

An integer matrix

totalReads

An integer matrix

Author(s)

Farhad Shokoohi <[email protected]>

See Also

SummarizedExperiment objects.

Examples

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)
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<-cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ1

BSDMCs object

Description

The BSDMCs object is an S4 class that represents differentially methylated CpG sites (DMCs) in BS-Seq Data.

Arguments

methReads

The matrix methReads contains the number of methylated reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

totalReads

The matrix totalReads contains the number of reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methLevels

The matrix methLevels contains the predicted methylation level spanning a CpG-site using Hidden Markov model. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methStates

The matrix methStates contains the state of methylation obtained from Hidden Markov model spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData. The value of state is stored in metadata, named Beta.

methVars

The matrix methVars contains the variances of the corresponding methLevels obtianed from MCMC.

Value

A BSDMCs-class object

Slots

methReads

An integer matrix

totalReads

An integer matrix

methLevels

A numeric matrix

methStates

An integer matrix

methVars

A double matrix

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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 <- cBSDMCs(rowRanges=r1,methReads=methc,totalReads=metht,
methLevels=methl,methStates=meths,methVars=methv,colData=cd1)
OBJ2

cBSData method

Description

Creates a BSData-class object

Usage

cBSData(
  methReads,
  totalReads,
  rowRanges,
  colData = DataFrame(row.names = colnames(methReads)),
  metadata = list(),
  ...
)

## S4 method for signature 'matrix,matrix,GRanges'
cBSData(
  methReads,
  totalReads,
  rowRanges,
  colData = DataFrame(row.names = colnames(methReads)),
  metadata = list(),
  ...
)

Arguments

methReads

The matrix methReads contains the number of methylated reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

totalReads

The matrix totalReads contains the number of reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

rowRanges

A GRanges or GRangesList object describing the ranges of interest. Names, if present, become the row names of the SummarizedExperiment object. The length of the GRanges or GRangesList must equal the number of rows of the matrices in assays. If rowRanges is missing, a SummarizedExperiment instance is returned.

colData

Object of class "DataFrame" containing information on variable values of the samples

metadata

An optional list of arbitrary content describing the overall experiment

...

other possible parameters

Details

The rows of a BSData 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.

Value

A BSData-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ1

cBSDMCs method

Description

Creates a BSDMCs-class object

Usage

cBSDMCs(
  methReads,
  totalReads,
  methLevels,
  methStates,
  methVars,
  rowRanges,
  colData = DataFrame(row.names = colnames(methReads)),
  metadata = list(),
  ...
)

## S4 method for signature 'matrix,matrix,matrix,matrix,matrix,GRanges'
cBSDMCs(
  methReads,
  totalReads,
  methLevels,
  methStates,
  methVars,
  rowRanges,
  colData = DataFrame(row.names = colnames(methReads)),
  metadata = list(),
  ...
)

Arguments

methReads

The matrix methReads contains the number of methylated reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

totalReads

The matrix totalReads contains the number of reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methLevels

The matrix methLevels contains the predicted methylation level spanning a CpG-site using Hidden Markov model. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methStates

The matrix methStates contains the state of methylation obtained from Hidden Markov model spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData. The value of state is stored in metadata, named Beta.

methVars

The matrix methVars contains the variances of the corresponding methLevels obtianed from MCMC.

rowRanges

A GRanges or GRangesList object describing the ranges of interest. Names, if present, become the row names of the SummarizedExperiment object. The length of the GRanges or GRangesList must equal the number of rows of the matrices in assays. If rowRanges is missing, a SummarizedExperiment instance is returned.

colData

Object of class "DataFrame" containing information on variable values of the samples

metadata

An optional list of arbitrary content describing the overall experiment

...

other possible parameters

Details

The rows of a BSDMCs 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.

Value

A BSDMCs-class

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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 <- cBSDMCs(rowRanges=r1,methReads=methc,totalReads=metht,
methLevels=methl,methStates=meths,methVars=methv,colData=cd1)
OBJ2

combine method

Description

combine two BSData-class or two BSDMCs-class

Usage

combine(obj1, obj2)

## S4 method for signature 'BSData,BSData'
combine(obj1, obj2)

## S4 method for signature 'BSDMCs,BSDMCs'
combine(obj1, obj2)

Arguments

obj1

A BSData-class or BSDMCs-class

obj2

A BSData-class or BSDMCs-class

Value

A BSData-class or BSDMCs-class

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc[,1:nc],totalReads=metht[,1:nc],
colData=cd1)
cd2 <- DataFrame(Group=rep('G2',each=nc),row.names=LETTERS[nc+1:nc])
OBJ2 <- cBSData(rowRanges=r1,methReads=methc[,nc+1:nc],totalReads=
metht[,nc+1:nc],colData=cd2)
OBJ3 <- combine(OBJ1, OBJ2)
OBJ3

data

Description

A part of BS-Seq data for three cell type: WGBS data were derived from whole blood collected on a cohort of healthy individuals from Sweden. Cell lines were separated into T-cells (19 samples), monocytes (13 samples) and B-cells (8 samples). Sequencing was performed on the Illumina HiSeq2000/2500 system for each of the 40 samples, separately. For illustration only 3 samples each containg 30,440 CpG sites around BLK gene are provided here. The whole data are analyzed in the cited paper.

Format

BED files

Details

The data is part of whole blood from Sweeden.

Author(s)

Farhad Shokoohi <[email protected]>

Source

Genomic Quebec


findDMCs method

Description

finds the DMCs after smoothing using HMM

Usage

findDMCs(
  object,
  formula,
  FDRthreshold,
  Methylthreshold,
  mc.cores,
  windowsize,
  weightfunction
)

## S4 method for signature 'BSDMCs'
findDMCs(
  object,
  formula,
  FDRthreshold,
  Methylthreshold,
  mc.cores,
  windowsize,
  weightfunction
)

Arguments

object

A BSData-class or BSDMCs-class object

formula

A formula

FDRthreshold

A numeric value

Methylthreshold

A positive numeric value; the default is 0.001

mc.cores

An integer greater than 0

windowsize

An integer value for partitioning data into windows of size windowsize.

weightfunction

A function to create weights using variance obtained form the MCMC algorithm

Value

BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ2 <- methHMEM(OBJ1, MaxK=2, mc.cores=2)
OBJ3 <- methHMMCMC(OBJ2, mc.cores=2)
OBJ4 <- findDMCs(OBJ3, mc.cores=2)
head(metadata(OBJ4)$DMCHMM)

manhattanDMCs method

Description

Creates a Manhattan plot based on the p-values obtained from findDMCs method

Usage

manhattanDMCs(
  object,
  col,
  chrlabs,
  suggestiveline,
  genomewideline,
  highlight,
  logp,
  annotatePval,
  annotateTop,
  ...
)

## S4 method for signature 'BSDMCs'
manhattanDMCs(
  object,
  col,
  chrlabs,
  suggestiveline,
  genomewideline,
  highlight,
  logp,
  annotatePval,
  annotateTop,
  ...
)

Arguments

object

A BSData-class or BSDMCs-class object

col

A character vector indicating which colors to alternate.

chrlabs

A character vector equal to the number of chromosomes specifying the chromosome labels (e.g., c(1:22, "X", "Y", "MT")).

suggestiveline

Where to draw a "suggestive" line. Default -log10(1e-5). Set to FALSE to disable.

genomewideline

Where to draw a "genome-wide sigificant" line. Default -log10(5e-8). Set to FALSE to disable.

highlight

A character vector of SNPs in your dataset to highlight. These SNPs should all be in your dataset.

logp

If TRUE, the -log10 of the p-value is plotted. It isn't very useful to plot raw p-values, but plotting the raw value could be useful for other genome-wide plots, for example, peak heights, bayes factors, test statistics, other "scores," etc.

annotatePval

If set, SNPs below this p-value will be annotated on the plot.

annotateTop

If TRUE, only annotates the top hit on each chromosome that is below the annotatePval threshold.

...

other possible parameters

Value

A Manhattan plot

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ2 <- methHMEM(OBJ1, MaxK=2, mc.cores=2)
OBJ3 <- methHMMCMC(OBJ2, mc.cores=2)
OBJ4 <- findDMCs(OBJ3, mc.cores=2)
manhattanDMCs(OBJ4)

methHMEM method

Description

Estimates the HMM methylation paths and the HMM order for each sample using the EM algorithm

Usage

methHMEM(object, MaxK, MaxEmiter, epsEM, useweight, mc.cores)

## S4 method for signature 'BSData'
methHMEM(object, MaxK, MaxEmiter, epsEM, useweight, mc.cores)

Arguments

object

A BSData-class or BSDMCs-class object

MaxK

An integer value

MaxEmiter

An integer value

epsEM

A positive numeric value

useweight

A logical value

mc.cores

An integer greater than 0

Value

BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ2 <- methHMEM(OBJ1, MaxK=2, mc.cores=2)
OBJ2

methHMMCMC method

Description

Estimates the HMM methylation paths and the HMM order for each sample using the MCMC algorithm

Usage

methHMMCMC(object, useweight, nburn, nthin, nsamp, mc.cores)

## S4 method for signature 'BSDMCs'
methHMMCMC(object, useweight, nburn, nthin, nsamp, mc.cores)

Arguments

object

A BSData-class or BSDMCs-class object

useweight

A logical value

nburn

An integer value

nthin

An integer value

nsamp

An integer value

mc.cores

An integer greater than 0

Value

BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ2 <- methHMEM(OBJ1, MaxK=2, mc.cores=2)
OBJ3 <- methHMMCMC(OBJ2, mc.cores=2)
OBJ3

methLevels method

Description

Returns methLevels stored in BSDMCs-class

Assigns methLevels to BSDMCs-class

Usage

methLevels(object)

methLevels(object) <- value

## S4 method for signature 'BSDMCs'
methLevels(object)

## S4 replacement method for signature 'BSDMCs,matrix'
methLevels(object) <- value

Arguments

object

A BSData-class or BSDMCs-class object

value

An integer matrix

Value

A matrix

A BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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 <- cBSDMCs(rowRanges=r1,methReads=methc,totalReads=metht,
methLevels=methl,methStates=meths,methVars=methv,colData=cd1)
methLevels(OBJ2)
methLevels(OBJ2) <- methl

methReads method

Description

Returns methReads stored in BSData-class

Assigns methReads to BSData-class

Returns methReads stored in BSDMCs-class

Assigns methReads to BSDMCs-class

Usage

methReads(object)

methReads(object) <- value

methReads(object)

methReads(object) <- value

## S4 method for signature 'BSData'
methReads(object)

## S4 replacement method for signature 'BSData,matrix'
methReads(object) <- value

## S4 method for signature 'BSDMCs'
methReads(object)

## S4 replacement method for signature 'BSDMCs,matrix'
methReads(object) <- value

Arguments

object

A BSData-class or BSDMCs-class object

value

An integer matrix

Value

A matrix

A BSData-class object

A matrix

A BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
methReads(OBJ1)
methReads(OBJ1) <- methc

methStates method

Description

Returns methStates stored in BSDMCs-class

Assigns methStates to BSDMCs-class

Usage

methStates(object)

methStates(object) <- value

## S4 method for signature 'BSDMCs'
methStates(object)

## S4 replacement method for signature 'BSDMCs,matrix'
methStates(object) <- value

Arguments

object

A BSData-class or BSDMCs-class object

value

An integer matrix

Value

A matrix

A BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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 <- cBSDMCs(rowRanges=r1,methReads=methc,totalReads=metht,
methLevels=methl,methStates=meths,methVars=methv,colData=cd1)
methStates(OBJ2)
methStates(OBJ2)<- meths

methVars method

Description

Returns methVars stored in BSDMCs-class

Assigns methVars to BSDMCs-class

Usage

methVars(object)

methVars(object) <- value

## S4 method for signature 'BSDMCs'
methVars(object)

## S4 replacement method for signature 'BSDMCs,matrix'
methVars(object) <- value

Arguments

object

A BSData-class or BSDMCs-class object

value

An integer matrix

Value

A matrix

A BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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 <- cBSDMCs(rowRanges=r1,methReads=methc,totalReads=metht,
methLevels=methl,methStates=meths,methVars=methv,colData=cd1)
methVars(OBJ2)
methVars(OBJ2)<- meths

params

Description

parameters name and their descriptions

Arguments

methReads

The matrix methReads contains the number of methylated reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

totalReads

The matrix totalReads contains the number of reads spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methLevels

The matrix methLevels contains the predicted methylation level spanning a CpG-site using Hidden Markov model. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData.

methVars

The matrix methVars contains the variances of the corresponding methLevels obtianed from MCMC.

methStates

The matrix methStates contains the state of methylation obtained from Hidden Markov model spanning a CpG-site. The rows represent the CpG sites in rowRanges and the columns represent the samples in colData. The value of state is stored in metadata, named Beta.

rowRanges

A GRanges or GRangesList object describing the ranges of interest. Names, if present, become the row names of the SummarizedExperiment object. The length of the GRanges or GRangesList must equal the number of rows of the matrices in assays. If rowRanges is missing, a SummarizedExperiment instance is returned.

colData

Object of class "DataFrame" containing information on variable values of the samples

metadata

An optional list of arbitrary content describing the overall experiment

object

A BSData-class or BSDMCs-class object

value

An integer matrix

obj1

A BSData-class or BSDMCs-class

obj2

A BSData-class or BSDMCs-class

files

A character list

file

A character

name

A character list

MaxK

An integer value

MaxEmiter

An integer value

epsEM

A positive numeric value

useweight

A logical value

mc.cores

An integer greater than 0

nburn

An integer value

nthin

An integer value

nsamp

An integer value

formula

A formula

FDRthreshold

A numeric value

Methylthreshold

A positive numeric value; the default is 0.001

weightfunction

A function to create weights using variance obtained form the MCMC algorithm

...

other possible parameters

col

A character vector indicating which colors to alternate.

chrlabs

A character vector equal to the number of chromosomes specifying the chromosome labels (e.g., c(1:22, "X", "Y", "MT")).

suggestiveline

Where to draw a "suggestive" line. Default -log10(1e-5). Set to FALSE to disable.

genomewideline

Where to draw a "genome-wide sigificant" line. Default -log10(5e-8). Set to FALSE to disable.

highlight

A character vector of SNPs in your dataset to highlight. These SNPs should all be in your dataset.

logp

If TRUE, the -log10 of the p-value is plotted. It isn't very useful to plot raw p-values, but plotting the raw value could be useful for other genome-wide plots, for example, peak heights, bayes factors, test statistics, other "scores," etc.

annotatePval

If set, SNPs below this p-value will be annotated on the plot.

annotateTop

If TRUE, only annotates the top hit on each chromosome that is below the annotatePval threshold.

windowsize

An integer value for partitioning data into windows of size windowsize.

Author(s)

Farhad Shokoohi <[email protected]>


qqDMCs method

Description

Creates a Q-Q plot based on the p-values obtained from findDMCs method

Usage

qqDMCs(object, ...)

## S4 method for signature 'BSDMCs'
qqDMCs(object, ...)

Arguments

object

A BSData-class or BSDMCs-class object

...

other possible parameters

Value

A QQ plot

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
OBJ2 <- methHMEM(OBJ1, MaxK=2, mc.cores=2)
OBJ3 <- methHMMCMC(OBJ2, mc.cores=2)
OBJ4 <- findDMCs(OBJ3, mc.cores=2)
qqDMCs(OBJ4)

readBismark method

Description

reads BS-Seq data

Usage

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)

Arguments

files

A character list

colData

Object of class "DataFrame" containing information on variable values of the samples

mc.cores

An integer greater than 0

Value

A BSData-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

fn <- list.files(system.file('extdata',package = 'DMCHMM'))
fn.f <- list.files(system.file('extdata',package='DMCHMM'), full.names=TRUE)
OBJ <- readBismark(fn.f, fn, mc.cores=2)
cdOBJ <- DataFrame(Cell = factor(c('BC', 'TC','Mono'),
labels = c('BC', 'TC', 'Mono')), row.names = c('BCU1568','BCU173','BCU551'))
colData(OBJ) <- cdOBJ
OBJ

totalReads method

Description

Returns totalReads stored in BSData-class

Assigns totalReads to BSData-class

Returns totalReads stored in BSDMCs-class

Assigns totalReads to BSDMCs-class

Usage

totalReads(object)

totalReads(object) <- value

totalReads(object)

totalReads(object) <- value

## S4 method for signature 'BSData'
totalReads(object)

## S4 replacement method for signature 'BSData,matrix'
totalReads(object) <- value

## S4 method for signature 'BSDMCs'
totalReads(object)

## S4 replacement method for signature 'BSDMCs,matrix'
totalReads(object) <- value

Arguments

object

A BSData-class or BSDMCs-class object

value

An integer matrix

Value

A matrix

A BSData-class object

A matrix

A BSDMCs-class object

Author(s)

Farhad Shokoohi <[email protected]>

Examples

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)
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 <- cBSData(rowRanges=r1,methReads=methc,totalReads=metht,colData=cd1)
totalReads(OBJ1)
totalReads(OBJ1) <- metht

writeBED method

Description

write BS-Seq data to BED files

Usage

writeBED(object, name, file)

## S4 method for signature 'BSData,character,character'
writeBED(object, name, file)

## S4 method for signature 'BSData,character,missing'
writeBED(object, name)

## S4 method for signature 'BSData,missing,character'
writeBED(object, file)

## S4 method for signature 'BSData,missing,missing'
writeBED(object)

## S4 method for signature 'BSDMCs,character,character'
writeBED(object, name, file)

## S4 method for signature 'BSDMCs,character,missing'
writeBED(object, name)

## S4 method for signature 'BSDMCs,missing,character'
writeBED(object, file)

## S4 method for signature 'BSDMCs,missing,missing'
writeBED(object)

Arguments

object

A BSData-class or BSDMCs-class object

name

A character list

file

A character

Value

BED files

Author(s)

Farhad Shokoohi <[email protected]>