RPeakComp: Find the overlap or differential peaks between two samples.

Description Usage Arguments Value Author(s) See Also Examples

Description

This function compares two peak file and report overlap or differential peaks according to the parameter "operation".

Usage

 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
atacpeakComp(
  atacProcPeak1,
  atacProcPeak2,
  bedInput1 = NULL,
  bedInput2 = NULL,
  bedOutput = NULL,
  olap.rate = 0.2,
  ...
)

## S4 method for signature 'ATACProc'
atacpeakComp(
  atacProcPeak1,
  atacProcPeak2,
  bedInput1 = NULL,
  bedInput2 = NULL,
  bedOutput = NULL,
  olap.rate = 0.2,
  ...
)

peakcomp(
  bedInput1 = NULL,
  bedInput2 = NULL,
  bedOutput = NULL,
  olap.rate = 0.2,
  ...
)

Arguments

atacProcPeak1

ATACProc-class object scalar. It has to be the return value of upstream process: atacPeakCalling.

atacProcPeak2

ATACProc-class object scalar. It has to be the return value of upstream process: atacPeakCalling.

bedInput1

Character scalar. Input peak file path. UCSC bed file is recommented. Other file should be able to import as GRanges objects through import.

bedInput2

Character scalar. Input peak file path. UCSC bed file is recommented. Other file should be able to import as GRanges objects through import.

bedOutput

The output file path for overlap peaks.

olap.rate

Overlap rate, if the overlap region between 2 peak is more than this rate of the short peak, these two peak are considered to be overlap and will be merged to a bigger peak. Default: 0.2. NOTICE: multi-peak will be merged together!

...

Additional arguments, currently unused.

Value

An invisible ATACProc-class object scalar for downstream analysis.

Author(s)

Wei Zhang

See Also

atacPeakCalling

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
library(R.utils)
p1bz <- system.file("extdata", "Example_peak1.bed.bz2", package="esATAC")
p2bz <- system.file("extdata", "Example_peak2.bed.bz2", package="esATAC")
## Not run: 
peak1_path <- as.vector(bunzip2(filename = p1bz,
destname = file.path(getwd(), "Example_peak1.bed"),
ext="bz2", FUN=bzfile, overwrite=TRUE , remove = FALSE))
peak2_path <- as.vector(bunzip2(filename = p2bz,
destname = file.path(getwd(), "Example_peak2.bed"),
ext="bz2", FUN=bzfile, overwrite=TRUE, remove = FALSE))
output <- peakcomp(bedInput1 = peak1_path, bedInput2 = peak2_path,
olap.rate = 0.1)

## End(Not run)

esATAC documentation built on Nov. 8, 2020, 6:58 p.m.