GInteractions-class: GInteractions class and constructors

Description Usage Arguments Value Description of the class Class construction Author(s) See Also Examples

Description

The GInteractions class stores pairwise genomic interactions, and is intended for use in data analysis from Hi-C or ChIA-PET experiments. Each row of the GInteractions corresponds to a pairwise interaction between two loci, with indexing to improve computational efficiency.

Usage

1
2
3
4
5
6
7
8
## S4 method for signature 'numeric,numeric,GRanges'
GInteractions(anchor1, anchor2, regions, metadata=list(), mode="normal", ...)

## S4 method for signature 'GRanges,GRanges,GenomicRanges_OR_missing'
GInteractions(anchor1, anchor2, regions, metadata=list(), mode="normal", ...)

## S4 method for signature 'missing,missing,GenomicRanges_OR_missing'
GInteractions(anchor1, anchor2, regions, metadata=list(), mode="normal", ...)

Arguments

anchor1, anchor2

Either a pair of numeric vectors containing indices to regions, or a pair of GRanges objects specifying the interacting loci. Lengths of both anchor1 and anchor2 must be equal.

regions

A GRanges object containing the coordinates of the interacting regions. This is only mandatory if anchor1 and anchor2 are numeric vectors.

metadata

An optional list of arbitrary content describing the overall experiment.

mode

A string indicating what type of GInteractions object should be constructed.

...

Optional metadata columns.

Value

For the constructors, a GInteractions (or StrictGInteractions, or ReverseStrictGInteractions) object is returned.

Description of the class

The GInteractions class inherits from the Vector class and has access to all of its data members and methods (e.g, metadata and elementMetadata - see Vector-class for more details). It also contains several additional slots:

anchor1:

An integer vector specifying the index of the first interacting region.

anchor2:

An integer vector specifying the index of the second interacting region.

regions:

A sorted GRanges object containing the coordinates of all interacting regions.

Each interaction is defined by the corresponding entries in the anchor1 and anchor2 slots, which point to the relevant coordinates in regions for each locus.

The StrictGInteractions class inherits from the GInteractions class and has the same slots. The only difference is that, for each interaction, anchor1 must be less than or equal to anchor2. This means that the first interacting region has a start position that is "lower" than the second interacting region. This condition is useful for comparing interactions within and between objects, as it ensures that redundant permutations of the same interaction are not being overlooked. However, it is not used by default as there may conceivably be instances where the order of interactions is informative. The ReverseStrictGInteractions class has the opposite behaviour, where anchor1 must be greater than or equal to anchor2.

Class construction

GInteractions objects can be constructed by specifying integer vectors to define the pairwise interactions in the GInteractions call. For entry x, the corresponding interaction is defined between regions[anchor1[x]] and regions[anchor2[x]]. Obviously, coordinates of all of the interacting loci must be specified in the regions argument. Any metadata in regions will be preserved. Note that regions will be resorted in the returned object, so the anchor indices may not be equal to the input anchor1 and anchor2.

Alternatively, GInteractions objects can be constructed by directly supplying the GRanges of the interacting loci to the GInteractions function. If regions is not specified, this will be constructed automatically from the two sets of supplied GRanges. If regions is supplied, exact matching will be performed to identify the indices in regions corresponding to the regions in the supplied anchor GRanges. Missing values are not tolerated and will cause an error to be raised. In both cases, any metadata in the input GRanges will be transferred to the mcols of the output GInteractions object.

All constructors will return a GInteractions object containing all of the specified information. Sorting of regions is performed automatically, with re-indexing of all anchor indices to preserve the correct pairings between regions. If mode="strict", a StrictGInterctions object is returned with anchor indices swapped such that anchor1 <= anchor2 for all interactions. If mode="reverse", a ReverseStrictGInterctions object is returned with anchor indices swapped such that anchor1 >= anchor2. If both anchors are missing, the constructor will return an empty GInteractions object.

Author(s)

Aaron Lun, with contributions from Malcolm Perry and Liz Ing-Simmons.

See Also

interaction-access, interaction-subset, interaction-compare, Vector-class

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
set.seed(1000)
N <- 30
all.starts <- sample(1000, N)
all.ends <- all.starts + round(runif(N, 5, 20))
all.regions <- GRanges(rep(c("chrA", "chrB"), c(N-10, 10)),
    IRanges(all.starts, all.ends))

Np <- 20
all.anchor1 <- sample(N, Np)
all.anchor2 <- sample(N, Np)
gi <- GInteractions(all.anchor1, all.anchor2, all.regions)

# Equivalent construction:
GInteractions(all.regions[all.anchor1], all.regions[all.anchor2])
GInteractions(all.regions[all.anchor1], all.regions[all.anchor2], 
    all.regions)

# Putting in metadata, elementMetadata
temp.gi <- gi
metadata(temp.gi)$name <- "My GI object"
mcols(temp.gi)$score <- runif(Np)

# Strict construction
sgi <- GInteractions(all.regions[all.anchor1], all.regions[all.anchor2], 
        all.regions, mode="strict")
rsgi <- GInteractions(all.regions[all.anchor1], all.regions[all.anchor2], 
        all.regions, mode="reverse")

InteractionSet documentation built on April 17, 2021, 6 p.m.