makeGRangesFromGFF: Make GRanges from a GFF/GTF file

Description Usage Arguments Details Value GFF/GTF specification Feature type Supported sources Ensembl GENCODE GENCODE vs. Ensembl RefSeq UCSC FlyBase WormBase Note See Also Examples

View source: R/makeGRangesFromGFF.R

Description

Make GRanges from a GFF/GTF file

Usage

1
2
3
4
5
6
makeGRangesFromGFF(
  file,
  level = c("genes", "transcripts"),
  ignoreVersion = TRUE,
  synonyms = FALSE
)

Arguments

file

character(1). File path.

level

character(1). Return as genes or transcripts.

ignoreVersion

logical(1). Ignore identifier (e.g. transcript, gene) versions. When applicable, the identifier containing version numbers will be stored in txIdVersion and geneIdVersion, and the variants without versions will be stored in txId, txIdNoVersion, geneId, and geneIdNoVersion.

synonyms

logical(1). Include gene synonyms. Queries the Ensembl web server, and is CPU intensive.

Details

Remote URLs and compressed files are supported.

Value

GRanges.

GFF/GTF specification

The GFF (General Feature Format) format consists of one line per feature, each containing 9 columns of data, plus optional track definition lines.

The GTF (General Transfer Format) format is identical to GFF version 2.

The UCSC website has detailed conventions on the GFF3 format, including the metadata columns.

Feature type

See also:

Supported sources

Currently makeGRangesFromGFF() supports genomes from these sources:

Ensembl

Note that makeGRangesFromEnsembl() offers native support for Ensembl genome builds and returns additional useful metadata that isn't defined inside a GFF/GTF file.

If you must load a GFF/GTF file directly, then use makeGRangesFromGFF().

Example URLs:

GENCODE

Example URLs:

GENCODE vs. Ensembl

Annotations available from Ensembl and GENCODE are very similar.

The GENCODE annotation is made by merging the manual gene annotation produced by the Ensembl-Havana team and the Ensembl-genebuild automated gene annotation. The GENCODE annotation is the default gene annotation displayed in the Ensembl browser. The GENCODE releases coincide with the Ensembl releases, although GENCODE can skip an Ensembl release if there is no update to the annotation with respect to the previous release. In practical terms, the GENCODE annotation is essentially identical to the Ensembl annotation.

However, GENCODE handles pseudoautosomal regions (PAR) differently than Ensembl. The Ensembl GTF file only includes this annotation once, for chromosome X. However, GENCODE GTF/GFF3 files include the annotation in the PAR regions of both chromosomes. You'll see these genes contain a "_PAR_Y" suffix.

Additionally, GENCODE GFF/GTF files import with a gene identifier containing a suffix, which differs slightly from the Ensembl GFF/GTF spec (e.g. GENCODE: ENSG00000000003.14; Ensembl: ENSG00000000003).

The GENCODE FAQ has additional details.

RefSeq

Refer to the current RefSeq spec for details.

Example URLs:

See also:

UCSC

Example URLs:

Related URLs:

FlyBase

Example URLs:

WormBase

Example URLs:

Note

Updated 2021-08-06.

See Also

Examples

 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
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
## Some examples here are commented because they are CPU-intensive and
## can cause CI timeouts.

## Ensembl ====
file <- pasteURL(
    "ftp.ensembl.org",
    "pub",
    "release-102",
    "gtf",
    "homo_sapiens",
    "Homo_sapiens.GRCh38.102.gtf.gz",
    protocol = "ftp"
)
genes <- makeGRangesFromGFF(
    file = file,
    level = "genes",
    ignoreVersion = FALSE
)
summary(genes)
## > transcripts <- makeGRangesFromGFF(
## >     file = file,
## >     level = "transcripts",
## >     ignoreVersion = FALSE
## > )
## > summary(transcripts)

## GENCODE ====
## > file <- pasteURL(
## >     "ftp.ebi.ac.uk",
## >     "pub",
## >     "databases",
## >     "gencode",
## >     "Gencode_human",
## >     "release_36",
## >     "gencode.v36.annotation.gtf.gz",
## >     protocol = "ftp"
## > )
## > genes <- makeGRangesFromGFF(file = file, level = "genes")
## > summary(genes)
## > transcripts <- makeGRangesFromGFF(file = file, level = "transcripts")
## > summary(transcripts)

## RefSeq ====
## > file <- pasteURL(
## >     "ftp.ncbi.nlm.nih.gov",
## >     "genomes",
## >     "refseq",
## >     "vertebrate_mammalian",
## >     "Homo_sapiens",
## >     "all_assembly_versions",
## >     "GCF_000001405.39_GRCh38.p13",
## >     "GCF_000001405.39_GRCh38.p13_genomic.gff.gz",
## >     protocol = "ftp"
## > )
## > genes <- makeGRangesFromGFF(file = file, level = "genes")
## > summary(genes)
## > transcripts <- makeGRangesFromGFF(file = file, level = "transcripts")
## > summary(transcripts)

## UCSC ====
## > file <- pasteURL(
## >     "hgdownload.soe.ucsc.edu",
## >     "goldenPath",
## >     "hg38",
## >     "bigZips",
## >     "genes",
## >     "hg38.ensGene.gtf.gz",
## >     protocol = "ftp"
## > )
## > genes <- makeGRangesFromGFF(file = file, level = "genes")
## > summary(genes)
## > transcripts <- makeGRangesFromGFF(file = file, level = "transcripts")
## > summary(transcripts)

acidgenomics/AcidGenomes documentation built on Sept. 16, 2021, 7:30 p.m.