Retrieve gene/transcript/exons annotations stored in an Ensembl based database package generated with the makeEnsembldbPackage function. Parameter filter enables to define filters to retrieve only specific data. Alternatively, a global filter might be added to the EnsDb object using the addFilter method.

# S4 method for EnsDb
exons(x, columns = listColumns(x,"exon"),
        filter = AnnotationFilterList(), order.by,
        order.type = "asc", return.type = "GRanges")

# S4 method for EnsDb
exonsBy(x, by = c("tx", "gene"),
        columns = listColumns(x, "exon"), filter =
        AnnotationFilterList(), use.names = FALSE)

# S4 method for EnsDb
intronsByTranscript(x, ..., use.names = FALSE)

# S4 method for EnsDb
exonsByOverlaps(x, ranges, maxgap = -1L, minoverlap = 0L,
        type = c("any", "start", "end"), columns = listColumns(x, "exon"),
        filter = AnnotationFilterList())

# S4 method for EnsDb
transcripts(x, columns = listColumns(x, "tx"),
        filter = AnnotationFilterList(), order.by, order.type = "asc",
        return.type = "GRanges")

# S4 method for EnsDb
transcriptsBy(x, by = c("gene", "exon"),
        columns = listColumns(x, "tx"), filter = AnnotationFilterList())

# S4 method for EnsDb
transcriptsByOverlaps(x, ranges, maxgap = -1L,
        minoverlap = 0L, type = c("any", "start", "end"),
        columns = listColumns(x, "tx"), filter = AnnotationFilterList())

# S4 method for EnsDb
promoters(x, upstream = 2000, downstream = 200,
        use.names = TRUE, ...)

# S4 method for EnsDb
genes(x, columns = c(listColumns(x, "gene"), "entrezid"),
        filter = AnnotationFilterList(), order.by, order.type = "asc",
        return.type = "GRanges")

# S4 method for EnsDb
cdsBy(x, by = c("tx", "gene"), columns = NULL,
        filter = AnnotationFilterList(), use.names = FALSE)

# S4 method for EnsDb
fiveUTRsByTranscript(x, columns = NULL,
        filter = AnnotationFilterList())

# S4 method for EnsDb
threeUTRsByTranscript(x, columns = NULL,
        filter = AnnotationFilterList())

# S4 method for GRangesList
toSAF(x, ...)

Arguments

(In alphabetic order)

...

For promoters: additional arguments to be passed to the transcripts method. For intronsByTranscript: additional arguments such as filter.

by

For exonsBy: wheter exons sould be fetched by genes or by transcripts; as in the corresponding function of the GenomicFeatures package. For transcriptsBy: whether transcripts should be fetched by genes or by exons; fetching transcripts by cds as supported by the transcriptsBy method in the GenomicFeatures package is currently not implemented. For cdsBy: whether cds should be fetched by transcript of by gene.

columns

Columns to be retrieved from the database tables.

Default values for genes are all columns from the gene database table, for exons and exonsBy the column names of the exon database table table and for transcript and transcriptBy the columns of the tx data base table (see details below for more information).

Note that any of the column names of the database tables can be submitted to any of the methods (use listTables or listColumns methods for a complete list of allowed column names).

For cdsBy: this argument is only supported for for by="tx".

downstream

For method promoters: the number of nucleotides downstream of the transcription start site that should be included in the promoter region.

filter

A filter describing which results to retrieve from the database. Can be a single object extending AnnotationFilter, an AnnotationFilterList object combining several such objects or a formula representing a filter expression (see examples below or AnnotationFilter for more details). Use the supportedFilters method to get an overview of supported filter classes and related fields.

maxgap

For exonsByOverlaps and transcriptsByOverlaps: see exonsByOverlaps in GenomicFeatures for more information.

minoverlap

For exonsByOverlaps and transcriptsByOverlaps: see exonsByOverlaps in GenomicFeatures for more information.

order.by

Character vector specifying the column(s) by which the result should be ordered. This can be either in the form of "gene_id, seq_name" or c("gene_id", "seq_name").

order.type

If the results should be ordered ascending (asc, default) or descending (desc).

ranges

For exonsByOverlaps and transcriptsByOverlaps: a GRanges object specifying the genomic regions.

return.type

Type of the returned object. Can be either "data.frame", "DataFrame" or "GRanges". In the latter case the return object will be a GRanges object with the GRanges specifying the chromosomal start and end coordinates of the feature (gene, transcript or exon, depending whether genes, transcripts or exons was called). All additional columns are added as metadata columns to the GRanges object.

type

For exonsByOverlaps and transcriptsByOverlaps: see exonsByOverlaps in GenomicFeatures for more information.

upstream

For method promoters: the number of nucleotides upstream of the transcription start site that should be included in the promoter region.

use.names

For cdsBy and exonsBy: only for by="gene": use the names of the genes instead of their IDs as names of the resulting GRangesList.

x

For toSAF a GRangesList object. For all other methods an EnsDb instance.

Methods and Functions

Note that many methods and functions from the GenomicFeatures package can also be used for EnsDb objects (such as exonicParts, intronicParts etc).

exons

Retrieve exon information from the database. Additional columns from transcripts or genes associated with the exons can be specified and are added to the respective exon annotation.

exonsBy

Retrieve exons grouped by transcript or by gene. This function returns a GRangesList as does the analogous function in the GenomicFeatures package. Using the columns parameter it is possible to determine which additional values should be retrieved from the database. These will be included in the GRanges object for the exons as metadata columns. The exons in the inner GRanges are ordered by the exon index within the transcript (if by="tx"), or increasingly by the chromosomal start position of the exon or decreasingly by the chromosomal end position of the exon depending whether the gene is encoded on the + or - strand (for by="gene"). The GRanges in the GRangesList will be ordered by the name of the gene or transcript.

intronsByTranscript

Retrieve introns by transcripts. Filters can also be passed to the function. For more information see the intronsByTranscript method in the GenomicFeatures package.

exonsByOverlaps

Retrieve exons overlapping specified genomic ranges. For more information see the exonsByOverlaps method in the GenomicFeatures package. The functionality is to some extent similar and redundant to the exons method in combination with GRangesFilter filter.

transcripts

Retrieve transcript information from the database. Additional columns from genes or exons associated with the transcripts can be specified and are added to the respective transcript annotation.

transcriptsBy

Retrieve transcripts grouped by gene or exon. This function returns a GRangesList as does the analogous function in the GenomicFeatures package. Using the columns parameter it is possible to determine which additional values should be retrieved from the database. These will be included in the GRanges object for the transcripts as metadata columns. The transcripts in the inner GRanges are ordered increasingly by the chromosomal start position of the transcript for genes encoded on the + strand and in a decreasing manner by the chromosomal end position of the transcript for genes encoded on the - strand. The GRanges in the GRangesList will be ordered by the name of the gene or exon.

transcriptsByOverlaps

Retrieve transcripts overlapping specified genomic ranges. For more information see transcriptsByOverlaps method in the GenomicFeatures package. The functionality is to some extent similar and redundant to the transcripts method in combination with GRangesFilter filter.

promoters

Retrieve promoter information from the database. Additional columns from genes or exons associated with the promoters can be specified and are added to the respective promoter annotation.

genes

Retrieve gene information from the database. Additional columns from transcripts or exons associated with the genes can be specified and are added to the respective gene annotation. Note that column "entrezid" is a list of Entrezgene identifiers to accomodate the potential 1:n mapping between Ensembl genes and Entrezgene IDs.

cdsBy

Returns the coding region grouped either by transcript or by gene. Each element in the GRangesList represents the cds for one transcript or gene, with the individual ranges corresponding to the coding part of its exons. For by="tx" additional annotation columns can be added to the individual GRanges (in addition to the default columns exon_id and exon_rank). Note that the GRangesList is sorted by its names.

fiveUTRsByTranscript

Returns the 5' untranslated region for protein coding transcripts.

threeUTRsByTranscript

Returns the 3' untranslated region for protein coding transcripts.

toSAF

Reformats a GRangesList object into a data.frame corresponding to a standard SAF (Simplified Annotation Format) file (i.e. with column names "GeneID", "Chr", "Start", "End" and "Strand"). Note: this method makes only sense on a GRangesList that groups features (exons, transcripts) by gene.

Details

A detailed description of all database tables and the associated attributes/column names is also given in the vignette of this package. An overview of the columns is given below:

gene_id

the Ensembl gene ID of the gene.

gene_name

the name of the gene (in most cases its official symbol).

entrezid

the NCBI Entrezgene ID of the gene. Note that this column contains a list of Entrezgene identifiers to accommodate the potential 1:n mapping between Ensembl genes and Entrezgene IDs.

gene_biotype

the biotype of the gene.

gene_seq_start

the start coordinate of the gene on the sequence (usually a chromosome).

gene_seq_end

the end coordinate of the gene.

seq_name

the name of the sequence the gene is encoded (usually a chromosome).

seq_strand

the strand on which the gene is encoded

seq_coord_system

the coordinate system of the sequence.

tx_id

the Ensembl transcript ID.

tx_biotype

the biotype of the transcript.

tx_seq_start

the chromosomal start coordinate of the transcript.

tx_seq_end

the chromosomal end coordinate of the transcript.

tx_cds_seq_start

the start coordinate of the coding region of the transcript (NULL for non-coding transcripts).

tx_cds_seq_end

the end coordinate of the coding region.

gc_content

the G and C nucleotide content of the transcript's sequence expressed as a percentage (i.e. between 0 and 100).

exon_id

the ID of the exon. In Ensembl, each exon specified by a unique chromosomal start and end position has its own ID. Thus, the same exon might be part of several transcripts.

exon_seq_start

the chromosomal start coordinate of the exon.

exon_seq_end

the chromosomal end coordinate of the exon.

exon_idx

the index of the exon in the transcript model. As noted above, an exon can be part of several transcripts and thus its position inside these transcript might differ.

Many EnsDb databases provide also protein related annotations. See listProteinColumns for more information.

Note

Ensembl defines genes not only on standard chromosomes, but also on patched chromosomes and chromosome variants. Thus it might be advisable to restrict the queries to just those chromosomes of interest (e.g. by specifying a SeqNameFilter(c(1:22, "X", "Y"))). In addition, also so called LRG genes (Locus Reference Genomic) are defined in Ensembl. Their gene id starts with LRG instead of ENS for Ensembl genes, thus, a filter can be applied to specifically select those genes or exclude those genes (see examples below).

Depending on the value of the global option "ucscChromosomeNames" (use getOption(ucscChromosomeNames, FALSE) to get its value or option(ucscChromosomeNames=TRUE) to change its value) the sequence/chromosome names of the returned GRanges objects or provided in the returned data.frame or DataFrame correspond to Ensembl chromosome names (if value is FALSE) or UCSC chromosome names (if TRUE). This ensures a better integration with the Gviz package, in which this option is set by default to TRUE.

Value

For exons, transcripts and genes, a data.frame, DataFrame

or a GRanges, depending on the value of the

return.type parameter. The result is ordered as specified by the parameter order.by or, if not provided, by seq_name

and chromosomal start coordinate, but NOT by any ordering of values in eventually submitted filter objects.

For exonsBy, transcriptsBy: a GRangesList, depending on the value of the

return.type parameter. The results are ordered by the value of the

by parameter.

For exonsByOverlaps and transcriptsByOverlaps: a

GRanges with the exons or transcripts overlapping the specified regions.

For toSAF: a data.frame with column names

"GeneID" (the group name from the GRangesList, i.e. the ID by which the GRanges are split), "Chr" (the seqnames from the GRanges), "Start" (the start coordinate),

"End" (the end coordinate) and "Strand" (the strand).

For cdsBy: a GRangesList with GRanges per either transcript or exon specifying the start and end coordinates of the coding region of the transcript or gene.

For fiveUTRsByTranscript: a GRangesList with

GRanges for each protein coding transcript representing the start and end coordinates of full or partial exons that constitute the 5' untranslated region of the transcript.

For threeUTRsByTranscript: a GRangesList with

GRanges for each protein coding transcript representing the start and end coordinates of full or partial exons that constitute the 3' untranslated region of the transcript.

Note

While it is possible to request values from a column "tx_name" (with the columns argument), no such column is present in the database. The returned values correspond to the ID of the transcripts.

Author

Johannes Rainer, Tim Triche

See also

supportedFilters to get an overview of supported filters. makeEnsembldbPackage, listColumns, lengthOf

addFilter for globally adding filters to an EnsDb object.

Examples


library(EnsDb.Hsapiens.v86)
edb <- EnsDb.Hsapiens.v86

