View source: R/dmpClustering.R
dmpClustering | R Documentation |
Given a 'pDMP' object carrying DMPs obtained in Methyl-IT
downstream analysis, function 'dmpClustering' build
clusters of DMPs, which can be further tested to identify differentially
methylated regions (DMRs) with countTest2
function.
dmpClustering(
dmps,
win.size = NULL,
step.size = NULL,
minNumDMPs = 1,
maxClustDist = NULL,
method = c("relaxed", "fixed.int"),
ignore.strand = TRUE,
verbose = FALSE
)
dmps |
An object from 'pDMP' class, which is returned
by |
win.size |
An integer. The size of the windows/intervals genomics.
Default: |
step.size |
Interval at which the regions/windows must be defined.
Default: |
minNumDMPs |
Minimum number of DMPs inside of each cluster.
Default: |
maxClustDist |
Clusters separated by a distance lesser than
'maxClustDist' positions are merged. Default: |
method |
Two different approaches are implemented to clustering DMPs:
|
ignore.strand |
Same as in
|
verbose |
if TRUE, prints the function log to stdout |
The number of DMPs reported in each cluster corresponds to the numbers of sites inside the cluster where DMPs were found in at least one of the samples (from control or from treatment). That is, dmpClustering is just a tool to locate regions with high density of DMPs from all the samples. It does not detect DMRs. It is assumed that only DMP coordinates are given in the 'dmps' object. That is, all the sites provided are considered in the computation.
A GRanges object carrying the coordinates of DMP clusters from all the samples and the number of DMPs on each of them.
Robersy Sanchez (https://github.com/genomaths).
## Creates a GRanges object carrying DMPs. Notice that only the DMP
## coordinates are needed.
gr <- GRanges(seqnames = Rle( c('chr1', 'chr2', 'chr3', 'chr4'),
c(5, 5, 5, 5)),
ranges = IRanges(start = 1:20, end = 1:20),
strand = rep(c('+', '-'), 10))
## Simple DMP clustering ignoring the DNA strand
dmpClustering(gr, win.size = 4, step.size = 4, minNumDMPs = 2,
method = "fixed.int")
## Now, the information on the DNA strand is included in the clustering
dmpClustering(dmps = gr, win.size = 4, step.size = 4, minNumDMPs = 2,
method = "fixed.int", ignore.strand = FALSE)
## Next, as before adding that clusters separated by a distance lesser than
## 'maxClustDist = 2' will be merged
dmpClustering(dmps = gr, win.size = 4, step.size = 4, minNumDMPs = 2,
method = "fixed.int", maxClustDist = 2, ignore.strand = FALSE)
## Finally, the relaxed approach. Notice that only two parameter values are
## needed
dmpClustering(gr, minNumDMPs = 2, maxClustDist = 2,
method = "relaxed", ignore.strand = FALSE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.