Reads pairwise ChIA-PET interaction from an input file.

parseLoopsTang(inFile, ...)

Arguments

inFile

input file with loops

...

additional arguments, that will be passed to GRanges functions.

Value

An GInteractions with loops from input file.

Details

It reads files with the following tab-delimited format:

chr124816035148161634chr12482306654823284827
chr77728466477285815chr777388242773889287
chr4128459961128460166chr41285083041285090824
chr124816035148161634chr124823066548232848

This file format was used for ChIA-PET interaction data by Tang et al. 2015 http://dx.doi.org/10.1016/j.cell.2015.11.024. The last column of input file is added as annotation column with colname "score".

Examples

exampleLoopTang2015File <- system.file("extdata", "ChIA-PET_GM12878_Tang2015.chr22_1-30000000.clusters.txt", package = "sevenC") gi <- parseLoopsTang(exampleLoopTang2015File)
#> INFO: Parse ChiA-PET interactions from file: /home/ibnsalem/project/sevenC/inst/extdata/ChIA-PET_GM12878_Tang2015.chr22_1-30000000.clusters.txt
#> Parsed with column specification: #> cols( #> X1 = col_character(), #> X2 = col_integer(), #> X3 = col_integer(), #> X4 = col_character(), #> X5 = col_integer(), #> X6 = col_integer(), #> X7 = col_integer() #> )
# read loops with custom seqinfo object: customSeqInfo <- Seqinfo(seqnames = c("chr1", "chr22"), seqlengths = c(10^8, 10^8), isCircular = c(FALSE, FALSE), genome = "custom") gi <- parseLoopsTang(exampleLoopTang2015File, seqinfo = customSeqInfo)
#> INFO: Parse ChiA-PET interactions from file: /home/ibnsalem/project/sevenC/inst/extdata/ChIA-PET_GM12878_Tang2015.chr22_1-30000000.clusters.txt
#> Parsed with column specification: #> cols( #> X1 = col_character(), #> X2 = col_integer(), #> X3 = col_integer(), #> X4 = col_character(), #> X5 = col_integer(), #> X6 = col_integer(), #> X7 = col_integer() #> )