makeGRangesFromGff: Make genomic ranges ('GRanges') from a GFF/GTF file

View source: R/makeGRangesFromGff.R

makeGRangesFromGffR Documentation

Make genomic ranges (GRanges) from a GFF/GTF file

Description

Make genomic ranges (GRanges) from a GFF/GTF file

Usage

makeGRangesFromGff(
  file,
  level = c("genes", "transcripts"),
  ignoreVersion = FALSE,
  extraMcols = TRUE
)

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.

extraMcols

logical(1). Include extra metadata columns (e.g. "broadClass").

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

  • CDS: CoDing Ssequence. A contiguous sequence that contains a genomic interval bounded by start and stop codons. CDS refers to the portion of a genomic DNA sequence that is translated, from the start codon to the stop codon.

  • exon: Genomic interval containing 5' UTR (five_prime_UTR), CDS, and 3' UTR (three_prime_UTR).

  • mRNA: Processed (spliced) mRNA transcript.

See also:

Supported sources

Currently makeGRangesFromGff() supports genomes from these sources:

  • Ensembl

  • GENCODE

  • RefSeq

  • UCSC

  • FlyBase

  • WormBase

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:

  • Ensembl Homo sapiens GRCh38.p13, release 108: GTF, GFF3

  • Ensembl Homo sapiens GRCh37, release 108 (87): GTF, GFF3

GENCODE

Example URLs:

  • GENCODE Homo sapiens GRCh38.p13, release 42: GTF, GFF3

  • GENCODE Homo sapiens GRCh37, release 42: GTF, GFF3

  • GENCODE Mus musculus GRCm39, release M31: GTF, GFF3

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:

  • RefSeq Homo sapiens GRCh38.p14 GTF, GFF3

See also:

  • RefSeq FAQ

  • https://ftp.ncbi.nih.gov/gene/DATA/gene2refseq.gz

UCSC

Example URLs:

Related URLs:

FlyBase

Example URLs:

  • FlyBase Drosophila melanogaster r6.49 GTF, GFF3

WormBase

Example URLs:

  • WormBase Caenorhabditis elegans WS287 GTF, GFF3

Note

Updated 2023-12-04.

See Also

  • rtracklayer::import().

  • GenomicFeatures::makeTxDbFromGRanges().

  • GenomicFeatures::makeTxDbFromGFF().

  • GenomicFeatures:::.make_splicings().

  • tximeta:::getRanges().

  • AnnotationDbi::columns().

  • GenomeInfoDb::GenomeDescription-class, which describes useful organism, commonName, providerVersion, provider, and releaseDate accessors.

Examples

## Ensembl ====
file <- AcidBase::pasteUrl(
    "ftp.ensembl.org",
    "pub",
    "release-108",
    "gtf",
    "homo_sapiens",
    "Homo_sapiens.GRCh38.108.gtf.gz",
    protocol = "ftp"
)
genes <- makeGRangesFromGff(file = file, level = "genes")
summary(genes)
## > transcripts <- makeGRangesFromGff(file = file, level = "transcripts")
## > summary(transcripts)

## GENCODE ====
## > file <- AcidBase::pasteUrl(
## >     "ftp.ebi.ac.uk",
## >     "pub",
## >     "databases",
## >     "gencode",
## >     "Gencode_human",
## >     "release_42",
## >     "gencode.v42.annotation.gtf.gz",
## >     protocol = "ftp"
## > )
## > genes <- makeGRangesFromGff(file = file, level = "genes")
## > summary(genes)
## > transcripts <- makeGRangesFromGff(file = file, level = "transcripts")
## > summary(transcripts)

## RefSeq ====
## > file <- AcidBase::pasteUrl(
## >     "ftp.ncbi.nlm.nih.gov",
## >     "genomes",
## >     "refseq",
## >     "vertebrate_mammalian",
## >     "Homo_sapiens",
## >     "all_assembly_versions",
## >     "GCF_000001405.40_GRCh38.p14",
## >     "GCF_000001405.40_GRCh38.p14_genomic.gff.gz",
## >     protocol = "ftp"
## > )
## > genes <- makeGRangesFromGff(file = file, level = "genes")
## > summary(genes)
## > transcripts <- makeGRangesFromGff(file = file, level = "transcripts")
## > summary(transcripts)

## UCSC ====
## > file <- AcidBase::pasteUrl(
## >     "hgdownload.soe.ucsc.edu",
## >     "goldenPath",
## >     "hg38",
## >     "bigZips",
## >     "genes",
## >     "hg38.knownGene.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 Dec. 10, 2023, 10:35 p.m.