Package 'annoLinker'

Title: Annotating genomic regions through chromatin interaction links
Description: Fast annotation of genomic peaks using DNA interaction data by constructing interaction networks with igraph, where peaks overlapping any node in a connected subgraph are annotated with all genes in that subgraph. The annotation evidence could be visualized as either a network graph or a genomic track integrated with gene annotation information.
Authors: Jianhong Ou [aut, cre] (ORCID: <https://orcid.org/0000-0002-8652-2488>), Kenneth Poss [aut, fnd]
Maintainer: Jianhong Ou <[email protected]>
License: GPL-3
Version: 1.1.0
Built: 2026-05-29 08:53:14 UTC
Source: https://github.com/bioc/annoLinker

Help Index


Annotating genomic regions through chromatin interaction links

Description

Fast annotation of genomic peaks using DNA interaction data by building interaction networks with igraph. Peaks overlapping any node in a connected subgraph are annotated with all genes in that subgraph.

Author(s)

Jianhong Ou

Maintainer: Jianhong Ou [email protected]

See Also

Useful links:


Annotate Peaks with DNA Interaction Networks Using Graph Clustering

Description

Fast annotation of genomic peaks using DNA interaction data by building interaction networks with igraph. Peaks overlapping any node in a connected subgraph are annotated with all genes in that subgraph.

Usage

annoLinker(
  peaks,
  annoData,
  interactions,
  bindingType = c("startSite", "body", "endSite"),
  bindingRegion = c(-5000, 5000),
  cluster_method = c("components", "louvain", "walktrap", "infomap"),
  maxgap = 0,
  interactionDistanceRange = c(10000, 1e+07),
  addEvidence = FALSE,
  parallel = FALSE,
  verbose = FALSE,
  ...
)

Arguments

peaks

GRanges object containing peak regions

annoData

annoGR or GRanges object with gene annotations

interactions

GInteractions, or Pairs object with interaction data (e.g., Hi-C, ChIA-PET)

bindingType

Character, one of "startSite", "body", or "endSite"

bindingRegion

Numeric vector of length 2 defining promoter window (e.g., c(-5000, 5000))

cluster_method

Character, clustering method: "components" (connected components), "louvain", "walktrap", or "infomap"

maxgap

Integer, bp to extend interaction anchors for overlap detection (default: 0)

interactionDistanceRange

Numeric vector of length 2 defining the minimal and maximal distance of interactions. This is used to make sure the annotations are not supper far away.

addEvidence

Logical, add evidence to the metadata or not.

parallel

Logical, use future_lapply to do parallel computing or not.

verbose

Logical, print the message or not

...

Parameters for cluster. see cluster_louvain, cluster_walktrap, and cluster_infomap.

Value

An annoLinkerResult object or NULL if no annotations found

Examples

## read the peaks and interactions
library(rtracklayer)
extPath <- system.file("extdata", package = "annoLinker")
peaks <- rtracklayer::import(file.path(extPath, "peaks.bed"))
interactions <- rtracklayer::import(file.path(extPath, "interaction.bedpe"))
library(TxDb.Drerio.UCSC.danRer10.refGene)
annoData <- genes(TxDb.Drerio.UCSC.danRer10.refGene)
anno <- annoLinker(peaks, annoData, interactions, verbose = TRUE)

Class "annoLinkerResult"

Description

An object of class "annoLinkerResult" represent the annotated peaks, which is a GRanges object with peaks annotated by gene clusters, and interaction graph, which is an igraph graph.

Usage

## S4 method for signature 'annoLinkerResult'
as.data.frame(x, row.names = NULL, optional = FALSE, ...)

anno_peaks(x)

## S4 method for signature 'annoLinkerResult'
anno_peaks(x)

anno_graph(x)

## S4 method for signature 'annoLinkerResult'
anno_graph(x)

anno_clusters(x)

## S4 method for signature 'annoLinkerResult'
anno_clusters(x)

anno_evidence(x, i)

## S4 method for signature 'annoLinkerResult'
anno_evidence(x, i)

anno_event(x, i)

## S4 method for signature 'annoLinkerResult'
anno_event(x, i)

anno_feature(x, i)

## S4 method for signature 'annoLinkerResult'
anno_feature(x, i)

anno_peakbin(x, i)

## S4 method for signature 'annoLinkerResult'
anno_peakbin(x, i)

anno_featurebin(x, i)

## S4 method for signature 'annoLinkerResult'
anno_featurebin(x, i)

## S4 method for signature 'annoLinkerResult'
length(x)

## S4 method for signature 'annoLinkerResult'
show(object)

## S4 method for signature 'annoLinkerResult'
head(x, ...)

Arguments

x, object

An annoLinkerResult object.

row.names, optional, ...

parameters used by as.data.frame

i

Numeric, index value.

Value

The object of 'annoLinkerResult', 'GRanges', 'igraph' or 'data.frame'

Objects from the Class

Objects can be created by calls of the form new("annoLinkerResult", annotated_peaks, graph, clusters).

Examples

library(igraph)
library(GenomicRanges)
new("annoLinkerResult",
    annotated_peaks = GRanges(),
    graph = make_empty_graph(),
    clusters = data.frame()
)

Plot interaction network for visualization

Description

Plot interaction network for visualization

Usage

plotEvidence(
  anno,
  event,
  output = c("graph", "htmlWidget", "trackPlot"),
  colors = c(peak = "darkgreen", feature = "brown", node = "tomato", background =
    "lightgray"),
  txdb,
  org
)

Arguments

anno

An object of annoLinkerResult output by annoLinker

event

Number to indicate the event to be plot

output

Output of the plot.

colors

Colors setting for the plot.

txdb, org

The TxDb and OrgDb object used for annotation plot.

Value

htmlWidget or plots.

Examples

anno <- readRDS(system.file("extdata", "sample_res.rds",
    package = "annoLinker"
))
library(org.Dr.eg.db)
library(TxDb.Drerio.UCSC.danRer10.refGene)
n <- 1 # length(anno$annotated_peaks$evidences)
plotEvidence(anno,
    event = n,
    output = "htmlWidget"
)
plotEvidence(anno,
    event = n,
    output = "trackPlot"
)