######   genes
##
## Get all genes encoded on chromosome Y
AllY <- genes(edb, filter = SeqNameFilter("Y"))
AllY
#> GRanges object with 523 ranges and 6 metadata columns:
#>                   seqnames            ranges strand |         gene_id
#>                      <Rle>         <IRanges>  <Rle> |     <character>
#>   ENSG00000251841        Y   2784749-2784853      + | ENSG00000251841
#>   ENSG00000184895        Y   2786855-2787699      - | ENSG00000184895
#>   ENSG00000237659        Y   2789827-2790328      + | ENSG00000237659
#>   ENSG00000232195        Y   2827982-2828218      + | ENSG00000232195
#>   ENSG00000129824        Y   2841486-2932000      + | ENSG00000129824
#>               ...      ...               ...    ... .             ...
#>   ENSG00000224240        Y 26549425-26549743      + | ENSG00000224240
#>   ENSG00000227629        Y 26586642-26591601      - | ENSG00000227629
#>   ENSG00000237917        Y 26594851-26634652      - | ENSG00000237917
#>   ENSG00000231514        Y 26626520-26627159      - | ENSG00000231514
#>   ENSG00000235857        Y 56855244-56855488      + | ENSG00000235857
#>                     gene_name           gene_biotype seq_coord_system
#>                   <character>            <character>      <character>
#>   ENSG00000251841  RNU6-1334P                  snRNA       chromosome
#>   ENSG00000184895         SRY         protein_coding       chromosome
#>   ENSG00000237659  RNASEH2CP1   processed_pseudogene       chromosome
#>   ENSG00000232195    TOMM22P2   processed_pseudogene       chromosome
#>   ENSG00000129824      RPS4Y1         protein_coding       chromosome
#>               ...         ...                    ...              ...
#>   ENSG00000224240     CYCSP49   processed_pseudogene       chromosome
#>   ENSG00000227629  SLC25A15P1 unprocessed_pseudogene       chromosome
#>   ENSG00000237917     PARP4P1 unprocessed_pseudogene       chromosome
#>   ENSG00000231514     FAM58CP   processed_pseudogene       chromosome
#>   ENSG00000235857     CTBP2P1   processed_pseudogene       chromosome
#>                        symbol entrezid
#>                   <character>   <list>
#>   ENSG00000251841  RNU6-1334P     <NA>
#>   ENSG00000184895         SRY     6736
#>   ENSG00000237659  RNASEH2CP1     <NA>
#>   ENSG00000232195    TOMM22P2     <NA>
#>   ENSG00000129824      RPS4Y1     6192
#>               ...         ...      ...
#>   ENSG00000224240     CYCSP49     <NA>
#>   ENSG00000227629  SLC25A15P1     <NA>
#>   ENSG00000237917     PARP4P1     <NA>
#>   ENSG00000231514     FAM58CP     <NA>
#>   ENSG00000235857     CTBP2P1     <NA>
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome

## Return the result as a DataFrame; also, we use a filter expression here
## to define which features to extract from the database.
AllY.granges <- genes(edb,
                      filter = ~ seq_name == "Y",
                      return.type="DataFrame")
AllY.granges
#> DataFrame with 523 rows and 10 columns
#>             gene_id   gene_name           gene_biotype gene_seq_start
#>         <character> <character>            <character>      <integer>
#> 1   ENSG00000251841  RNU6-1334P                  snRNA        2784749
#> 2   ENSG00000184895         SRY         protein_coding        2786855
#> 3   ENSG00000237659  RNASEH2CP1   processed_pseudogene        2789827
#> 4   ENSG00000232195    TOMM22P2   processed_pseudogene        2827982
#> 5   ENSG00000129824      RPS4Y1         protein_coding        2841486
#> ...             ...         ...                    ...            ...
#> 519 ENSG00000224240     CYCSP49   processed_pseudogene       26549425
#> 520 ENSG00000227629  SLC25A15P1 unprocessed_pseudogene       26586642
#> 521 ENSG00000237917     PARP4P1 unprocessed_pseudogene       26594851
#> 522 ENSG00000231514     FAM58CP   processed_pseudogene       26626520
#> 523 ENSG00000235857     CTBP2P1   processed_pseudogene       56855244
#>     gene_seq_end    seq_name seq_strand seq_coord_system      symbol entrezid
#>        <integer> <character>  <integer>      <character> <character>   <list>
#> 1        2784853           Y          1       chromosome  RNU6-1334P       NA
#> 2        2787699           Y         -1       chromosome         SRY     6736
#> 3        2790328           Y          1       chromosome  RNASEH2CP1       NA
#> 4        2828218           Y          1       chromosome    TOMM22P2       NA
#> 5        2932000           Y          1       chromosome      RPS4Y1     6192
#> ...          ...         ...        ...              ...         ...      ...
#> 519     26549743           Y          1       chromosome     CYCSP49       NA
#> 520     26591601           Y         -1       chromosome  SLC25A15P1       NA
#> 521     26634652           Y         -1       chromosome     PARP4P1       NA
#> 522     26627159           Y         -1       chromosome     FAM58CP       NA
#> 523     56855488           Y          1       chromosome     CTBP2P1       NA

## Include all transcripts of the gene and their chromosomal
## coordinates, sort by chrom start of transcripts and return as
## GRanges.
AllY.granges.tx <- genes(edb,
                         filter = SeqNameFilter("Y"),
                         columns = c("gene_id", "seq_name",
                                     "seq_strand", "tx_id", "tx_biotype",
                                     "tx_seq_start", "tx_seq_end"),
                         order.by = "tx_seq_start")
