object |
A profileplyr object |
annotation_subset |
If specific annotations (from ChIPseeker package) are desired, specify them here in a character vector. Can be one or any combination of “Promoter”, “Exon”, “Intron”, “Downstream”, “Distal Intergenic”, “3p UTR”, or “5p UTR”. This argument is optional and all annotation types will be included if argument is left out |
TxDb |
This must be either a TxDb object, a character string that is a path to a GTF file, or character string indicating genome if one of the following - “hg19”, “hg38”, “mm9”, “mm10”. |
annoDb |
The annotation package to be used. If the ‘TxDb’ argument is set to “hg19”, “hg38”, “mm9”, or “mm10” this will automatically be set and this can be left as NULL. |
tssRegion |
This needs to be a vector of two numbers that will define promoter regions. The first number must be negative, while the second number must be positive. Default values are c(-3000, 3000) |
changeGroupToAnnotation |
If the grouping should be changed to the annotations (typically when the ranges will be exported for visualization based on this annotation), this should be TRUE. The default if FALSE, which will keep the grouping that existed before annotating the object. This is typical if the output will be used for finding overlaps with gene lists in the ‘groupBy’ function. . |
heatmap_grouping |
Only relevant if ‘keepAnnotationAsGroup’ is set to TRUE. This argument needs to be either “group”, or “annotation”. This will determine how the ranges are grouped in the resulting object. Default is heatmap_grouping = “Group”. If there are no groups in the deepTools matrix that was used in the function, this argument is unnecessary |
… |
pass argument to ‘annotatePeak’ function from the ChIPseeker package |
Example
proplyr_object_subset <- annotateRanges(proplyr_object_subset,
TxDb = "hg19")
## >> preparing features information... 2019-08-13 16:42:15
## >> identifying nearest features... 2019-08-13 16:42:16
## >> calculating distance from peak to TSS... 2019-08-13 16:42:17
## >> assigning genomic annotation... 2019-08-13 16:42:17
## >> adding gene annotation... 2019-08-13 16:42:33
## >> assigning chromosome lengths 2019-08-13 16:42:33
## >> done... 2019-08-13 16:42:33
rowRanges(proplyr_object_subset)[1:3]
## GRanges object with 3 ranges and 18 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <numeric>
## giID1 chr21 34696445-34697205 + | <NA> 0
## giID10 chr21 36207782-36209962 + | <NA> 0
## giID100 chr7 150734773-150737527 + | <NA> 0
## sgGroup giID names
## <factor> <factor> <factor>
## giID1 K27ac_top10_HUES64.bed giID1 <NA>
## giID10 K27ac_top10_HUES64.bed giID10 <NA>
## giID100 K27ac_top10_HUES64.bed giID100 <NA>
## annotation geneChr geneStart
## <character> <integer> <integer>
## giID1 Promoter (<=1kb) 21 34696734
## giID10 Promoter (<=1kb) 21 36169710
## giID100 Exon (uc003wik.4/11194, exon 10 of 16) 7 150745379
## geneEnd geneLength geneStrand geneId transcriptId
## <integer> <integer> <integer> <character> <character>
## giID1 34732128 35395 1 3454 uc011adv.2
## giID10 36206898 37189 2 861 uc002yul.1
## giID100 150749843 4465 1 9311 uc003win.3
## distanceToTSS ENSEMBL SYMBOL
## <numeric> <character> <character>
## giID1 0 ENSG00000142166 IFNAR1
## giID10 -884 ENSG00000159216 RUNX1
## giID100 -7852 ENSG00000213199 ASIC3
## GENENAME annotation_short
## <character> <ordered>
## giID1 interferon alpha and beta receptor subunit 1 Promoter
## giID10 RUNX family transcription factor 1 Promoter
## giID100 acid sensing ion channel subunit 3 Exon
## -------
## seqinfo: 24 sequences from an unspecified genome; no seqlengths
Visualize annotation of ranges
- Here we visualize the genomic annotation of each range by using the column name ‘annotation_short’ from the range metadata in the the ‘extra_annotation_columns’ argument
generateEnrichedHeatmap(proplyr_object_subset,
extra_annotation_columns = c("annotation_short"),
color_by_sample_group = "chip",
ylim = "chip",
show_heatmap_legend = rep(FALSE, 6)
)

Visualize annotation of ranges after clustering (using ‘clusterRanges’ function)
library(dplyr)
annotateRanges(proplyr_object_subset,
TxDb = "hg19") %>%
clusterRanges(kmeans_k = 3) %>%
generateEnrichedHeatmap(extra_annotation_columns = c("annotation_short"),
color_by_sample_group = "chip",
ylim = "chip",
show_heatmap_legend = rep(FALSE, 6)
)
## >> preparing features information... 2019-08-13 15:53:44
## >> identifying nearest features... 2019-08-13 15:53:44
## >> calculating distance from peak to TSS... 2019-08-13 15:53:45
## >> assigning genomic annotation... 2019-08-13 15:53:45
## >> adding gene annotation... 2019-08-13 15:53:46
## >> assigning chromosome lengths 2019-08-13 15:53:46
## >> done... 2019-08-13 15:53:46

NOTE: The ‘SYMBOL’ column of this object contains the gene annotation for that region, and importantly this can be used with the ‘groupBy’ function within profileplyr to group the ranges by gene sets.