getNearestFeature: Get nearest annotation boundary for a position range.

Description Usage Arguments Value Note See Also Examples

View source: R/hiAnnotator.R

Description

Given a query object, the function retrieves the nearest feature and its properties from a subject and then appends them as new columns within the query object. When used in genomic context, the function can be used to retrieve the nearest gene 5' or 3' end relative to genomic position of interest.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
getNearestFeature(
  sites.rd,
  features.rd,
  colnam = NULL,
  side = "either",
  feature.colnam = NULL,
  dists.only = FALSE,
  parallel = FALSE,
  relativeTo = "subject"
)

Arguments

sites.rd

GRanges object to be used as the query.

features.rd

GRanges object to be used as the subject or the annotation table.

colnam

column name to be added to sites.rd for the newly calculated annotation...serves a core!

side

boundary of annotation to use to calculate the nearest distance. Options are '5p','3p', 'either'(default), or 'midpoint'.

feature.colnam

column name from features.rd to be used for retrieving the nearest feature name. By default this is NULL assuming that features.rd has a column that includes the word 'name' somewhere in it.

dists.only

flag to return distances only. If this is TRUE, then 'feature.colnam' is not required and only distance to the nearest feature will be returned. By default this is FALSE.

parallel

use parallel backend to perform calculation with foreach. Defaults to FALSE. If no parallel backend is registered, then a serial version of foreach is ran using registerDoSEQ.

relativeTo

calculate distance relative to query or subject. Default is 'subject'. This essentially means whether to use query or subject as the anchor point to get distance from!

Value

a GRanges object with new annotation columns appended at the end of sites.rd.

Note

See Also

makeGRanges, getFeatureCounts, getSitesInFeature, get2NearestFeature.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
# Convert a dataframe to GRanges object
data(sites)
alldata.rd <- makeGRanges(sites, soloStart = TRUE)

data(genes)
genes.rd <- makeGRanges(genes)

nearestGenes <- getNearestFeature(alldata.rd, genes.rd, "NearestGene")
nearestGenes
nearestGenes <- getNearestFeature(alldata.rd, genes.rd, "NearestGene",
side = "5p")
nearestGenes
## Not run: 
nearestGenes <- getNearestFeature(alldata.rd, genes.rd, "NearestGene",
side = "3p")
nearestGenes
nearestGenes <- getNearestFeature(alldata.rd, genes.rd, "NearestGene",
side = "midpoint")
## Parallel version of getNearestFeature
nearestGenes <- getNearestFeature(alldata.rd, genes.rd, "NearestGene",
parallel = TRUE)
nearestGenes

## End(Not run)

malnirav/hiAnnotator documentation built on July 29, 2021, 6:32 a.m.