AllY.granges.tx
#> GRanges object with 740 ranges and 5 metadata columns:
#>                   seqnames            ranges strand |         gene_id
#>                      <Rle>         <IRanges>  <Rle> |     <character>
#>   ENSG00000251841        Y   2784749-2784853      + | ENSG00000251841
#>   ENSG00000184895        Y   2786855-2787699      - | ENSG00000184895
#>   ENSG00000237659        Y   2789827-2790328      + | ENSG00000237659
#>   ENSG00000232195        Y   2827982-2828218      + | ENSG00000232195
#>   ENSG00000129824        Y   2841486-2932000      + | ENSG00000129824
#>               ...      ...               ...    ... .             ...
#>   ENSG00000224240        Y 26549425-26549743      + | ENSG00000224240
#>   ENSG00000227629        Y 26586642-26591601      - | ENSG00000227629
#>   ENSG00000237917        Y 26594851-26634652      - | ENSG00000237917
#>   ENSG00000231514        Y 26626520-26627159      - | ENSG00000231514
#>   ENSG00000235857        Y 56855244-56855488      + | ENSG00000235857
#>                             tx_id             tx_biotype tx_seq_start
#>                       <character>            <character>    <integer>
#>   ENSG00000251841 ENST00000516032                  snRNA      2784749
#>   ENSG00000184895 ENST00000383070         protein_coding      2786855
#>   ENSG00000237659 ENST00000454281   processed_pseudogene      2789827
#>   ENSG00000232195 ENST00000430735   processed_pseudogene      2827982
#>   ENSG00000129824 ENST00000250784         protein_coding      2841486
#>               ...             ...                    ...          ...
#>   ENSG00000224240 ENST00000420810   processed_pseudogene     26549425
#>   ENSG00000227629 ENST00000456738 unprocessed_pseudogene     26586642
#>   ENSG00000237917 ENST00000435945 unprocessed_pseudogene     26594851
#>   ENSG00000231514 ENST00000435741   processed_pseudogene     26626520
#>   ENSG00000235857 ENST00000431853   processed_pseudogene     56855244
#>                   tx_seq_end
#>                    <integer>
#>   ENSG00000251841    2784853
#>   ENSG00000184895    2787699
#>   ENSG00000237659    2790328
#>   ENSG00000232195    2828218
#>   ENSG00000129824    2867268
#>               ...        ...
#>   ENSG00000224240   26549743
#>   ENSG00000227629   26591601
#>   ENSG00000237917   26634652
#>   ENSG00000231514   26627159
#>   ENSG00000235857   56855488
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome



######   transcripts
##
## Get all transcripts of a gene
Tx <- transcripts(edb,
                  filter = GeneIdFilter("ENSG00000184895"),
                  order.by = "tx_seq_start")
Tx
#> GRanges object with 1 range and 6 metadata columns:
#>                   seqnames          ranges strand |           tx_id
#>                      <Rle>       <IRanges>  <Rle> |     <character>
#>   ENST00000383070        Y 2786855-2787699      - | ENST00000383070
#>                       tx_biotype tx_cds_seq_start tx_cds_seq_end
#>                      <character>        <integer>      <integer>
#>   ENST00000383070 protein_coding          2786989        2787603
#>                           gene_id         tx_name
#>                       <character>     <character>
#>   ENST00000383070 ENSG00000184895 ENST00000383070
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome

## Get all transcripts of two genes along with some information on the
## gene and transcript
Tx <- transcripts(edb,
                  filter = GeneIdFilter(c("ENSG00000184895",
                                          "ENSG00000092377")),
                  columns = c("gene_id", "gene_seq_start", "gene_seq_end",
                              "gene_biotype", "tx_biotype"))
Tx
#> GRanges object with 4 ranges and 6 metadata columns:
#>                   seqnames          ranges strand |         gene_id
#>                      <Rle>       <IRanges>  <Rle> |     <character>
#>   ENST00000383070        Y 2786855-2787699      - | ENSG00000184895
#>   ENST00000383032        Y 6910686-7091683      + | ENSG00000092377
#>   ENST00000355162        Y 6910686-7091683      + | ENSG00000092377
#>   ENST00000346432        Y 6910686-7091683      + | ENSG00000092377
#>                   gene_seq_start gene_seq_end   gene_biotype     tx_biotype
#>                        <integer>    <integer>    <character>    <character>
#>   ENST00000383070        2786855      2787699 protein_coding protein_coding
#>   ENST00000383032        6910686      7091683 protein_coding protein_coding
#>   ENST00000355162        6910686      7091683 protein_coding protein_coding
#>   ENST00000346432        6910686      7091683 protein_coding protein_coding
#>                             tx_id
#>                       <character>
#>   ENST00000383070 ENST00000383070
#>   ENST00000383032 ENST00000383032
#>   ENST00000355162 ENST00000355162
#>   ENST00000346432 ENST00000346432
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome

