mergeDMRsIteratively: Merge DMRs iteratively

Description Usage Arguments Value Author(s) See Also Examples

View source: R/computeDMRs.R

Description

This function takes a list of DMRs and attempts to merge DMRs while keeping the new DMRs statistically significant.

Usage

1
2
3
4
mergeDMRsIteratively(DMRs, minGap, respectSigns = TRUE, methylationData1,
  methylationData2, context = "CG", minProportionDifference = 0.4,
  minReadsPerCytosine = 4, pValueThreshold = 0.01, test = "fisher",
  alternative = "two.sided", cores = 1)

Arguments

DMRs

the list of DMRs as a GRanges object; e.g. see computeDMRs

minGap

DMRs separated by a gap of at least minGap are not merged.

respectSigns

logical value indicating whether to respect the sign when joining DMRs.

methylationData1

the methylation data in condition 1 (see methylationDataList).

methylationData2

the methylation data in condition 2 (see methylationDataList).

context

the context in which the DMRs are computed ("CG", "CHG" or "CHH").

minProportionDifference

two adjacent DMRs are merged only if the difference in methylation proportion of the new DMR is higher than minProportionDifference.

minReadsPerCytosine

two adjacent DMRs are merged only if the number of reads per cytosine of the new DMR is higher than minReadsPerCytosine.

pValueThreshold

two adjacent DMRs are merged only if the p-value of the new DMR (see test below) is lower than pValueThreshold. Note that we adjust the p-values using the Benjamini and Hochberg's method to control the false discovery rate.

test

the statistical test used to call DMRs ("fisher" for Fisher's exact test or "score" for Score test).

alternative

indicates the alternative hypothesis and must be one of "two.sided", "greater" or "less".

cores

the number of cores used to compute the DMRs.

Value

the reduced list of DMRs as a GRanges object; e.g. see computeDMRs

Author(s)

Nicolae Radu Zabet

See Also

filterDMRs, computeDMRs, analyseReadsInsideRegionsForCondition and DMRsNoiseFilterCG

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
# load the methylation data
data(methylationDataList)

#load the DMRs in CG context they were computed with minGap = 200
data(DMRsNoiseFilterCG)


#merge the DMRs
DMRsNoiseFilterCGLarger <- mergeDMRsIteratively(DMRsNoiseFilterCG[1:100],
                           minGap = 500, respectSigns = TRUE,
                           methylationDataList[["WT"]],
                           methylationDataList[["met1-3"]],
                           context = "CG", minProportionDifference=0.4,
                           minReadsPerCytosine = 1, pValueThreshold=0.01,
                           test="score",alternative = "two.sided")


## Not run: 
#set genomic coordinates where to compute DMRs
regions <- GRanges(seqnames = Rle("Chr3"), ranges = IRanges(1,1E5))

# compute DMRs and remove gaps smaller than 200 bp
DMRsNoiseFilterCG200 <- computeDMRs(methylationDataList[["WT"]],
                       methylationDataList[["met1-3"]], regions = regions,
                       context = "CG", method = "noise_filter",
                       windowSize = 100, kernelFunction = "triangular",
                       test = "score", pValueThreshold = 0.01,
                       minCytosinesCount = 1, minProportionDifference = 0.4,
                       minGap = 200, minSize = 0, minReadsPerCytosine = 1,
                       cores = 1)

DMRsNoiseFilterCG0 <- computeDMRs(methylationDataList[["WT"]],
                       methylationDataList[["met1-3"]], regions = regions,
                       context = "CG", method = "noise_filter",
                       windowSize = 100, kernelFunction = "triangular",
                       test = "score", pValueThreshold = 0.01,
                       minCytosinesCount = 1, minProportionDifference = 0.4,
                       minGap = 0, minSize = 0, minReadsPerCytosine = 1,
                       cores = 1)
DMRsNoiseFilterCG0Merged200 <- mergeDMRsIteratively(DMRsNoiseFilterCG0,
                             minGap = 200, respectSigns = TRUE,
                             methylationDataList[["WT"]],
                             methylationDataList[["met1-3"]],
                             context = "CG", minProportionDifference=0.4,
                             minReadsPerCytosine = 1, pValueThreshold=0.01,
                             test="score",alternative = "two.sided")

#check that all newley computed DMRs are identical
print(all(DMRsNoiseFilterCG200 == DMRsNoiseFilterCG0Merged200))


## End(Not run)

DMRcaller documentation built on Nov. 8, 2020, 5:26 p.m.