assignTxID {CAGEfightR} | R Documentation |
Annotate a set of ranges in a GRanges object with transcript IDs based on their genic context. All overlapping transcripts are returned. Transcripts are obtained from a TxDb object, or can manually supplied as a GRanges.
assignTxID(object, txModels, ...) ## S4 method for signature 'GenomicRanges,GenomicRanges' assignTxID(object, txModels, outputColumn = "txID", swap = NULL) ## S4 method for signature 'RangedSummarizedExperiment,GenomicRanges' assignTxID(object, txModels, ...) ## S4 method for signature 'GenomicRanges,TxDb' assignTxID( object, txModels, outputColumn = "txID", swap = NULL, upstream = 1000, downstream = 0 ) ## S4 method for signature 'RangedSummarizedExperiment,TxDb' assignTxID(object, txModels, ...)
object |
GRanges or RangedSummarizedExperiment: Ranges to be annotated. |
txModels |
TxDb or GRanges: Transcript models via a TxDb, or manually specified as a GRanges. |
... |
additional arguments passed to methods. |
outputColumn |
character: Name of column to hold txID values. |
swap |
character or NULL: If not NULL, use another set of ranges contained in object to calculate overlaps, for example peaks in the thick column. |
upstream |
integer: Distance to extend annotated promoter upstream. |
downstream |
integer: Distance to extend annotated promoter downstream. |
object with txID added as a column in rowData (or mcols)
Other Annotation functions:
assignGeneID()
,
assignMissingID()
,
assignTxType()
data(exampleUnidirectional) # Obtain transcript models from a TxDb-object: library(TxDb.Mmusculus.UCSC.mm9.knownGene) txdb <- TxDb.Mmusculus.UCSC.mm9.knownGene # Assign txIDs assignTxID(exampleUnidirectional, txModels=txdb, outputColumn='geneID') # Assign txIDs using only TC peaks: assignTxID(exampleUnidirectional, txModels=txdb, outputColumn='geneID', swap='thick')