View source: R/makeVennDiagram.R
makeVennDiagram | R Documentation |
Make Venn Diagram from two or more peak ranges, Also calculate p-value to determine whether those peaks overlap significantly.
makeVennDiagram(
Peaks,
NameOfPeaks,
maxgap = -1L,
minoverlap = 0L,
totalTest,
by = c("region", "feature", "base"),
ignore.strand = TRUE,
connectedPeaks = c("min", "merge", "keepAll", "keepFirstListConsistent"),
method = c("hyperG", "permutation"),
TxDb,
plot = TRUE,
...
)
Peaks |
A list of peaks in GRanges format: See example below. |
NameOfPeaks |
Character vector to specify the name of Peaks, e.g., c("TF1", "TF2"). This will be used as label in the Venn Diagram. |
maxgap , minoverlap |
Used in the internal call to
|
totalTest |
Numeric value to specify the total number of tests performed to obtain the list of peaks. It should be much larger than the number of peaks in the largest peak set. |
by |
"region", "feature" or "base", default = "region". "feature" means using feature field in the GRanges for calculating overlap, "region" means using chromosome range for calculating overlap, and "base" means calculating overlap in nucleotide level. |
ignore.strand |
Logical: when set to TRUE, the strand information is ignored in the overlap calculations. |
connectedPeaks |
If multiple peaks involved in overlapping in several groups, set it to "merge" will count it as only 1, while set it to "min" will count it as the minimal involved peaks in any connected peak group. "keepAll" will show all the orginal counts for each list while the final counts will be same as "min". "keepFirstListConsistent" will keep the counts consistent with first list. |
method |
method to be used for p value calculation. hyperG means hypergeometric test and permutation means peakPermTest. |
TxDb |
An object of TxDb. |
plot |
logical. If TRUE (default), a venn diagram is plotted. |
... |
Additional arguments to be passed to venn.diagram. |
For customized graph options, please see venn.diagram in VennDiagram package.
A p.value is calculated by hypergeometric test or permutation test to determine whether the overlaps of peaks or features are significant.
Lihua Julie Zhu, Jianhong Ou
findOverlapsOfPeaks, venn.diagram, peakPermTest
if (interactive()){
peaks1 <- GRanges(seqnames=c("1", "2", "3"),
IRanges(start=c(967654, 2010897, 2496704),
end=c(967754, 2010997, 2496804),
names=c("Site1", "Site2", "Site3")),
strand="+",
feature=c("a","b","f"))
peaks2 = GRanges(seqnames=c("1", "2", "3", "1", "2"),
IRanges(start = c(967659, 2010898,2496700,
3075866,3123260),
end = c(967869, 2011108, 2496920,
3076166, 3123470),
names = c("t1", "t2", "t3", "t4", "t5")),
strand = c("+", "+", "-", "-", "+"),
feature=c("a","b","c","d","a"))
makeVennDiagram(list(peaks1, peaks2), NameOfPeaks=c("TF1", "TF2"),
totalTest=100,scaled=FALSE, euler.d=FALSE,
fill=c("#009E73", "#F0E442"), # circle fill color
col=c("#D55E00", "#0072B2"), #circle border color
cat.col=c("#D55E00", "#0072B2"))
makeVennDiagram(list(peaks1, peaks2), NameOfPeaks=c("TF1", "TF2"),
totalTest=100,
fill=c("#009E73", "#F0E442"), # circle fill color
col=c("#D55E00", "#0072B2"), #circle border color
cat.col=c("#D55E00", "#0072B2"))
###### 4-way diagram using annotated feature instead of chromosome ranges
makeVennDiagram(list(peaks1, peaks2, peaks1, peaks2),
NameOfPeaks=c("TF1", "TF2","TF3", "TF4"),
totalTest=100, by="feature",
main = "Venn Diagram for 4 peak lists",
fill=c(1,2,3,4))
}
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.