The linkSet object is a container for storing gene-enhancer interactions.
Slots
nameBaitA character vector of the bait names.
anchor1A integer vector of the first anchor indices.
anchor2A integer vector of the second anchor indices.
regionsA GenomicRanges object of the regions.
NAMESA character vector of the region names.
elementMetadataA DataFrame of the element metadata.
Examples
showClass("linkSet") # shows the known subclasses
#> Class "linkSet" [package "linkSet"]
#>
#> Slots:
#>
#> Name: nameBait anchor1
#> Class: character integer
#>
#> Name: anchor2 regions
#> Class: integer GenomicRanges_OR_missing
#>
#> Name: NAMES elementMetadata
#> Class: character_OR_NULL DataFrame
#>
#> Name: metadata
#> Class: list
#>
#> Extends:
#> Class "Vector", directly
#> Class "Annotated", by class "Vector", distance 2
#> Class "vector_OR_Vector", by class "Vector", distance 2
set.seed(7000)
N <- 40
all.starts <- round(runif(N, 1, 100))
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))
genes = c(rep("SP7",4),rep("ASPN",10),rep("XBP1",6))
Np <- 20
all.anchor1 <- sample(N, Np)
gr1 <- all.regions[all.anchor1]
gr1$symbol <- genes
all.anchor2 <- setdiff(1:40,all.anchor1)
gr2 <- all.regions[all.anchor2]
x <- linkSet(gr1, gr2,specificCol = "symbol")
x
#> linkSet object with 20 interactions and 1 metadata column:
#> bait seqnames_oe ranges_oe | anchor1.symbol
#> <character> <Rle> <IRanges> | <character>
#> [1] SP7 --- chrA 47-58 | SP7
#> [2] SP7 --- chrA 41-54 | SP7
#> [3] SP7 --- chrA 20-33 | SP7
#> [4] SP7 --- chrA 87-104 | SP7
#> [5] ASPN --- chrA 59-76 | ASPN
#> ... ... ... ... ... . ...
#> [16] XBP1 --- chrB 54-61 | XBP1
#> [17] XBP1 --- chrB 1-18 | XBP1
#> [18] XBP1 --- chrB 83-91 | XBP1
#> [19] XBP1 --- chrB 33-45 | XBP1
#> [20] XBP1 --- chrB 96-114 | XBP1
#> -------
#> regions: 40 ranges and 0 metadata columns
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
x2 <- linkSet(genes, gr2)
x2
#> linkSet object with 20 interactions and 0 metadata columns:
#> bait seqnames_oe ranges_oe
#> <character> <Rle> <IRanges>
#> [1] SP7 --- chrA 47-58
#> [2] SP7 --- chrA 41-54
#> [3] SP7 --- chrA 20-33
#> [4] SP7 --- chrA 87-104
#> [5] ASPN --- chrA 59-76
#> ... ... ... ... ...
#> [16] XBP1 --- chrB 54-61
#> [17] XBP1 --- chrB 1-18
#> [18] XBP1 --- chrB 83-91
#> [19] XBP1 --- chrB 33-45
#> [20] XBP1 --- chrB 96-114
#> -------
#> regions: 20 ranges and 0 metadata columns
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
