Description Usage Arguments Details Value Author(s) Examples
Plots enhancer-promoter interactions
1 2 3 4 5 6 | plotInteractions(interactions, rangeGenes = NULL, selectGene = NULL,
selectRegulatoryRegion = NULL, filters = "hgnc_symbol",
benchInteractions = NULL, statistics = NULL, coloring = NULL,
interactionsNaming = "Enh~Prom Interactions",
benchmarkNaming = "Benchmark", sizes = 0.3, cex.title = 1,
plotAllAnchors1 = FALSE)
|
interactions |
A GInteractions object or a list of GInteractions
objects. Each GInteractions object has an anchor1 which corresponds to
regulatory region (enahncer), and anchor2 which should correspond to the TSS
of a gene (that info is used to plot promoters, thus please input TSS, not
gene coordinated). A minimum necessary info that should be present as
meta-data is info about gene names stored in the column "name2". Additionaly,
statistics (as "pval" column), and colors of the loops (as "color" column)
could be added as meta-data. Read more about arguments statistics and
coloring.
Such object can be produced by modelling [ |
rangeGenes |
(default NULL) By default, function retrieves
automatically gene annotations using TxDb.Hsapiens.UCSC.hg19.knownGene
package for genes present in the interactions object. Importantly,
hgnc_symbol or ensembl_gene_id are stored as the name2 column in the
interactions object.
However, by using this argument, one can define and import gene annotations
in formats such as TxDb, GRanges, GRangesList, data.frame,
character.In the case you want to import your own gene annotations, please
take a look at |
selectGene |
(default NULL) a character vector of gene name symbols (eg "FTO","IRX3", etc.) which user wants to plot. |
selectRegulatoryRegion |
(default NULL) GRanges object or a character vector which stores info about the region of interest which user wants to plot (format: "chr16:53112601-53114200"). This region does not necessarilly need to be equal to the regulatory regions reported in the interactions input objects, whereas it only needs to overlap some regulatory regions. In addition, to provide the context, a function plots not only regulatory regions queried, but all that other that are associated with genes which are in the first place associated with queried regions. Primary interactions (those who overlap with queried region) are plotted red, whereas secondary interactions (additional interactions associated with identified genes) are plotted in grey. |
filters |
default "hgnc_symbol"; other option "ensembl_gene_id".
One should define if the gene id's in the name2 column of the
interactions objects are ensembl gene ids or hgnc symbols.
In summary, filters define a restriction on the query for
|
benchInteractions |
(default NULL) A GInteractions object or a list of GInteractions objects storing info about interaction regions in the genome (for example, regions from the literature obtained using chromatin conformation capture related methods. |
statistics |
(default NULL) Column name where info about the heights of loops is stored. If NULL, then all interaction loops have an equal height Otherwise, loops of different height are plotted based on the info stored in the corresponding column, eg "pval" or "qval". Height is calculated as -log10(pval). |
coloring |
(default NULL) Information about colors used to plot enhancer promoter interactions. If it is NULL, then all interactions will be plotted red. Otherwise, a character, e.g. a name of the column of the interactions object. In selected column, info about loop color is stored. For example, interactions for different genes can be colored differently, or statistically significant/insignificant interactions can be colored differently. User pre-defines this column by itself. In addition, if selectRegulatoryRegion!=NULL, then primary selected interactions are colored red, wheres secondary interactions (additional interactions associated with identified genes) are plotted in grey |
interactionsNaming |
(character string) Track name for the interaction set |
benchmarkNaming |
(character string) Track name for the benchmark interaction set |
sizes |
(default 0.3) Default size of tracks. In the case that many tracks are plotted this argument should be adjusted (lower height of tracks) such that all tracks get plotted. |
cex.title |
(default NULL) Numeric. This argument controls the size of text which describes each track. If, many tracks are plotted, then text might disapear, but then cex.title argument should be set to smaller value, until text appears again |
plotAllAnchors1 |
(default FALSE). If interactions argument is a list of GInteractions objects, whether to plot all anchor1 or not, eg whether to plot separate track(s) for enhancer regions defined in different tracks, or only in one track |
Plots InteractionTrack, GeneRegionTrack, InteractionsTracks, AnnotationTrack, IdeogramTrack and GenomeAxisTrack for queried (list of) GInteraction object(s) [which contain info about enhancer-gene associations], and (list of) benchmark datasets.
User can select genes or regions of interest, and these regions will be selected from the input GInteractions object. Check description of arguments for this function for specific cases (coloring sheme). Names of the tracks correspond to the name of the list (in the case that list is used as an input), or one can define name for individual objects. E track (stands for enhancers) indicates locations of enhancer regions, wereas P track (stands for promoters) plots TSS locations for input genes.
Take a look at the plot :-) Gene annotations are retieved from TxDb.Hsapiens.UCSC.hg19.knownGene, and exons are plotted (no info about transcript names is used). Entrez -> hgcn symbol are obtained using biomaRt.
Nice plot.
Inga Patarcic
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 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 | library(Gviz)
library(GenomicInteractions)
library(GenomicFeatures)
library(GenomicRanges)
library(InteractionSet)
library(reg2gene)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(biomaRt)
# Creating an example GInteractions set:
enhancers <- GRanges(rep("chr16",6),
IRanges(c(53112601,55531601,53777201,
53778801,54084001,53946467),
c(53114200,55533250, 53778800,
53780400, 54084400 ,53947933)))
genes <- GRanges(rep("chr16",6),
IRanges(c(53737874, 54964773, 54320676,
53737874, 54964773, 54320676),
c(53737874, 54964773, 54320676,
53737874, 54964773, 54320676)))
GenomeInteractions <- GInteractions(enhancers,genes)
GenomeInteractions$name2 <- c("FTO","IRX5","IRX3")
GenomeInteractions$pval <- c(0.20857403, 0.72856090, 0.03586015,
0.32663439, 0.32534945, 0.03994488)
GenomeInteractions$color <- c("red","blue","grey")
# simplest example
# if no info about coloring and height of loops provided, all loops are equal
plotInteractions(interactions = GenomeInteractions)
# provide info about statistics * height of loops
plotInteractions(interactions = GenomeInteractions,
statistics ="pval",
coloring = "color")
# Specific gene can be individually plotted
plotInteractions(interactions = GenomeInteractions,
selectGene = "FTO")
# This should be equal to the example where genes are not selected
# if one wants to plot regulatory region
plotInteractions(interactions = GenomeInteractions,
selectRegulatoryRegion = "chr16:53112601-53114200")
# NOTE: red interactions correspond to the enhancer regions that
# overlap queried region, whereas grey interaction is interactions
#represent in the dataset and associated with a gene with which the
# queried regulatory region interacts (provides info about neighbourhood
# for queried interaction)
# and function plots not regulatory regions queried, but all that
# are present in the interactions and overlap that region
plotInteractions(interactions = GenomeInteractions,
selectRegulatoryRegion = "chr16:53112601-53778800",
coloring = "color",
statistics = "pval")
# more than one interaction track can be plotted, if input is list
GenomeInteractionsList <- list(GenomeInteractions,GenomeInteractions[1:3])
names(GenomeInteractionsList) <- c("EnhP1","EnhP2")
# names from the list are taken as an argument for Enh&Inter plots
plotInteractions(interactions = GenomeInteractionsList,
coloring = "color",
statistics = "pval",
interactionsNaming = "EP Interactions",
cex.title = 0.4,
sizes = 0.4)
# additionaly, benchmark data can be plotted as well
plotInteractions(interactions = GenomeInteractions,
coloring = "color",
statistics = "pval",
benchInteractions = GenomeInteractions[1:3])
# add a list of bemchmark data
benchInteractions = list(GenomeInteractions[1:3],GenomeInteractions[4:5],
GenomeInteractions[6])
names(benchInteractions) <- c("HiC","ChIA-PET","GTEx")
plotInteractions(interactions = GenomeInteractions,
coloring = "color",
statistics = "pval",
benchInteractions = benchInteractions,
# interactionsNaming = "EP Interactions",
cex.title = 0.4,
sizes = 0.4)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.