Finding the nearest genomic range neighbor in a TxDb
The distance
methods for TxDb objects and subclasses.
## S4 method for signature 'GenomicRanges,TxDb' distance(x, y, ignore.strand=FALSE, ..., id, type=c("gene", "tx", "exon", "cds"))
x |
The query GenomicRanges instance. |
y |
For |
id |
A |
type |
A |
ignore.strand |
A |
... |
Additional arguments for methods. |
distance:
Returns the distance for each range in x
to the range
extracted from the TxDb object y
. Values in
id
are matched to one of ‘gene_id’, ‘tx_id’,
‘exon_id’ or ‘cds_id’ identifiers in the TxDb
and the corresponding ranges are extracted. The type
argument
specifies which identifier is represented in id
. The extracted
ranges are used in the distance calculation with the ranges in
x
.
The method returns NA
values when the genomic region defined
by id
cannot be collapsed into a single range (e.g.,
when a gene spans multiple chromosomes) or if the id
is not found in y
.
The behavior of distance
with respect to zero-width ranges
has changed in Bioconductor 2.12. See the man page ?distance
in IRanges for details.
For distance
, an integer vector of distances between the ranges
in x
and y
.
Valerie Obenchain <vobencha@fhcrc.org>
nearest-methods man page in IRanges.
nearest-methods man page in GenomicRanges.
## ----------------------------------------------------------- ## distance() ## ----------------------------------------------------------- library(TxDb.Dmelanogaster.UCSC.dm3.ensGene) txdb <- TxDb.Dmelanogaster.UCSC.dm3.ensGene gr <- GRanges(c("chr2L", "chr2R"), IRanges(c(100000, 200000), width=100)) distance(gr, txdb, id=c("FBgn0259717", "FBgn0261501"), type="gene") distance(gr, txdb, id=c("10000", "23000"), type="cds") ## The id's must be in the appropriate order with respect to 'x'. distance(gr, txdb, id=c("4", "4097"), type="tx") ## 'id' "4" is on chr2L and "4097" is on chr2R. transcripts(txdb, filter=list(tx_id=c("4", "4097"))) ## If we reverse the 'id' the chromosomes are incompatable with gr. distance(gr, txdb, id=c("4097", "4"), type="tx") ## distance() compares each 'x' to the corresponding 'y'. ## If an 'id' is not found in the TxDb 'y' will not ## be the same lenth as 'x' and an error is thrown. ## Not run: distance(gr, txdb, id=c("FBgn0000008", "INVALID"), type="gene") ## will fail ## End(Not run)
Please choose more modern alternatives, such as Google Chrome or Mozilla Firefox.