annotateTranscripts: Annotate transcripts

View source: R/matchGenes.R

annotateTranscriptsR Documentation

Annotate transcripts

Description

Annotate transcripts

Usage

annotateTranscripts(txdb, annotationPackage = NULL, by = c("tx","gene"), codingOnly=FALSE, verbose = TRUE, requireAnnotation = FALSE, mappingInfo = NULL, simplifyGeneID = FALSE)

Arguments

txdb

A TxDb database object such as TxDb.Hsapiens.UCSC.hg19.knownGene

annotationPackage

An annotation data package from which to obtain gene/transcript annotation. For example org.Hs.eg.db. If none is provided the function tries to infer it from organism(txdb) and if it can't it proceeds without annotation unless requireAnnotation = TRUE.

by

Should we create a GRanges of transcripts (tx) or genes (gene).

codingOnly

Should we exclude all the non-coding transcripts.

verbose

logical value. If 'TRUE', it writes out some messages indicating progress. If 'FALSE' nothing should be printed.

requireAnnotation

logical value. If 'TRUE' function will stop if no annotation package is successfully loaded.

mappingInfo

a named list with elements 'column', 'keytype' and 'multiVals'. If specified this information will be used with mapIds when mapping the gene ids using annotationPackage. This is useful when working with a txdb object from ENSEMBL or GENCODE among other databases.

simplifyGeneID

logical value. If 'TRUE', gene ids will be shortened to before a dot is present in the id. This is useful for changing GENCODE gene ids to ENSEMBL ids.

Details

This function prepares a GRanges for the matchGenes function. It adds information and in particular adds exons information to each gene/transcript.

Value

A GRanges object with an attribute description set to annotatedTranscripts. The following columns are added. seqinfo is the information returned by seqinfo, CSS is the coding region start, CSE is the coding region end, Tx is the transcript ID used in TxDb, Entrez is the Entrez ID, Gene is the gene symbol, Refseq is the RefSeq annotation, Nexons is the number of exons, Exons is an IRanges with the exon information.

Author(s)

Harris Jaffee and Rafael A. Irizarry. 'mappingInfo' and 'simplifyGeneID' contributed by Leonardo Collado-Torres.

See Also

matchGenes

Examples

## Not run: 
    library("TxDb.Hsapiens.UCSC.hg19.knownGene")
    genes <- annotateTranscripts(TxDb.Hsapiens.UCSC.hg19.knownGene)
    
    ##and to avoid guessing the annotation package:
    genes <- annotateTranscripts(TxDb.Hsapiens.UCSC.hg19.knownGene,annotation="org.Hs.eg.db")
    

## End(Not run)

rafalab/bumphunter documentation built on March 20, 2024, 6:22 a.m.