annotateRegions: Assign genomic states to regions

View source: R/annotateRegions.R

annotateRegionsR Documentation

Assign genomic states to regions

Description

This function takes the regions found in calculatePvalues and assigns them genomic states contructed with makeGenomicState. The main workhorse functions are countOverlaps and findOverlaps.

Usage

annotateRegions(regions, genomicState, annotate = TRUE, ...)

Arguments

regions

The ⁠$regions⁠ output from calculatePvalues.

genomicState

A GRanges object created with makeGenomicState. It can be either the genomicState$fullGenome or genomicState$codingGenome component.

annotate

If TRUE then the regions are annotated by the genomic state. Otherwise, only the overlaps between the regions and the genomic states are computed.

...

Arguments passed to other methods and/or advanced arguments. Advanced arguments:

verbose

If TRUE basic status updates will be printed along the way.

ignore.strand

Passed on to findOverlaps-methods and countOverlaps. Default: TRUE.

Passed to extendedMapSeqlevels, countOverlaps and findOverlaps-methods.

Details

You might want to specify arguments such as minoverlap to control how the overlaps are determined. See findOverlaps for further details.

Value

A list with elements countTable and annotationList (only if annotate=TRUE).

countTable

This is a data.frame with the number of overlaps from the regions vs the genomic states with one type per column. For example, if fullOrCoding='full' then the columns are exon, intergenic and intron.

annotationList

This is a GRangesList with the genomic states that overlapped with the regions. The names of this GRangesList correspond to the region index in regions.

Author(s)

Andrew Jaffe, Leonardo Collado-Torres

See Also

makeGenomicState, calculatePvalues

Examples

## Annotate regions, first two only
annotatedRegions <- annotateRegions(
    regions = genomeRegions$regions[1:2],
    genomicState = genomicState$fullGenome, minoverlap = 1
)
annotatedRegions

lcolladotor/derfinder documentation built on May 10, 2023, 11:07 p.m.