GenomicInteractions: The GenomicInteractions class

Description Constructors Getters Setters Comparison methods Author(s) See Also Examples

Description

The GenomicInteractions class is derived from the Pairs class and is dedicated to representing interactions between pairs of genomic intervals. Each genomic interval is referred to as an “anchor”, and each interaction is defined between a pair of anchors. We will refer to the two anchors in a pair as the “first” and “second” anchors (arbitrarily defined).

All first anchors in a GenomicInteractions instance are represented by one GRangesFactor object, while all second anchors are represented by another GRangesFactor object of the same length. This improves efficiency by avoiding storage and processing of redundant copies of the same genomic coordinates and metadata; especially for interaction data, where one interval may be involved in multiple interactions.

It is worth clarifying the terminology at this point. The “anchors” refer to the genomic intervals themselves; the “anchor type” refers to the anchor ordering, i.e., first or second; the “anchor levels” are the universe of possible anchors across all interactions for a given anchor type; and the “anchor indices” are integer indices to the anchor levels, such that subsetting the anchor levels by the anchor indices yields the anchors.

Constructors

GenomicInteractions(anchor1, anchor2, regions, ..., metadata=list()) will create a GenomicInteractions object, given:

  1. A GRangesFactor object in each of anchor1 and anchor2, and missing regions.

  2. A GRanges object in each of anchor1 and anchor2, and missing regions.

  3. A GRanges object in each of anchor1 and anchor2, and a GRanges in regions that contains a superset of intervals in anchor1 and anchor2.

  4. An integer vector in each of anchor1 and anchor2, and a GRanges in regions. Anchor regions are defined by using integers to index regions.

Any arguments in ... should be vector-like and are added to the element-wise metadata of the output object. For options 1-3, any metadata in anchor1 or anchor2 are moved to the element-wise metadata, with the anchor of origin prepended to the name.

Option 2 has an additional common=TRUE argument, where the anchor levels for both the first and second anchors are defined as the union of the regions used in both anchors. This is provided for backwards compatibility with the old GenomicInteractions behaviour, and provides identical regions regardless of type (see below).

The metadata argument is expected to be a list-like object, to be used as the metadata slot of the output object.

Getters

In the following code snippets, x is a GenomicInteractions object:

anchors(x, type="both", id=FALSE):

Returns the anchor regions. If type="both", it returns a Pairs of length 2, where the first and second entries contain a GRangesFactor of the first and second anchors, respectively.

A GRangesFactor of only the first anchor regions can be returned directly with type=1 or type="first". The same applies for the second anchors with type=2 or type="second".

id=TRUE will return the integer indices pointing to the anchor regions in regions(x). If type="both", this will be a DataFrame of integer vectors, with one column per partner. Otherwise, the appropriate integer vector is directly returned for the partner specified by type.

regions(x, type=NA):

Returns the universe of all regions corresponding to anchor region type. Setting type=1 or "first" returns a GRanges of all first anchors, while setting type=2 or "second" does the same for all second anchors. Setting type="both" returns a List of length 2 containing the regions used in the first and second anchors.

type=NA behaves like type=1 with a deprecation warning. This is intended to encourage users to actually specify a value of type in calls to regions(x), as this option was not available in previous versions.

Note that type=1 and type=2 may not yield the same result! Anchor indices are only directly comparable if the universe of regions is the same for both anchors.

first(x):

A synonym for anchors(x, 1).

second(x):

A synonym for anchors(x, 2).

x$name:

Retrieves the column named name from mcols(x), yielding some vector-like object with the same length as x.

All getter methods applicable to Vector objects can also be used, e.g., mcols, metadata.

Setters

In the following code snippets, x is a GenomicInteractions object:

anchors(x, type="both", id=FALSE) <- value:

Replaces the anchor regions specified by type with value. If type="both", value should be a Pairs of two GRanges, where the first and second entries are to be used as the replacement first and second anchor regions, respectively. If type=1 or "first", value should be a GRanges or GRangesFactor to replace the first anchors; same for type=2 or "second" for the second anchors.

id=TRUE will replace the integer indices pointing to the anchor regions in regions(x). If type="both", value should be a DataFrame of integer vectors, with one column per partner. Otherwise, it should be an integer vector to replace the indices for the partner specified by type.

regions(x, type=NA) <- value:

Replaces the anchor levels in x. If type=1 or "first", value should be a GRanges to replace the first anchor levels; same for type=2 or "second" for the second anchor levels. If type="both", value should be a List of length 2 to replace both levels simultaneously.

type=NA is provided for backwards compatibility, and will replace the region sets for both anchors with the GRanges value. This will trigger a deprecation warning - users wanting this behaviour should use type="both" with the replacement value set to List(value, value).

first(x) <- value:

A synonym for anchors(x, 1) <- value.

second(x) <- value:

A synonym for anchors(x, 2) <- value.

x$name <- value:

Sets the column named name in mcols(x) with value, where value is recycled to have the same length as x.

Again, all setter methods applicable to Vector objects can also be used here, e.g., mcols<-, metadata<-.

Comparison methods

Methods like sort, match and pcompare are supported via inheritance from the Pairs class.

One entry X of a GenomicInteractions object is considered “greater” than another entry Y if the first anchor of X is “greater” than that of Y; if these are equal, then if the second anchor of X is greater than the the second anchor of Y. Comparisons between anchors are defined according to comparisons between GRanges, see ?"GenomicRanges-comparison".

Author(s)

Aaron Lun

See Also

GRangesFactor, used to represent the anchor regions.

Vector, from which this class is derived.

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
anchor1 <- GRanges(c("chr1", "chr1", "chr1", "chr1"), 
    IRanges(c(10, 20, 30, 20), width=5))
anchor1 <- sample(anchor1, 20, replace=TRUE)

anchor2 <- GRanges(c("chr1", "chr1", "chr1", "chr2"), 
    IRanges(c(100, 200, 300, 50), width=5))
anchor2 <- sample(anchor2, 20, replace=TRUE)

test <- GenomicInteractions(anchor1, anchor2)
test

# Alternative methods of construction:
combined <- sort(unique(c(anchor1, anchor2)))
m1 <- match(anchor1, combined)
m2 <- match(anchor1, combined)
test2 <- GenomicInteractions(m1, m2, combined)
test2

# Getting
anchors(test, 1)
anchors(test, 2) 
regions(test, 1)
regions(test, 2)

# Setting
anchors(test, 1) <- rev(anchors(test, 1))
test

first.regions <- regions(test, type=1)
regions(test, type=1) <- resize(first.regions, 1000)
test

LTLA/GenomicInteractions documentation built on June 24, 2019, 2:09 p.m.