Annotate a GRanges
object with
sets of reference GRanges
or
GRangesList
objects, with
respect to overlaps and nearest neighbors.
The GRanges
object to annotate.
Named list
with GRanges
or GRangesList
object(s).
For each list element, a logical vector "X.hasOverlap" will be added to the
mcols
of the result, with TRUE
for each tile that overlaps
any region in that element. "X" is obtained from names(hasOverlap)
.
Named list
with GRanges
or GRangesList
object(s).
For each list element, a numeric vector "X.fracOverlap" will be added to the
mcols
of the result, with a value between 0 and 1 giving the fraction
of bases in a tile that overlaps with any region in that element. "X" is
obtained from names(fracOverlap)
.
Named list
with GRanges
or GRangesList
object(s).
For each list element, two numeric vectors "X.numOverlapWithin" and
"X.numOverlapAny" will be added to the mcols
of the result, giving
the number of ranges in that element that are fully contained within
a tile, or that overlap with a tile in any way, respectively. "X" is
obtained from names(numOverlap)
.
Named list
with GRanges
or GRangesList
object(s).
For each list element, two numeric vectors "X.nearestName" and
"X.nearestDistance" will be added to the mcols
of the result, giving
the name and distance of the nearest range in that element for each tile. "X" is
obtained from names(nearest)
, and the values of "X.nearestName" from
names(nearest$X)
. If multiple nearest ranges are at the same
distance from a tile, an arbitrary one is reported in "X.nearestName".
Logical scalar passed to
findOverlaps
when searching for overlaps
between x
and reference regions.
A GRanges
similar to x
, with
annotations added to its metadata columns (mcols
).
getGenomicTiles
that uses this function,
findOverlaps
and
nearest
in package GenomicRanges used
internally.
library(GenomicRanges)
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> Loading required package: generics
#>
#> Attaching package: ‘generics’
#> The following objects are masked from ‘package:base’:
#>
#> as.difftime, as.factor, as.ordered, intersect, is.element, setdiff,
#> setequal, union
#>
#> Attaching package: ‘BiocGenerics’
#> The following objects are masked from ‘package:stats’:
#>
#> IQR, mad, sd, var, xtabs
#> The following objects are masked from ‘package:base’:
#>
#> Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#> as.data.frame, basename, cbind, colnames, dirname, do.call,
#> duplicated, eval, evalq, get, grep, grepl, is.unsorted, lapply,
#> mapply, match, mget, order, paste, pmax, pmax.int, pmin, pmin.int,
#> rank, rbind, rownames, sapply, saveRDS, table, tapply, unique,
#> unsplit, which.max, which.min
#> Loading required package: S4Vectors
#>
#> Attaching package: ‘S4Vectors’
#> The following object is masked from ‘package:utils’:
#>
#> findMatches
#> The following objects are masked from ‘package:base’:
#>
#> I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
x <- GRanges("chr1", IRanges(c(1, 12), width = 10))
tss <- GRanges("chr1", IRanges(c(1, 10, 30), width = 1,
names = paste0("t", 1:3)))
blacklist <- GRanges("chr1", IRanges(20, width = 5))
annotateRegions(x, hasOverlap = list(Blacklist = blacklist),
fracOverlap = list(Blacklist = blacklist),
numOverlap = list(TSS = tss),
nearest = list(TSS = tss))
#> GRanges object with 2 ranges and 6 metadata columns:
#> seqnames ranges strand | Blacklist.hasOverlap Blacklist.fracOverlap
#> <Rle> <IRanges> <Rle> | <logical> <numeric>
#> [1] chr1 1-10 * | FALSE 0.0
#> [2] chr1 12-21 * | TRUE 0.2
#> TSS.numOverlapWithin TSS.numOverlapAny TSS.nearestName
#> <integer> <integer> <character>
#> [1] 2 2 t2
#> [2] 0 0 t2
#> TSS.nearestDistance
#> <integer>
#> [1] 0
#> [2] 1
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths