Title: | A R/Bioconductor package to detect the interaction nodes from HiC/HiChIP/HiCAR data |
---|---|
Description: | The GenomicInteractionNodes package can import interactions from bedpe file and define the interaction nodes, the genomic interaction sites with multiple interaction loops. The interaction nodes is a binding platform regulates one or multiple genes. The detected interaction nodes will be annotated for downstream validation. |
Authors: | Jianhong Ou [aut, cre], Yarui Diao [fnd] |
Maintainer: | Jianhong Ou <[email protected]> |
License: | file LICENSE |
Version: | 1.11.0 |
Built: | 2024-11-29 07:59:18 UTC |
Source: | https://github.com/bioc/GenomicInteractionNodes |
Assigne gene id and gene symbols to node regions by interacted.
annotateNodes( node_regions, txdb, orgDb, upstream = 2000, downstream = 500, ... )
annotateNodes( node_regions, txdb, orgDb, upstream = 2000, downstream = 500, ... )
node_regions |
GRanges object represent regions interacted with nodes. |
txdb |
An object of TxDb to extract gene information |
orgDb |
An object of OrgDb to extract gene symbols |
upstream , downstream
|
An integer(1) value indicating the number of bases upstream or downstream from the transcription start site. For additional details see promoters. |
... |
parameter can be passed to genes |
GRanges object with gene_id and symbols metadata.
library(TxDb.Hsapiens.UCSC.hg19.knownGene) ## for human hg19 library(org.Hs.eg.db) ## used to convert gene_id to gene_symbol set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene) annotateNodes(node_regions, TxDb.Hsapiens.UCSC.hg19.knownGene, org.Hs.eg.db)
library(TxDb.Hsapiens.UCSC.hg19.knownGene) ## for human hg19 library(org.Hs.eg.db) ## used to convert gene_id to gene_symbol set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene) annotateNodes(node_regions, TxDb.Hsapiens.UCSC.hg19.knownGene, org.Hs.eg.db)
Generate a list of random nodes used for example or test.
createRandomNodes( txdb, seq = "chr22", size = 500, upstream = 500, downstream = 500, maxDist = 1e+06, wid = 5000 )
createRandomNodes( txdb, seq = "chr22", size = 500, upstream = 500, downstream = 500, maxDist = 1e+06, wid = 5000 )
txdb |
An TxDb object. |
seq |
seqlevels to be kept. |
size |
the length of regions involved in nodes |
upstream , downstream
|
upstream or downstream for promoters |
maxDist |
maximal distance from promoters |
wid |
regions width. |
An GRanges object with comp_id metadata.
library(TxDb.Hsapiens.UCSC.hg19.knownGene) set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(TxDb.Hsapiens.UCSC.hg19.knownGene) set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene)
Define the interaction node from input Pairs.
detectNodes(interaction, pval_cutoff = 0.05, ...)
detectNodes(interaction, pval_cutoff = 0.05, ...)
interaction |
An object of Pairs to represent interactions. |
pval_cutoff |
Cutoff P value for interaction node by poisson distribution |
... |
Not used. |
A list of interaction nodes with elements: node_connection, Pairs object represent interactions interacted with nodes; nodes, GRanges object represent regions involved in nodes; node_regions, GRanges object represent regions interacted with nodes.
library(rtracklayer) p <- system.file("extdata", "WT.2.bedpe", package = "GenomicInteractionNodes") interactions <- import(con=p, format="bedpe") nodes <- detectNodes(interactions)
library(rtracklayer) p <- system.file("extdata", "WT.2.bedpe", package = "GenomicInteractionNodes") interactions <- import(con=p, format="bedpe") nodes <- detectNodes(interactions)
GO enrichment analysis for nodes
enrichmentAnalysis( node_regions, orgDb, onto = c("BP", "CC", "MF"), minGeneNum = 3, evidence = list(Experimental_evidence_codes = c("EXP", "IDA", "IPI", "IMP", "IGI", "IEP", "HTP", "HDA", "HMP", "HGI", "HEP"), `Phylogenetically-inferred_annotations` = c("IBA", "IBD", "IKR", "IRD"), Computational_analysis_evidence_codes = c("ISS", "ISO", "ISA", "ISM", "IGC", "RCA"), Author_statement_evidence_codes = c("TAS", "NAS"), Curator_statement_evidence_codes = c("IC", "ND"), Electronic_annotation_evidence_code = c("IEA")), ... )
enrichmentAnalysis( node_regions, orgDb, onto = c("BP", "CC", "MF"), minGeneNum = 3, evidence = list(Experimental_evidence_codes = c("EXP", "IDA", "IPI", "IMP", "IGI", "IEP", "HTP", "HDA", "HMP", "HGI", "HEP"), `Phylogenetically-inferred_annotations` = c("IBA", "IBD", "IKR", "IRD"), Computational_analysis_evidence_codes = c("ISS", "ISO", "ISA", "ISM", "IGC", "RCA"), Author_statement_evidence_codes = c("TAS", "NAS"), Curator_statement_evidence_codes = c("IC", "ND"), Electronic_annotation_evidence_code = c("IEA")), ... )
node_regions |
GRanges object represent regions interacted with nodes. The object must be annotated by annotateNodes with comp_id and gene_id in the metadata. |
orgDb |
An object of OrgDb to extract gene symbols. |
onto |
Ontology category. |
minGeneNum |
An integer(1) value indicating the minimal number of gene to start the enrichment analysis. If total gene counts is smaller than the 'minGeneNum', the NULL will be returned. |
evidence |
The acceptable evidence code. |
... |
Not used. |
A list with element enriched and enriched_in_compound. Or NULL if total counts of gene is smaller than 'minGeneNum'.
library(TxDb.Hsapiens.UCSC.hg19.knownGene) ## for human hg19 library(org.Hs.eg.db) ## used to convert gene_id to gene_symbol library(GO.db) set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene) node_regions <- annotateNodes(node_regions, TxDb.Hsapiens.UCSC.hg19.knownGene, org.Hs.eg.db) enr <- enrichmentAnalysis(node_regions, org.Hs.eg.db, onto="BP")
library(TxDb.Hsapiens.UCSC.hg19.knownGene) ## for human hg19 library(org.Hs.eg.db) ## used to convert gene_id to gene_symbol library(GO.db) set.seed(123) node_regions <- createRandomNodes(TxDb.Hsapiens.UCSC.hg19.knownGene) node_regions <- annotateNodes(node_regions, TxDb.Hsapiens.UCSC.hg19.knownGene, org.Hs.eg.db) enr <- enrichmentAnalysis(node_regions, org.Hs.eg.db, onto="BP")