---
title: "Obtaining and Utilizing TxDb Objects"
author:
- name: Marc Carlson
- name: Patrick Aboyoun
- name: Hervé Pagès
- name: Seth Falcon
- name: Martin Morgan
date: "`r format(Sys.Date(), '%A, %B %d, %Y')`"
package: GenomicFeatures
vignette: |
  %\VignetteIndexEntry{Obtaining and Utilizing TxDb Objects}
  %\VignetteEngine{knitr::rmarkdown}
  %\VignetteEncoding{UTF-8}
  %\VignetteDepends{GenomicFeatures,TxDb.Hsapiens.UCSC.hg19.knownGene,BSgenome.Hsapiens.UCSC.hg19}
  %\VignetteKeywords{annotation, database}
output:
  BiocStyle::html_document
---



# Introduction

The `GenomicFeatures` package implements the `TxDb` container for
storing transcript metadata for a given organism. A `TxDb` object
stores the genomic positions of the 5' and 3' untranslated regions
(UTRs), protein coding sequences (CDSs), and exons for a set of mRNA
transcripts. The genomic positions are stored and reported with respect
to a given genome assembly. `TxDb` objects have numerous accessors
functions to allow such features to be retrieved individually or grouped
together in a way that reflects the underlying biology.

All `TxDb` objects are backed by a SQLite database that stores
the genomic positions and relationships between pre-processed mRNA
transcripts, exons, protein coding sequences, and their related
gene identifiers.



# Installing the `GenomicFeatures` package

Install the package with:

```{r installgenomicfeatures, eval=FALSE}
if (!require("BiocManager", quietly=TRUE))
    install.packages("BiocManager")

BiocManager::install("GenomicFeatures")
```

Then load it with:

```{r loadgenomicfeatures}
suppressPackageStartupMessages(library(GenomicFeatures))
```



# Obtaining a `TxDb` object

There are three ways that users can obtain a `TxDb` object.

One way is to use the `loadDb` method to load the object directly
from an appropriate `.sqlite` database file.

Here we are loading a previously created `TxDb` object
based on UCSC known gene data.  This database only contains a small
subset of the possible annotations for human and is only included to
demonstrate and test the functionality of the
`GenomicFeatures` package as a demonstration.

```{r loadDb}
samplefile <- system.file("extdata", "hg19_knownGene_sample.sqlite",
                          package="GenomicFeatures")
txdb <- loadDb(samplefile)
txdb
```

In this case, the `TxDb` object has been returned by
the `loadDb` method.

More commonly however, we expect that users will just load a
TxDb annotation package like this:

```{r loadPackage}
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
hg19_txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene  # shorthand (for convenience)
hg19_txdb
```

Loading the package like this will also create a `TxDb`
object, and by default that object will have the same name as the
package itself.

Finally, the third way to obtain a `TxDb` object is to use one of
the numerous tools defined in the `txdbmaker` package. `txdbmaker`
provides a set of tools for making `TxDb` objects from genomic
annotations from various sources (e.g. UCSC, Ensembl, and GFF files).
See the vignette in the `txdbmaker` package for more information.




# Retrieving Data from a `TxDb` object


## Pre-filtering data based on Chromosomes

It is possible to filter the data that is returned from a
`TxDb` object based on it's chromosome.  This can be a
useful way to limit the things that are returned if you are only
interested in studying a handful of chromosomes.

To determine which chromosomes are currently active, use the
`seqlevels` method.  For example:

```{r seqlevels}
head(seqlevels(hg19_txdb))
```

Will tell you all the chromosomes that are active for the
TxDb.Hsapiens.UCSC.hg19.knownGene `TxDb` object (by
default it will be all of them).

If you then wanted to only set Chromosome 1 to be active you could do
it like this:
```{r seqlevels2}
seqlevels(hg19_txdb) <- "chr1"
```

So if you ran this, then from this point on in your R session only
chromosome 1 would be consulted when you call the various retrieval
methods...  If you need to reset back to the original seqlevels (i.e.
to the seqlevels stored in the db), then set the seqlevels to
`seqlevels0(hg19_txdb)`.
```{r seqlevels3}
seqlevels(hg19_txdb) <- seqlevels0(hg19_txdb)
```

**Exercise:**
Use `seqlevels` to set only chromsome 15 to be active.  BTW,
the rest of this vignette will assume you have succeeded at this.