######   promoters
##
## Get the bona-fide promoters (2k up- to 200nt downstream of TSS)
promoters(edb, filter = GeneIdFilter(c("ENSG00000184895",
                                       "ENSG00000092377")))
#> GRanges object with 4 ranges and 6 metadata columns:
#>                   seqnames          ranges strand |           tx_id
#>                      <Rle>       <IRanges>  <Rle> |     <character>
#>   ENST00000383070        Y 2787500-2789699      - | ENST00000383070
#>   ENST00000383032        Y 6908686-6910885      + | ENST00000383032
#>   ENST00000355162        Y 6908686-6910885      + | ENST00000355162
#>   ENST00000346432        Y 6908686-6910885      + | ENST00000346432
#>                       tx_biotype tx_cds_seq_start tx_cds_seq_end
#>                      <character>        <integer>      <integer>
#>   ENST00000383070 protein_coding          2786989        2787603
#>   ENST00000383032 protein_coding          7025085        7091492
#>   ENST00000355162 protein_coding          7025085        7091492
#>   ENST00000346432 protein_coding          7025085        7091492
#>                           gene_id         tx_name
#>                       <character>     <character>
#>   ENST00000383070 ENSG00000184895 ENST00000383070
#>   ENST00000383032 ENSG00000092377 ENST00000383032
#>   ENST00000355162 ENSG00000092377 ENST00000355162
#>   ENST00000346432 ENSG00000092377 ENST00000346432
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome

######   exons
##
## Get all exons of protein coding transcript for the gene ENSG00000184895
Exon <- exons(edb,
              filter = ~ gene_id == "ENSG00000184895" &
                  tx_biotype == "protein_coding",
              columns = c("gene_id", "gene_seq_start", "gene_seq_end",
                          "tx_biotype", "gene_biotype"))
Exon
#> GRanges object with 1 range and 6 metadata columns:
#>                   seqnames          ranges strand |         gene_id
#>                      <Rle>       <IRanges>  <Rle> |     <character>
#>   ENSE00001494622        Y 2786855-2787699      - | ENSG00000184895
#>                   gene_seq_start gene_seq_end     tx_biotype   gene_biotype
#>                        <integer>    <integer>    <character>    <character>
#>   ENSE00001494622        2786855      2787699 protein_coding protein_coding
#>                           exon_id
#>                       <character>
#>   ENSE00001494622 ENSE00001494622
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome



#####    exonsBy
##
## Get all exons for transcripts encoded on chromosomes X and Y.
ETx <- exonsBy(edb, by = "tx",
               filter = SeqNameFilter(c("X", "Y")))
ETx
#> GRangesList object of length 6829:
#> $ENST00000014935
#> GRanges object with 9 ranges and 2 metadata columns:
#>       seqnames              ranges strand |         exon_id exon_rank
#>          <Rle>           <IRanges>  <Rle> |     <character> <integer>
#>   [1]        X 154409418-154410039      - | ENSE00001050512         1
#>   [2]        X 154409112-154409196      - | ENSE00001050508         2
#>   [3]        X 154405434-154405655      - | ENSE00001450968         3
#>   [4]        X 154404995-154405083      - | ENSE00000678418         4
#>   [5]        X 154404828-154404914      - | ENSE00003790806         5
#>   [6]        X 154403522-154403622      - | ENSE00000678413         6
#>   [7]        X 154403269-154403381      - | ENSE00003549633         7
#>   [8]        X 154402942-154403190      - | ENSE00000678409         8
#>   [9]        X 154401758-154402841      - | ENSE00001200781         9
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> ...
#> <6828 more elements>
## Get all exons for genes encoded on chromosome 1 to 22, X and Y and
## include additional annotation columns in the result
EGenes <- exonsBy(edb, by = "gene",
                  filter = SeqNameFilter(c("X", "Y")),
                  columns = c("gene_biotype", "gene_name"))
EGenes
#> GRangesList object of length 2922:
#> $ENSG00000000003
#> GRanges object with 20 ranges and 3 metadata columns:
#>        seqnames              ranges strand |   gene_biotype   gene_name
#>           <Rle>           <IRanges>  <Rle> |    <character> <character>
#>    [1]        X 100639945-100639991      - | protein_coding      TSPAN6
#>    [2]        X 100636793-100637104      - | protein_coding      TSPAN6
#>    [3]        X 100636608-100636806      - | protein_coding      TSPAN6
#>    [4]        X 100636191-100636689      - | protein_coding      TSPAN6
#>    [5]        X 100635558-100635746      - | protein_coding      TSPAN6
#>    ...      ...                 ...    ... .            ...         ...
#>   [16]        X 100632541-100632568      - | protein_coding      TSPAN6
#>   [17]        X 100632063-100632068      - | protein_coding      TSPAN6
#>   [18]        X 100630759-100630866      - | protein_coding      TSPAN6
#>   [19]        X 100628670-100629986      - | protein_coding      TSPAN6
#>   [20]        X 100627109-100629986      - | protein_coding      TSPAN6
#>                exon_id
#>            <character>
#>    [1] ENSE00001828996
#>    [2] ENSE00001863395
#>    [3] ENSE00001855382
#>    [4] ENSE00001886883
#>    [5] ENSE00003662440
#>    ...             ...
#>   [16] ENSE00001849132
#>   [17] ENSE00003731560
#>   [18] ENSE00000868868
#>   [19] ENSE00001459322
#>   [20] ENSE00003730948
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> ...
#> <2921 more elements>

