Description Usage Arguments Value See Also Examples
Annotate a set of ranges in a GRanges object with transcript type (txType) based on their genic context. Transcripts are obtained from a TxDb object, but can alternatively be specified manually as a GRangesList.
| 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 | assignTxType(object, txModels, ...)
## S4 method for signature 'GenomicRanges,GenomicRangesList'
assignTxType(
  object,
  txModels,
  outputColumn = "txType",
  swap = NULL,
  noOverlap = "intergenic"
)
## S4 method for signature 'RangedSummarizedExperiment,GenomicRangesList'
assignTxType(object, txModels, ...)
## S4 method for signature 'GenomicRanges,TxDb'
assignTxType(
  object,
  txModels,
  outputColumn = "txType",
  swap = NULL,
  tssUpstream = 100,
  tssDownstream = 100,
  proximalUpstream = 1000,
  detailedAntisense = FALSE
)
## S4 method for signature 'RangedSummarizedExperiment,TxDb'
assignTxType(object, txModels, ...)
 | 
| object | GRanges or RangedSummarizedExperiment: Ranges to be annotated. | 
| txModels | TxDb or GRangesList: Transcript models via a TxDb, or manually specified as a GRangesList. | 
| ... | additional arguments passed to methods. | 
| outputColumn | character: Name of column to hold txType 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. | 
| noOverlap | character: In case categories are manually supplied with as a GRangesList, what to call regions with no overlap. | 
| tssUpstream | integer: Distance to extend annotated promoter upstream. | 
| tssDownstream | integer: Distance to extend annotated promoter downstream. | 
| proximalUpstream | integer: Maximum distance upstream of promoter to be considered proximal. | 
| detailedAntisense | logical: Wether to mirror all txType categories in the antisense direction (TRUE) or lump them all together (FALSE). | 
object with txType added as factor column in rowData (or mcols)
Other Annotation functions: 
assignGeneID(),
assignMissingID(),
assignTxID()
| 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 | ## Not run: 
data(exampleUnidirectional)
# Obtain transcript models from a TxDb-object:
library(TxDb.Mmusculus.UCSC.mm9.knownGene)
txdb <- TxDb.Mmusculus.UCSC.mm9.knownGene
# Assign txIDs
assignTxType(exampleUnidirectional,
             txModels=txdb)
# Assign txIDs using only TC peaks:
exampleUnidirectional <- assignTxType(exampleUnidirectional,
                                      txModels=txdb,
                                      swap='thick')
# The 'promoter' and 'proximal' category boundaries can be changed:
assignTxType(exampleUnidirectional,
             txModels=txdb,
             swap='thick',
             tssUpstream=50,
             tssDownstream=50,
             proximalUpstream=100)
# Annotation using complete antisense categories:
exampleUnidirectional <- assignTxType(exampleUnidirectional,
                                    txModels=txdb,
                                    outputColumn='txTypeExtended',
                                    swap='thick',
                                    detailedAntisense=TRUE)
# The output is always a factor added as a column:
summary(rowRanges(exampleUnidirectional)$txType)
summary(rowRanges(exampleUnidirectional)$txTypeExtended)
# To avoid using a TxDb-object, a GRangesList can be supplied:
custom_hierarchy <- GRangesList(promoters=granges(promoters(txdb)),
                                exons=granges(exons(txdb)))
assignTxType(exampleUnidirectional,
             txModels=custom_hierarchy,
             outputColumn='customType',
             swap='thick',
             noOverlap = 'intergenic')
## End(Not run)
 | 
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.