The Bioc.gff package provides support for the General Feature Format (GFF) file format, which is widely used for representing genomic features and annotations. This package allows users to read, write, and manipulate GFF versions 1, 2 (GTF), and 3 in R, making it easier to work with genomic data.
While the rtracklayer
package offers robust GFF support, it is a large package with many
features beyond file import. Bioc.gff fills a specific
niche in the Bioconductor ecosystem by providing a lightweight, focused
solution with minimal dependencies. This modularity is beneficial for
developers of other packages who need to parse GFF files without
inheriting the extensive dependency footprint of
rtracklayer. For the end user, it offers a simple and
direct interface for GFF manipulation.
Note that much of the code in the package is ported and adapted from
the rtracklayer
Bioconductor package. The intention is that the GFF functionality
residing in rtracklayer
will be removed from that package in favor of using
Bioc.gff, thereby reducing its size and complexity.
You can install the Bioc.gff package from Bioconductor
using the following:
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> Loading required package: generics
#>
#> Attaching package: 'generics'
#> The following objects are masked from 'package:base':
#>
#> as.difftime, as.factor, as.ordered, intersect, is.element, setdiff,
#> setequal, union
#>
#> Attaching package: 'BiocGenerics'
#> The following objects are masked from 'package:stats':
#>
#> IQR, mad, sd, var, xtabs
#> The following objects are masked from 'package:base':
#>
#> Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#> as.data.frame, basename, cbind, colnames, dirname, do.call,
#> duplicated, eval, evalq, get, grep, grepl, is.unsorted, lapply,
#> mapply, match, mget, order, paste, pmax, pmax.int, pmin, pmin.int,
#> rank, rbind, rownames, sapply, saveRDS, table, tapply, unique,
#> unsplit, which.max, which.min
#> Loading required package: S4Vectors
#>
#> Attaching package: 'S4Vectors'
#> The following object is masked from 'package:utils':
#>
#> findMatches
#> The following objects are masked from 'package:base':
#>
#> I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: Seqinfo
The Bioc.gff package provides functions to read and
write GFF files, as well as to manipulate GFF data structures. The
following sections provide some examples of how to use the package.
You can read GFF files using the readGFF function. This
function supports GFF versions 1, 2, and 3. The following example
demonstrates how to read a GFF file:
import functionThe import function deduces that the input string is a
GFF file type and calls the appropriate method for reading the GFF
file.
#> GRanges object with 31 ranges and 10 metadata columns:
#> seqnames ranges strand | source type score phase
#> <Rle> <IRanges> <Rle> | <factor> <factor> <numeric> <integer>
#> [1] chr10 92828-95504 - | rtracklayer gene 5 <NA>
#> [2] chr10 92828-95178 - | rtracklayer mRNA NA <NA>
#> [3] chr10 92828-95504 - | rtracklayer mRNA NA <NA>
#> [4] chr10 92828-94054 - | rtracklayer exon NA <NA>
#> [5] chr10 92997-94054 - | rtracklayer CDS NA <NA>
#> ... ... ... ... . ... ... ... ...
#> [27] chr12 89675-89827 + | rtracklayer CDS NA <NA>
#> [28] chr12 90587-90655 + | rtracklayer exon NA <NA>
#> [29] chr12 90587-90655 + | rtracklayer CDS NA <NA>
#> [30] chr12 90796-91263 + | rtracklayer exon NA <NA>
#> [31] chr12 90796-91263 * | rtracklayer CDS NA <NA>
#> ID Name geneName Alias genome
#> <character> <character> <character> <list> <character>
#> [1] GeneID:347688 TUBB8 tubulin, beta 8 FLJ40100,TUBB8 hg19
#> [2] 873 TUBB8 <NA> <NA>
#> [3] 872 TUBB8 <NA> <NA>
#> [4] <NA> <NA> <NA> <NA>
#> [5] <NA> <NA> <NA> <NA>
#> ... ... ... ... ... ...
#> [27] <NA> <NA> <NA> <NA>
#> [28] <NA> <NA> <NA> <NA>
#> [29] <NA> <NA> <NA> <NA>
#> [30] <NA> <NA> <NA> <NA>
#> [31] <NA> <NA> <NA> <NA>
#> Parent
#> <list>
#> [1]
#> [2] GeneID:347688
#> [3] GeneID:347688
#> [4] 872,873
#> [5] 872,873
#> ... ...
#> [27] 4644
#> [28] 4644
#> [29] 4644
#> [30] 4644
#> [31] 4644
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
Here the output object is a GRanges class, which
contains genomic coordinates and associated metadata. The
GRanges object indicates the ranges with the
seqnames, ranges, and strand,
columns. The associated metadata is stored in the mcols
slot, which (in this example) includes attributes like
source, type, phase,
ID, Name, geneName,
Alias and genome.
You can also selectively import ranges from the GFF file using the
which argument. This allows you to filter the data based on
specific criteria, such as chromosome or strand. For example, to import
only the ranges on chromosome 1:
#> GRanges object with 5 ranges and 10 metadata columns:
#> seqnames ranges strand | source type score phase
#> <Rle> <IRanges> <Rle> | <factor> <factor> <numeric> <integer>
#> [1] chr10 92828-95504 - | rtracklayer gene 5 <NA>
#> [2] chr10 92828-95178 - | rtracklayer mRNA NA <NA>
#> [3] chr10 92828-95504 - | rtracklayer mRNA NA <NA>
#> [4] chr10 92828-94054 - | rtracklayer exon NA <NA>
#> [5] chr10 92997-94054 - | rtracklayer CDS NA <NA>
#> ID Name geneName Alias genome
#> <character> <character> <character> <list> <character>
#> [1] GeneID:347688 TUBB8 tubulin, beta 8 FLJ40100,TUBB8 hg19
#> [2] 873 TUBB8 <NA> <NA>
#> [3] 872 TUBB8 <NA> <NA>
#> [4] <NA> <NA> <NA> <NA>
#> [5] <NA> <NA> <NA> <NA>
#> Parent
#> <list>
#> [1]
#> [2] GeneID:347688
#> [3] GeneID:347688
#> [4] 872,873
#> [5] 872,873
#> -------
#> seqinfo: 298 sequences (2 circular) from hg19 genome
Note that you can indicate the genome build using the
genome argument, which is useful for ensuring that the
imported data is compatible with other genomic data you may be working
with.
readGFF functionAs an alternative, you can use the readGFF function
directly, which is more explicit about the GFF version:
#> DataFrame with 31 rows and 14 columns
#> seqid source type start end score strand
#> <factor> <factor> <factor> <integer> <integer> <numeric> <character>
#> 1 chr10 rtracklayer gene 92828 95504 5 -
#> 2 chr10 rtracklayer mRNA 92828 95178 NA -
#> 3 chr10 rtracklayer mRNA 92828 95504 NA -
#> 4 chr10 rtracklayer exon 92828 94054 NA -
#> 5 chr10 rtracklayer CDS 92997 94054 NA -
#> ... ... ... ... ... ... ... ...
#> 27 chr12 rtracklayer CDS 89675 89827 NA +
#> 28 chr12 rtracklayer exon 90587 90655 NA +
#> 29 chr12 rtracklayer CDS 90587 90655 NA +
#> 30 chr12 rtracklayer exon 90796 91263 NA +
#> 31 chr12 rtracklayer CDS 90796 91263 NA *
#> phase ID Name geneName Alias
#> <integer> <character> <character> <character> <list>
#> 1 NA GeneID:347688 TUBB8 tubulin, beta 8 FLJ40100,TUBB8
#> 2 NA 873 TUBB8 NA
#> 3 NA 872 TUBB8 NA
#> 4 NA NA NA NA
#> 5 NA NA NA NA
#> ... ... ... ... ... ...
#> 27 NA NA NA NA
#> 28 NA NA NA NA
#> 29 NA NA NA NA
#> 30 NA NA NA NA
#> 31 NA NA NA NA
#> genome Parent
#> <character> <list>
#> 1 hg19
#> 2 NA GeneID:347688
#> 3 NA GeneID:347688
#> 4 NA 872,873
#> 5 NA 872,873
#> ... ... ...
#> 27 NA 4644
#> 28 NA 4644
#> 29 NA 4644
#> 30 NA 4644
#> 31 NA 4644
This function returns a DataFrame object, which contains
the same information as the GRanges object but in a tabular
format. Note that one can use the makeGRangesFromDataFrame
function to convert the DataFrame returned by
readGFF into a GRanges object, which is often
more convenient for downstream analysis:
#> GRanges object with 31 ranges and 10 metadata columns:
#> seqnames ranges strand | source type score phase
#> <Rle> <IRanges> <Rle> | <factor> <factor> <numeric> <integer>
#> [1] chr10 92828-95504 - | rtracklayer gene 5 <NA>
#> [2] chr10 92828-95178 - | rtracklayer mRNA NA <NA>
#> [3] chr10 92828-95504 - | rtracklayer mRNA NA <NA>
#> [4] chr10 92828-94054 - | rtracklayer exon NA <NA>
#> [5] chr10 92997-94054 - | rtracklayer CDS NA <NA>
#> ... ... ... ... . ... ... ... ...
#> [27] chr12 89675-89827 + | rtracklayer CDS NA <NA>
#> [28] chr12 90587-90655 + | rtracklayer exon NA <NA>
#> [29] chr12 90587-90655 + | rtracklayer CDS NA <NA>
#> [30] chr12 90796-91263 + | rtracklayer exon NA <NA>
#> [31] chr12 90796-91263 * | rtracklayer CDS NA <NA>
#> ID Name geneName Alias genome
#> <character> <character> <character> <list> <character>
#> [1] GeneID:347688 TUBB8 tubulin, beta 8 FLJ40100,TUBB8 hg19
#> [2] 873 TUBB8 <NA> <NA>
#> [3] 872 TUBB8 <NA> <NA>
#> [4] <NA> <NA> <NA> <NA>
#> [5] <NA> <NA> <NA> <NA>
#> ... ... ... ... ... ...
#> [27] <NA> <NA> <NA> <NA>
#> [28] <NA> <NA> <NA> <NA>
#> [29] <NA> <NA> <NA> <NA>
#> [30] <NA> <NA> <NA> <NA>
#> [31] <NA> <NA> <NA> <NA>
#> Parent
#> <list>
#> [1]
#> [2] GeneID:347688
#> [3] GeneID:347688
#> [4] 872,873
#> [5] 872,873
#> ... ...
#> [27] 4644
#> [28] 4644
#> [29] 4644
#> [30] 4644
#> [31] 4644
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
You can also read GFF files from remote URLs. The import
function can handle URLs directly, allowing you to read GFF files hosted
on remote servers. In this example, we read a GFF3 file from the miRBase
database:
To show the example output in the vignette, a more advanced approach
is used below. This is done to avoid repeated downloads of the remote
file when the vignette is re-built on the Bioconductor Build System
(BBS). We use the BiocFileCache package to cache the file
locally.
#> Loading required package: dbplyr
bfc <- BiocFileCache::BiocFileCache()
remote_gff_url <- "https://www.mirbase.org/download/hsa.gff3"
bquery <- bfcquery(bfc, remote_gff_url, "rname", exact = TRUE)
if (!nrow(bquery))
bfcadd(x = bfc, rname = remote_gff_url, rtype = "web", download = TRUE)#> BFC1
#> "/github/home/.cache/R/BiocFileCache/150711816bda_hsa.gff3"
gff_local <- bfcrpath(
bfc, rnames = remote_gff_url, exact = TRUE, download = FALSE, rtype = "web"
)Finally, the relevant function remains the same for reading a GFF
file i.e., via import():
#> GRanges object with 4801 ranges and 8 metadata columns:
#> seqnames ranges strand | source type
#> <Rle> <IRanges> <Rle> | <factor> <factor>
#> [1] chr1 17369-17436 - | NA miRNA_primary_transcript
#> [2] chr1 17409-17431 - | NA miRNA
#> [3] chr1 17369-17391 - | NA miRNA
#> [4] chr1 30366-30503 + | NA miRNA_primary_transcript
#> [5] chr1 30438-30458 + | NA miRNA
#> ... ... ... ... . ... ...
#> [4797] chrY 2609229-2609252 + | NA miRNA
#> [4798] chrY 4606120-4606228 + | NA miRNA_primary_transcript
#> [4799] chrY 4606140-4606158 + | NA miRNA
#> [4800] chrY 13479177-13479266 + | NA miRNA_primary_transcript
#> [4801] chrY 13479189-13479214 + | NA miRNA
#> score phase ID Alias Name
#> <numeric> <integer> <character> <list> <character>
#> [1] NA <NA> MI0022705 MI0022705 hsa-mir-6859-1
#> [2] NA <NA> MIMAT0027618 MIMAT0027618 hsa-miR-6859-5p
#> [3] NA <NA> MIMAT0027619 MIMAT0027619 hsa-miR-6859-3p
#> [4] NA <NA> MI0006363 MI0006363 hsa-mir-1302-2
#> [5] NA <NA> MIMAT0005890 MIMAT0005890 hsa-miR-1302
#> ... ... ... ... ... ...
#> [4797] NA <NA> MIMAT0023714_1 MIMAT0023714 hsa-miR-6089
#> [4798] NA <NA> MI0032313 MI0032313 hsa-mir-9985
#> [4799] NA <NA> MIMAT0039763 MIMAT0039763 hsa-miR-9985
#> [4800] NA <NA> MI0039722 MI0039722 hsa-mir-12120
#> [4801] NA <NA> MIMAT0049014 MIMAT0049014 hsa-miR-12120
#> Derives_from
#> <character>
#> [1] <NA>
#> [2] MI0022705
#> [3] MI0022705
#> [4] <NA>
#> [5] MI0006363
#> ... ...
#> [4797] MI0023563
#> [4798] <NA>
#> [4799] MI0032313
#> [4800] <NA>
#> [4801] MI0039722
#> -------
#> seqinfo: 24 sequences from an unspecified genome; no seqlengths
While TxDb objects are highly efficient for querying
transcript annotations within R, it is often necessary to export this
data to a standard, portable format for use with external tools or for
sharing with collaborators. The GFF3 format is a widely accepted
standard for this purpose. Converting a TxDb or a
derivative object (like a GRangesList of exons) into a
GFF3-compatible GRanges object allows for easy export. This
is particularly useful for visualizing annotations in genome browsers
like IGV or UCSC, or for input into downstream analysis pipelines that
expect GFF3 files.
#> Loading required package: AnnotationDbi
#> Loading required package: Biobase
#> Welcome to Bioconductor
#>
#> Vignettes contain introductory material; view with
#> 'browseVignettes()'. To cite Bioconductor, see
#> 'citation("Biobase")', and for packages 'citation("pkgname")'.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
exonsBy(txdb, by = "tx") |>
asGFF()#> GRanges object with 2517401 ranges and 7 metadata columns:
#> seqnames ranges strand | type ID Name
#> <Rle> <IRanges> <Rle> | <character> <character> <character>
#> [1] chr1 10370-10582 + | mRNA mRNA1 1
#> [2] chr1 11121-14413 + | mRNA mRNA2 2
#> [3] chr1 11125-14405 + | mRNA mRNA3 3
#> [4] chr1 11410-14413 + | mRNA mRNA4 4
#> [5] chr1 11411-14413 + | mRNA mRNA5 5
#> ... ... ... ... . ... ... ...
#> [2517397] chrMT 5826-5891 - | exon exon2135410 <NA>
#> [2517398] chrMT 7446-7514 - | exon exon2135411 <NA>
#> [2517399] chrMT 14149-14673 - | exon exon2135412 <NA>
#> [2517400] chrMT 14674-14742 - | exon exon2135413 <NA>
#> [2517401] chrMT 15956-16023 - | exon exon2135414 <NA>
#> exon_id exon_name exon_rank Parent
#> <integer> <character> <integer> <character>
#> [1] <NA> <NA> <NA> <NA>
#> [2] <NA> <NA> <NA> <NA>
#> [3] <NA> <NA> <NA> <NA>
#> [4] <NA> <NA> <NA> <NA>
#> [5] <NA> <NA> <NA> <NA>
#> ... ... ... ... ...
#> [2517397] 891109 <NA> 1 mRNA381983
#> [2517398] 891110 <NA> 1 mRNA381984
#> [2517399] 891111 <NA> 1 mRNA381985
#> [2517400] 891112 <NA> 1 mRNA381986
#> [2517401] 891113 <NA> 1 mRNA381987
#> -------
#> seqinfo: 298 sequences (2 circular) from hg19 genome
You can convert a GFF object to a TxDb
object using the makeTxDbFromGFF in the txdbmaker
package. This is useful for creating a transcript database from GFF
annotations, which can then be used for various genomic analyses.
library(txdbmaker)
txdb <- makeTxDbFromGFF(
file = gff_local,
format = "gff3",
dataSource = "https://www.mirbase.org/download/hsa.gff3",
organism = "Homo sapiens",
taxonomyId = 9606
)#> Import genomic features from the file as a GRanges object ... OK
#> Prepare the 'metadata' data frame ... OK
#> Make the TxDb object ...
#> Warning in .extract_transcripts_from_GRanges(tx_IDX, gr, mcols0$type, mcols0$ID, : the transcript names ("tx_name" column in the TxDb object) imported from the
#> "Name" attribute are not unique
#> Warning in .makeTxDb_normarg_chrominfo(chrominfo): genome version information
#> is not available for this TxDb object
#> OK
genome <- grepv("genome-build-id", readLines(gff_local)) |>
strsplit("# genome-build-id:\\s+") |>
unlist() |>
tail(1L)
genome(txdb) <- genome
txdb#> TxDb object:
#> # Db type: TxDb
#> # Supporting package: GenomicFeatures
#> # Data source: https://www.mirbase.org/download/hsa.gff3
#> # Organism: Homo sapiens
#> # Taxonomy ID: 9606
#> # Genome: NA
#> # Nb of transcripts: 4801
#> # Db created by: txdbmaker package from Bioconductor
#> # Creation time: 2026-01-23 03:26:39 +0000 (Fri, 23 Jan 2026)
#> # txdbmaker version at creation time: 1.7.3
#> # RSQLite version at creation time: 2.4.5
#> # DBSCHEMAVERSION: 1.2
#> R version 4.5.2 (2025-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] txdbmaker_1.7.3
#> [2] TxDb.Hsapiens.UCSC.hg19.knownGene_3.22.1
#> [3] GenomicFeatures_1.63.1
#> [4] AnnotationDbi_1.73.0
#> [5] Biobase_2.71.0
#> [6] BiocFileCache_3.1.0
#> [7] dbplyr_2.5.1
#> [8] GenomicRanges_1.63.1
#> [9] Seqinfo_1.1.0
#> [10] IRanges_2.45.0
#> [11] S4Vectors_0.49.0
#> [12] BiocGenerics_0.57.0
#> [13] generics_0.1.4
#> [14] Bioc.gff_1.1.0
#> [15] BiocStyle_2.39.0
#>
#> loaded via a namespace (and not attached):
#> [1] tidyselect_1.2.1 dplyr_1.1.4
#> [3] blob_1.3.0 filelock_1.0.3
#> [5] Biostrings_2.79.4 bitops_1.0-9
#> [7] fastmap_1.2.0 RCurl_1.98-1.17
#> [9] GenomicAlignments_1.47.0 XML_3.99-0.20
#> [11] digest_0.6.39 lifecycle_1.0.5
#> [13] KEGGREST_1.51.1 RSQLite_2.4.5
#> [15] magrittr_2.0.4 compiler_4.5.2
#> [17] rlang_1.1.7 sass_0.4.10
#> [19] progress_1.2.3 tools_4.5.2
#> [21] yaml_2.3.12 rtracklayer_1.71.3
#> [23] knitr_1.51 prettyunits_1.2.0
#> [25] S4Arrays_1.11.1 bit_4.6.0
#> [27] curl_7.0.0 DelayedArray_0.37.0
#> [29] abind_1.4-8 BiocParallel_1.45.0
#> [31] withr_3.0.2 purrr_1.2.1
#> [33] sys_3.4.3 grid_4.5.2
#> [35] biomaRt_2.67.1 SummarizedExperiment_1.41.0
#> [37] cli_3.6.5 rmarkdown_2.30
#> [39] crayon_1.5.3 httr_1.4.7
#> [41] rjson_0.2.23 BiocBaseUtils_1.13.0
#> [43] DBI_1.2.3 cachem_1.1.0
#> [45] stringr_1.6.0 parallel_4.5.2
#> [47] BiocManager_1.30.27 XVector_0.51.0
#> [49] restfulr_0.0.16 matrixStats_1.5.0
#> [51] vctrs_0.7.0 Matrix_1.7-4
#> [53] jsonlite_2.0.0 hms_1.1.4
#> [55] bit64_4.6.0-1 maketools_1.3.2
#> [57] jquerylib_0.1.4 glue_1.8.0
#> [59] codetools_0.2-20 stringi_1.8.7
#> [61] GenomeInfoDb_1.47.2 BiocIO_1.21.0
#> [63] UCSC.utils_1.7.1 tibble_3.3.1
#> [65] pillar_1.11.1 rappdirs_0.3.4
#> [67] htmltools_0.5.9 GenomeInfoDbData_1.2.15
#> [69] R6_2.6.1 httr2_1.2.2
#> [71] evaluate_1.0.5 lattice_0.22-7
#> [73] png_0.1-8 Rsamtools_2.27.0
#> [75] cigarillo_1.1.0 memoise_2.0.1
#> [77] bslib_0.9.0 SparseArray_1.11.10
#> [79] xfun_0.56 MatrixGenerics_1.23.0
#> [81] buildtools_1.0.0 pkgconfig_2.0.3