## Note that this might also contain "LRG" genes.
length(grep(names(EGenes), pattern="LRG"))
#> [1] 40

## to fetch just Ensemblgenes, use an GeneIdFilter with value
## "ENS%" and condition "like"
eg <- exonsBy(edb, by = "gene",
              filter = AnnotationFilterList(SeqNameFilter(c("X", "Y")),
                                            GeneIdFilter("ENS", "startsWith")),
              columns = c("gene_biotype", "gene_name"))
eg
#> GRangesList object of length 2882:
#> $ENSG00000000003
#> GRanges object with 20 ranges and 4 metadata columns:
#>        seqnames              ranges strand |   gene_biotype   gene_name
#>           <Rle>           <IRanges>  <Rle> |    <character> <character>
#>    [1]        X 100639945-100639991      - | protein_coding      TSPAN6
#>    [2]        X 100636793-100637104      - | protein_coding      TSPAN6
#>    [3]        X 100636608-100636806      - | protein_coding      TSPAN6
#>    [4]        X 100636191-100636689      - | protein_coding      TSPAN6
#>    [5]        X 100635558-100635746      - | protein_coding      TSPAN6
#>    ...      ...                 ...    ... .            ...         ...
#>   [16]        X 100632541-100632568      - | protein_coding      TSPAN6
#>   [17]        X 100632063-100632068      - | protein_coding      TSPAN6
#>   [18]        X 100630759-100630866      - | protein_coding      TSPAN6
#>   [19]        X 100628670-100629986      - | protein_coding      TSPAN6
#>   [20]        X 100627109-100629986      - | protein_coding      TSPAN6
#>                exon_id         gene_id
#>            <character>     <character>
#>    [1] ENSE00001828996 ENSG00000000003
#>    [2] ENSE00001863395 ENSG00000000003
#>    [3] ENSE00001855382 ENSG00000000003
#>    [4] ENSE00001886883 ENSG00000000003
#>    [5] ENSE00003662440 ENSG00000000003
#>    ...             ...             ...
#>   [16] ENSE00001849132 ENSG00000000003
#>   [17] ENSE00003731560 ENSG00000000003
#>   [18] ENSE00000868868 ENSG00000000003
#>   [19] ENSE00001459322 ENSG00000000003
#>   [20] ENSE00003730948 ENSG00000000003
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> ...
#> <2881 more elements>
length(grep(names(eg), pattern="LRG"))
#> [1] 0

#####    transcriptsBy
##
TGenes <- transcriptsBy(edb, by = "gene",
                        filter = SeqNameFilter(c("X", "Y")))
