| Title: | High-Performance GenBank I/O using 'gb-io' |
|---|---|
| Description: | Interfaces with the 'gb-io' Rust crate to provide fast reading and writing of GenBank files. Supports parsing of sequences, features, and metadata into R-friendly data structures. |
| Authors: | Richard Stöckl [aut, cre] (ORCID: <https://orcid.org/0000-0002-0451-0652>) |
| Maintainer: | Richard Stöckl <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.4.0 |
| Built: | 2026-05-27 15:13:23 UTC |
| Source: | https://github.com/richardstoeckl/rgbio |
Reads one or more GenBank records and returns selected components in either Bioconductor, tidy, or base format.
read_gbk( file, format = "bioconductor", sequences = TRUE, features = TRUE, metadata = TRUE, records = NULL, validate = TRUE )read_gbk( file, format = "bioconductor", sequences = TRUE, features = TRUE, metadata = TRUE, records = NULL, validate = TRUE )
file |
Character path to a GenBank file. |
format |
Output format. One of "bioconductor", "tidy", or "base". |
sequences |
Logical; include sequence data. |
features |
Logical; include feature annotations. |
metadata |
Logical; include record metadata. |
records |
Integer indices or character accession/name selectors. |
validate |
Logical; validate parsed records. |
A variable object based on selected components. Returns either a
single object or a named list with sequences, features, and/or metadata.
Deprecated compatibility wrapper for read_gbk().
read_genbank(file)read_genbank(file)
file |
Path to a GenBank file. |
Parsed records in legacy list format.
Writes GenBank records from sequence, feature, and metadata components.
write_gbk( file, sequences, features = NULL, metadata = NULL, append = FALSE, validate = TRUE, line_width = 80 )write_gbk( file, sequences, features = NULL, metadata = NULL, append = FALSE, validate = TRUE, line_width = 80 )
file |
Character path to output file. |
sequences |
DNAStringSet or named character vector. |
features |
GRanges with |
metadata |
DataFrame, data.frame, or list with record metadata. |
append |
Logical; append to file. |
validate |
Logical; validate inputs. |
line_width |
Integer sequence line width. |
Logical TRUE on success.
Deprecated compatibility wrapper for write_gbk().
write_genbank(file, sequence, features, metadata = list())write_genbank(file, sequence, features, metadata = list())
file |
Path to output file. |
sequence |
Sequence string. |
features |
Feature table with columns key, location, qualifiers. |
metadata |
Metadata list. |
Logical TRUE on success.