**Solution:**
```{r seqlevels4}
seqlevels(hg19_txdb) <- "chr15"
seqlevels(hg19_txdb)
```


## Retrieving data using the `select()` method

The `TxDb` objects inherit from `AnnotationDb`
objects (just as the `ChipDb` and `OrgDb` objects do).
One of the implications of this relationship is that these object
ought to be used in similar ways to each other.  Therefore we have
written supporting `columns`, `keytypes`, `keys`
and `select` methods for `TxDb` objects.

These methods can be a useful way of extracting data from a
`TxDb` object.  And they are used in the same way that
they would be used to extract information about a `ChipDb` or
an `OrgDb` object.  Here is a simple example of how to find the
UCSC transcript names that match with a set of gene IDs.

```{r selectExample}
keys <- c("100033416", "100033417", "100033420")
columns(hg19_txdb)
keytypes(hg19_txdb)
select(hg19_txdb, keys = keys, columns="TXNAME", keytype="GENEID")
```

**Exercise:**
For the genes in the example above, find the chromosome and strand
information that will go with each of the transcript names.

**Solution:**
```{r selectExercise}
columns(hg19_txdb)
cols <- c("TXNAME", "TXSTRAND", "TXCHROM")
select(hg19_txdb, keys=keys, columns=cols, keytype="GENEID")
```


## Methods for returning `GRanges` objects

Retrieving data with select is useful, but sometimes it is more
convenient to extract the result as a `GRanges` object.  This is
often the case when you are doing counting or specialized overlap
operations downstream.  For these use cases there is another family of
methods available.

Perhaps the most common operations for a `TxDb` object
is to retrieve the genomic coordinates or *ranges* for exons,
transcripts or coding sequences.  The functions
`transcripts`, `exons`, and `cds` return
the coordinate information as a `GRanges` object.

As an example, all transcripts present in a `TxDb` object
can be obtained as follows:

```{r transcripts1}
GR <- transcripts(hg19_txdb)
GR[1:3]
```

The `transcripts` function returns a `GRanges` class
object.  You can learn a lot more about the manipulation of these
objects by reading the `GenomicRanges` introductory
vignette.  The `show` method for a `GRanges` object
will display the ranges, seqnames (a chromosome or a contig), and
strand on the left side and then present related metadata on the right
side.

The `strand` function is used to obtain the strand
information from the transcripts.  The sum of the Lengths of the
`Rle` object that `strand` returns is equal to the
length of the `GRanges` object.

```{r transcripts2}
tx_strand <- strand(GR)
tx_strand
sum(runLength(tx_strand))
length(GR)
```

In addition, the `transcripts` function can also be used to
retrieve a subset of the transcripts available such as those on the
`+`-strand of chromosome 1.

```{r transcripts3}
GR <- transcripts(hg19_txdb, filter=list(tx_chrom = "chr15", tx_strand = "+"))
length(GR)
unique(strand(GR))
```

The `exons` and `cds` functions can also be used
in a similar fashion to retrive genomic coordinates for exons and
coding sequences.

The `promoters` function computes a `GRanges` object
that spans the promoter region around the transcription start site
for the transcripts in a `TxDb` object.  The `upstream`
and `downstream` arguments define the number of bases upstream
and downstream from the transcription start site that make up the
promoter region.

```{r transcripts4}
PR <- promoters(hg19_txdb, upstream=2000, downstream=400)
PR
```

A similar function (`terminators`) is provided to compute the
terminator region around the transcription end site for the
transcripts in a `TxDb` object.

**Exercise:**
Use `exons` to retrieve all the exons from chromosome 15.
How does the length of this compare to the value returned by
`transcripts`?

**Solution:**
```{r exonsExer1}
EX <- exons(hg19_txdb)
EX[1:4]
length(EX)
length(GR)
```


## Working with Grouped Features

Often one is interested in how particular genomic features relate to
each other, and not just their genomic positions.  For example, it might
be of interest to group transcripts by gene or to group exons by transcript.
Such groupings are supported by the `transcriptsBy`,
`exonsBy`, and `cdsBy` functions.

The following call can be used to group transcripts by genes:

```{r transcriptsBy}
GRList <- transcriptsBy(hg19_txdb, by = "gene")
length(GRList)
names(GRList)[10:13]
GRList[11:12]
```

The `transcriptsBy` function returns a `GRangesList`
class object.  As with `GRanges` objects, you can learn more
about these objects by reading the `GenomicRanges`
introductory vignette.  The `show` method for a
`GRangesList` object will display as a list of `GRanges`
objects.  And, at the bottom the seqinfo will be displayed once for
the entire list.

For each of these three functions, there is a limited set of options
that can be passed into the `by` argument to allow grouping.
For the `transcriptsBy` function, you can group by gene,
exon or cds, whereas for the `exonsBy` and `cdsBy`
functions can only be grouped by transcript (tx) or gene.

So as a further example, to extract all the exons for each transcript
you can call:

```{r exonsBy}
GRList <- exonsBy(hg19_txdb, by = "tx")
length(GRList)
names(GRList)[10:13]
GRList[[12]]
```

As you can see, the `GRangesList` objects returned from each
function contain genomic positions and identifiers grouped into a
list-like object according to the type of feature specified in the `by`
argument. The object returned can then be used by functions like
`findOverlaps` to contextualize alignments from
high-throughput sequencing.

The identifiers used to label the `GRanges` objects depend upon
the data source used to create the `TxDb` object.  So
the list identifiers will not always be Entrez Gene IDs, as they were
in the first example.  Furthermore, some data sources do not provide a
unique identifier for all features.  In this situation, the group
label will be a synthetic ID created by `GenomicFeatures` to
keep the relations between features consistent in the database this
was the case in the 2nd example.  Even though the results will
sometimes have to come back to you as synthetic IDs, you can still
always retrieve the original IDs.

**Exercise:**
Starting with the tx_ids that are the names of the GRList object we
just made, use `select` to retrieve that matching transcript
names.  Remember that the list used a `by` argument = "tx", so
the list is grouped by transcript IDs.

**Solution:**
```{r internalID}
GRList <- exonsBy(hg19_txdb, by = "tx")
tx_ids <- names(GRList)
head(select(hg19_txdb, keys=tx_ids, columns="TXNAME", keytype="TXID"))
```

Finally, the order of the results in a `GRangesList` object can
vary with the way in which things were grouped. In most cases the
grouped elements of the `GRangesList` object will be listed in
the order that they occurred along the chromosome.  However, when
exons or CDS parts are grouped by transcript, they will instead be
grouped according to their position along the transcript itself.
This is important because alternative splicing can mean that the
order along the transcript can be different from that along the
chromosome.


## Predefined grouping functions

The `intronsByTranscript`, `fiveUTRsByTranscript`
and `threeUTRsByTranscript` are convenience functions that
provide behavior equivalent to the grouping functions, but in
prespecified form. These functions return a `GRangesList`
object grouped by transcript for introns, 5' UTR's, and 3' UTR's,
respectively.  Below are examples of how you can call these methods.

```{r introns-UTRs}
length(intronsByTranscript(hg19_txdb))
length(fiveUTRsByTranscript(hg19_txdb))
length(threeUTRsByTranscript(hg19_txdb))
```



# Getting the actual sequence data

The `GenomicFeatures` package also provides functions for
converting from ranges to actual sequence (when paired
with an appropriate `BSgenome` package).

```{r extract}
suppressPackageStartupMessages(library(BSgenome.Hsapiens.UCSC.hg19))
genome <- BSgenome.Hsapiens.UCSC.hg19  # shorthand (for convenience)
tx_seqs1 <- extractTranscriptSeqs(genome, hg19_txdb, use.names=TRUE)
```

And, once these sequences have been extracted, you can translate them
into proteins with `translate`:

```{r translate1}
suppressWarnings(translate(tx_seqs1))
```

**Exercise:**
But of course this is not a meaningful translation, because the call
to `extractTranscriptSeqs` will have extracted all
the transcribed regions of the genome regardless of whether or not
they are translated. Look at the manual page for
`extractTranscriptSeqs` and see how you can use cdsBy
to only translate only the coding regions.

**Solution:**
```{r betterTranslation}
cds_seqs <- extractTranscriptSeqs(Hsapiens,
                                  cdsBy(hg19_txdb, by="tx", use.names=TRUE))
translate(cds_seqs)
```



# Session Information

```{r SessionInfo, echo=FALSE}
sessionInfo()
```