TGenes
#> GRangesList object of length 2922:
#> $ENSG00000000003
#> GRanges object with 5 ranges and 6 metadata columns:
#>       seqnames              ranges strand |           tx_id
#>          <Rle>           <IRanges>  <Rle> |     <character>
#>   [1]        X 100633442-100639991      - | ENST00000494424
#>   [2]        X 100627109-100637104      - | ENST00000612152
#>   [3]        X 100632063-100637104      - | ENST00000614008
#>   [4]        X 100628670-100636806      - | ENST00000373020
#>   [5]        X 100632541-100636689      - | ENST00000496771
#>                 tx_biotype tx_cds_seq_start tx_cds_seq_end         gene_id
#>                <character>        <integer>      <integer>     <character>
#>   [1] processed_transcript             <NA>           <NA> ENSG00000000003
#>   [2]       protein_coding        100630798      100635569 ENSG00000000003
#>   [3]       protein_coding        100632063      100635569 ENSG00000000003
#>   [4]       protein_coding        100630798      100636694 ENSG00000000003
#>   [5] processed_transcript             <NA>           <NA> ENSG00000000003
#>               tx_name
#>           <character>
#>   [1] ENST00000494424
#>   [2] ENST00000612152
#>   [3] ENST00000614008
#>   [4] ENST00000373020
#>   [5] ENST00000496771
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> $ENSG00000000005
#> GRanges object with 2 ranges and 6 metadata columns:
#>       seqnames              ranges strand |           tx_id
#>          <Rle>           <IRanges>  <Rle> |     <character>
#>   [1]        X 100584802-100599885      + | ENST00000373031
#>   [2]        X 100593624-100597531      + | ENST00000485971
#>                 tx_biotype tx_cds_seq_start tx_cds_seq_end         gene_id
#>                <character>        <integer>      <integer>     <character>
#>   [1]       protein_coding        100585019      100599717 ENSG00000000005
#>   [2] processed_transcript             <NA>           <NA> ENSG00000000005
#>               tx_name
#>           <character>
#>   [1] ENST00000373031
#>   [2] ENST00000485971
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> $ENSG00000001497
#> GRanges object with 5 ranges and 6 metadata columns:
#>       seqnames            ranges strand |           tx_id
#>          <Rle>         <IRanges>  <Rle> |     <character>
#>   [1]        X 65512583-65534775      - | ENST00000484069
#>   [2]        X 65512582-65534756      - | ENST00000374811
#>   [3]        X 65512583-65534756      - | ENST00000374804
#>   [4]        X 65512582-65534754      - | ENST00000374807
#>   [5]        X 65520429-65523617      - | ENST00000469091
#>                   tx_biotype tx_cds_seq_start tx_cds_seq_end         gene_id
#>                  <character>        <integer>      <integer>     <character>
#>   [1] nonsense_mediated_de..         65525021       65534715 ENSG00000001497
#>   [2]         protein_coding         65512775       65534715 ENSG00000001497
#>   [3]         protein_coding         65512775       65534715 ENSG00000001497
#>   [4]         protein_coding         65512775       65534715 ENSG00000001497
#>   [5]         protein_coding         65520655       65523617 ENSG00000001497
#>               tx_name
#>           <character>
#>   [1] ENST00000484069
#>   [2] ENST00000374811
#>   [3] ENST00000374804
#>   [4] ENST00000374807
#>   [5] ENST00000469091
#>   -------
#>   seqinfo: 2 sequences from GRCh38 genome
#> 
#> ...
#> <2919 more elements>

## convert this to a SAF formatted data.frame that can be used by the
## featureCounts function from the Rsubreader package.
head(toSAF(TGenes))
#>            GeneID Chr     Start       End Strand
#> 1 ENSG00000000003   X 100633442 100639991      -
#> 2 ENSG00000000003   X 100627109 100637104      -
#> 3 ENSG00000000003   X 100632063 100637104      -
#> 4 ENSG00000000003   X 100628670 100636806      -
#> 5 ENSG00000000003   X 100632541 100636689      -
#> 6 ENSG00000000005   X 100584802 100599885      +


#####   transcriptsByOverlaps
##
ir <- IRanges(start = c(2654890, 2709520, 28111770),
              end = c(2654900, 2709550, 28111790))
gr <- GRanges(rep("Y", length(ir)), ir)

## Retrieve all transcripts overlapping any of the regions.
txs <- transcriptsByOverlaps(edb, gr)
txs
#> GRanges object with 0 ranges and 6 metadata columns:
#>    seqnames    ranges strand |       tx_id  tx_biotype tx_cds_seq_start
#>       <Rle> <IRanges>  <Rle> | <character> <character>        <integer>
#>    tx_cds_seq_end     gene_id     tx_name
#>         <integer> <character> <character>
#>   -------
#>   seqinfo: 1 sequence from GRCh38 genome

## Alternatively, use a GRangesFilter
grf <- GRangesFilter(gr, type = "any")
txs <- transcripts(edb, filter = grf)
txs
#> GRanges object with 0 ranges and 6 metadata columns:
#>    seqnames    ranges strand |       tx_id  tx_biotype tx_cds_seq_start
#>       <Rle> <IRanges>  <Rle> | <character> <character>        <integer>
#>    tx_cds_seq_end     gene_id     tx_name
#>         <integer> <character> <character>
#>   -------
#>   seqinfo: no sequences


####    cdsBy
## Get the coding region for all transcripts on chromosome Y.
## Specifying also additional annotation columns (in addition to the default
## exon_id and exon_rank).
cds <- cdsBy(edb, by = "tx", filter = SeqNameFilter("Y"),
             columns = c("tx_biotype", "gene_name"))

####    the 5' untranslated regions:
fUTRs <- fiveUTRsByTranscript(edb, filter = SeqNameFilter("Y"))

####    the 3' untranslated regions with additional column gene_name.
tUTRs <- threeUTRsByTranscript(edb, filter = SeqNameFilter("Y"),
                               columns = "gene_name")