TxDb-class: TxDb objects

TxDb-classR Documentation

TxDb objects

Description

The TxDb class is a container for storing transcript annotations.

Methods

In the code snippets below, x is a TxDb object.

metadata(x):

Return x's metadata in a data frame.

seqlevels0(x):

Get the sequence levels originally in x. This ignores any change the user might have made to the sequence levels with the seqlevels setter.

seqlevels(x), seqlevels(x) <- value:

Get or set the sequence levels in x.

seqinfo(x), seqinfo(x) <- value:

Get or set the information about the underlying sequences. Note that, for now, the setter only supports replacement of the sequence names, i.e., except for their sequence names (accessed with seqnames(value) and seqnames(seqinfo(x)), respectively), Seqinfo objects value (supplied) and seqinfo(x) (current) must be identical.

isActiveSeq(x):

Return the currently active sequences for this txdb object as a named logical vector. Only active sequences will be tapped when using the supplied accessor methods. Inactive sequences will be ignored. By default, all available sequences will be active.

isActiveSeq(x) <- value:

Allows the user to change which sequences will be actively accessed by the accessor methods by altering the contents of this named logical vector.

seqlevelsStyle(x), seqlevelsStyle(x) <- value:

Get or set the seqname style for x. See the seqlevelsStyle generic getter and setter in the GenomeInfoDb package for more information.

as.list(x):

Dump the entire db into a list of data frames, say txdb_dump, that can then be used to recreate the original db with do.call(makeTxDb, txdb_dump) with no loss of information (except possibly for some of the metadata). Note that the transcripts are dumped in the same order in all the data frames.

Author(s)

Hervé Pagès, Marc Carlson

See Also

  • makeTxDbFromUCSC, makeTxDbFromBiomart, and makeTxDbFromEnsembl, for making a TxDb object from online resources.

  • makeTxDbFromGRanges and makeTxDbFromGFF for making a TxDb object from a GRanges object, or from a GFF or GTF file.

  • saveDb and loadDb in the AnnotationDbi package for saving and loading a TxDb object as an SQLite file.

  • transcripts, transcriptsBy, and transcriptsByOverlaps, for extracting genomic feature locations from a TxDb-like object.

  • transcriptLengths for extracting the transcript lengths (and other metrics) from a TxDb object.

  • select-methods for how to use the simple "select" interface to extract information from a TxDb object.

  • The Seqinfo class in the GenomeInfoDb package.

Examples

txdb_file <- system.file("extdata", "Biomart_Ensembl_sample.sqlite",
                         package="GenomicFeatures")
txdb <- loadDb(txdb_file)
txdb

## Use of seqinfo():
seqlevelsStyle(txdb)
seqinfo(txdb)
seqlevels(txdb)
seqlengths(txdb)  # shortcut for 'seqlengths(seqinfo(txdb))'
isCircular(txdb)  # shortcut for 'isCircular(seqinfo(txdb))'
names(which(isCircular(txdb)))

## You can set user-supplied seqlevels on 'txdb' to restrict any further
## operations to a subset of chromosomes:
seqlevels(txdb) <- c("Y", "6")
## Then you can restore the seqlevels stored in the db:
seqlevels(txdb) <- seqlevels0(txdb)

## Use of as.list():
txdb_dump <- as.list(txdb)
txdb_dump
txdb1 <- do.call(makeTxDb, txdb_dump)
stopifnot(identical(as.list(txdb1), txdb_dump))

Bioconductor/GenomicFeatures documentation built on Nov. 7, 2024, 4:25 a.m.