bedbaser 1.0.3
bedbaser is an R API client for BEDbase that provides access to the bedhost API and includes convenience functions, such as to create GRanges and GRangesList objects.
Install bedbaser using BiocManager.
if (!"BiocManager" %in% rownames(installed.packages())) {
install.packages("BiocManager")
}
BiocManager::install("bedbaser")
Load the package and create a BEDbase instance, optionally setting the cache
to cache_path
. If cache_path
is not set, bedbaser will
choose the default location.
library(bedbaser)
bedbase <- BEDbase(tempdir())
## 125977 BED files available.
bedbaser can use the same cache as
geniml’s BBClient by setting the
cache_path
to the same location. It will create the following structure:
cache_path
bedfiles
a/f/afile.bed.gz
bedsets
a/s/aset.txt
bedbaser includes convenience functions prefixed with bb_ to
facilitate finding BED files, exploring their metadata, downloading files, and
creating GRanges
objects.
Use bbs_stats()
to display the total available BED files, BEDsets, and
genomes. Set detailed
to TRUE
to display the type of BED formats and genomes
available.
Use bb_list_beds()
and bb_list_bedsets()
to browse available resources in
BEDbase. Both functions display the id and names of BED files and BEDsets. An
id can be used to access a specific resource.
bb_list_beds(bedbase)
## # A tibble: 1,000 × 52
## name genome_alias bed_compliance data_format compliant_columns
## <chr> <chr> <chr> <chr> <chr>
## 1 Plasma B cells - T… "hg38" bed6+4 encode_nar… 6
## 2 Young_Daughter 3y3A "mm10" bed5+5 bed_like_rs 5
## 3 DU01: Damage H2AK1… "" bed6+0 ucsc_bed 6
## 4 total RNA at diagn… "hg19" bed4+2 bed_like 4
## 5 encode_16468 "hg38" bed6+4 encode_nar… 6
## 6 Human colon cancer… "hg19" bed6+3 encode_bro… 6
## 7 FGFR1(RA)_ChIP_Seq "mm10" bed6+4 encode_nar… 6
## 8 encode_4362 "hg38" bed6+4 encode_nar… 6
## 9 TF ChIP-seq from H… "hg38" bed6+4 encode_nar… 6
## 10 H3K4me3 Wnt3KO Chi… "mm9" bed6+4 encode_nar… 6
## # ℹ 990 more rows
## # ℹ 47 more variables: non_compliant_columns <chr>, id <chr>,
## # description <chr>, submission_date <chr>, last_update_date <chr>,
## # is_universe <chr>, license_id <chr>, annotation.organism <chr>,
## # annotation.species_id <chr>, annotation.genotype <chr>,
## # annotation.phenotype <chr>, annotation.description <chr>,
## # annotation.cell_type <chr>, annotation.cell_line <chr>, …
bb_list_bedsets(bedbase)
## # A tibble: 1,000 × 9
## # Groups: id, name, md5sum, submission_date, last_update_date, description,
## # author, source [1,000]
## id name md5sum submission_date last_update_date description bed_ids
## <chr> <chr> <chr> <chr> <chr> <chr> <list>
## 1 encode_ba… enco… 462ea… 2025-04-25T01:… 2025-04-25T01:4… "Encode pr… <tibble>
## 2 encode_ba… enco… f6f15… 2025-04-24T23:… 2025-04-24T23:3… "Encode pr… <tibble>
## 3 encode_ba… enco… cf9f5… 2025-04-27T03:… 2025-04-27T03:1… "Encode pr… <tibble>
## 4 encode_ba… enco… 50cba… 2025-04-30T05:… 2025-04-30T05:3… "Encode pr… <tibble>
## 5 encode_ba… enco… 6054b… 2025-04-28T22:… 2025-04-28T22:2… "Encode pr… <tibble>
## 6 encode_ba… enco… 59f42… 2025-04-29T06:… 2025-04-29T06:2… "Encode pr… <tibble>
## 7 encode_ba… enco… 1cd93… 2025-05-02T20:… 2025-05-02T20:5… "Encode pr… <tibble>
## 8 excludera… excl… f6826… 2025-05-05T18:… 2025-05-05T18:3… "Exclude r… <tibble>
## 9 gse100000 gse1… 3da66… 2025-06-01T14:… 2025-06-01T14:5… "Data from… <tibble>
## 10 gse100302 gse1… a8df0… 2025-05-23T03:… 2025-05-23T03:4… "Data from… <tibble>
## # ℹ 990 more rows
## # ℹ 2 more variables: author <chr>, source <chr>
Use bb_metadata()
to learn more about a BED or BEDset associated with an id.
ex_bed <- bb_example(bedbase, "bed")
md <- bb_metadata(bedbase, ex_bed$id)
head(md)
## $name
## [1] "E1864"
##
## $genome_alias
## [1] "hg38"
##
## $genome_digest
## NULL
##
## $bed_compliance
## [1] "bed4+1"
##
## $data_format
## [1] "bed_like"
##
## $compliant_columns
## [1] 4
Use bb_beds_in_bedset()
to display the id of BEDs in a BEDset.
bb_beds_in_bedset(bedbase, "excluderanges")
## # A tibble: 81 × 32
## name genome_alias genome_digest bed_compliance data_format compliant_columns
## <chr> <chr> <chr> <chr> <chr> <chr>
## 1 mm10… mm10 0f10d83b1050… bed4+1 bed_like 4
## 2 hg38… hg38 2230c535660f… bed4+1 bed_like 4
## 3 T2T.… t2t-chm13 <NA> bed4+8 bed_like 4
## 4 TAIR… tair10 <NA> bed4+2 bed_like 4
## 5 mm9.… mm9 <NA> bed4+1 bed_like 4
## 6 T2T.… t2t-chm13 <NA> bed4+1 bed_like 4
## 7 danR… danrer10 <NA> bed4+1 bed_like 4
## 8 mm39… mm39 <NA> bed4+7 bed_like 4
## 9 mm9.… mm9 <NA> bed4+4 bed_like 4
## 10 hg19… hg19 baa91c8f6e27… bed4+7 bed_like 4
## # ℹ 71 more rows
## # ℹ 26 more variables: non_compliant_columns <chr>, id <chr>,
## # description <chr>, submission_date <chr>, last_update_date <chr>,
## # is_universe <chr>, license_id <chr>, annotation.species_name <chr>,
## # annotation.species_id <chr>, annotation.genotype <chr>,
## # annotation.phenotype <chr>, annotation.description <chr>,
## # annotation.cell_type <chr>, annotation.cell_line <chr>, …
Search for BED files by keywords. bb_bed_text_search()
returns all BED files
scored against a keyword query.
bb_bed_text_search(bedbase, "cancer", limit = 10)
## # A tibble: 10 × 52
## id payload.species_name payload.species_id payload.genotype
## <chr> <chr> <chr> <chr>
## 1 b33dca39-9e67-576a-… Homo sapiens 9606 ""
## 2 43b3e61f-5e0b-b8fe-… Homo sapiens 9606 ""
## 3 bb612399-36ab-6125-… Homo sapiens 9606 ""
## 4 0f9848ff-51dc-5e7f-… Homo sapiens 9606 ""
## 5 910a5753-a15b-5afc-… Homo sapiens 9606 ""
## 6 b85718c4-511d-85e4-… Homo sapiens 9606 ""
## 7 5fd3b113-2e75-23ea-… Homo sapiens 9606 ""
## 8 5bd3a925-7dc2-1386-… Homo sapiens 9606 ""
## 9 835c0159-780f-ca54-… Homo sapiens 9606 ""
## 10 9c273213-8d8a-f64f-… Homo sapiens 9606 ""
## # ℹ 48 more variables: payload.phenotype <chr>, payload.description <chr>,
## # payload.cell_type <chr>, payload.cell_line <chr>, payload.tissue <chr>,
## # payload.library_source <chr>, payload.assay <chr>, payload.antibody <chr>,
## # payload.target <chr>, payload.treatment <chr>,
## # payload.global_sample_id <chr>, payload.global_experiment_id <chr>,
## # payload.original_file_name <chr>, score <chr>, metadata.name <chr>,
## # metadata.genome_alias <chr>, metadata.genome_digest <chr>, …
Create a GRanges object with a BED id with bb_to_granges
, which
downloads and imports a BED file using rtracklayer.
ex_bed <- bb_example(bedbase, "bed")
# Allow bedbaser to assign column names and types
bb_to_granges(bedbase, ex_bed$id, quietly = FALSE)
## Assigning column names and types.
##
## Attaching package: 'Biostrings'
##
##
## The following object is masked from 'package:base':
##
## strsplit
## GRanges object with 4493 ranges and 2 metadata columns:
## seqnames ranges strand | name V5
## <Rle> <IRanges> <Rle> | <character> <numeric>
## [1] chr1 750001-800000 * | <NA> 2.48932
## [2] chr1 2150001-2200000 * | <NA> 1.42802
## [3] chr1 2650001-2700000 * | <NA> 1.21170
## [4] chr1 2800001-2850000 * | <NA> 1.86796
## [5] chr1 3450001-3500000 * | <NA> 2.82448
## ... ... ... ... . ... ...
## [4489] chrX 152650001-152700000 * | <NA> 1.67777
## [4490] chrX 152950001-153000000 * | <NA> 3.18726
## [4491] chrX 153500001-153550000 * | <NA> 1.44397
## [4492] chrX 154550001-154600000 * | <NA> 2.68313
## [4493] chrX 154900001-154950000 * | <NA> 1.96125
## -------
## seqinfo: 711 sequences (1 circular) from hg38 genome
For BEDX+Y formats, a named list with column types may be passed through
extra_cols
if the column name and type are known. Otherwise, bb_to_granges
guesses the column types and assigns column names.
# Manually assign column name and type using `extra_cols`
bb_to_granges(bedbase, ex_bed$id, extra_cols = c("column_name" = "character"))
bb_to_granges
automatically assigns the column names and types for broad peak
and narrow peak files.
bed_id <- "bbad85f21962bb8d972444f7f9a3a932"
md <- bb_metadata(bedbase, bed_id)
head(md)
## $name
## [1] "PM_137_NPC_CTCF_ChIP"
##
## $genome_alias
## [1] "hg38"
##
## $genome_digest
## [1] "2230c535660fb4774114bfa966a62f823fdb6d21acf138d4"
##
## $bed_compliance
## [1] "bed6+4"
##
## $data_format
## [1] "encode_narrowpeak_rs"
##
## $compliant_columns
## [1] 6
bb_to_granges(bedbase, bed_id)
## GRanges object with 26210 ranges and 6 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <numeric>
## [1] chr1 869762-870077 * | 111-11-DSP-NPC-CTCF-.. 587
## [2] chr1 904638-904908 * | 111-11-DSP-NPC-CTCF-.. 848
## [3] chr1 921139-921331 * | 111-11-DSP-NPC-CTCF-.. 177
## [4] chr1 939191-939364 * | 111-11-DSP-NPC-CTCF-.. 139
## [5] chr1 976105-976282 * | 111-11-DSP-NPC-CTCF-.. 185
## ... ... ... ... . ... ...
## [26206] chrY 18445992-18446211 * | 111-11-DSP-NPC-CTCF-.. 203
## [26207] chrY 18608331-18608547 * | 111-11-DSP-NPC-CTCF-.. 203
## [26208] chrY 18669820-18670062 * | 111-11-DSP-NPC-CTCF-.. 244
## [26209] chrY 18997783-18997956 * | 111-11-DSP-NPC-CTCF-.. 191
## [26210] chrY 19433165-19433380 * | 111-11-DSP-NPC-CTCF-.. 275
## signalValue pValue qValue peak
## <numeric> <numeric> <numeric> <integer>
## [1] 20.94161 58.7971 54.9321 152
## [2] 30.90682 84.8282 80.3102 118
## [3] 9.62671 17.7065 14.8446 69
## [4] 8.10671 13.9033 11.1352 49
## [5] 9.26375 18.5796 15.6985 129
## ... ... ... ... ...
## [26206] 10.64005 20.3549 17.4328 106
## [26207] 8.00064 20.3991 17.4753 149
## [26208] 12.16006 24.4764 21.4585 119
## [26209] 8.97342 19.1163 16.2230 69
## [26210] 12.21130 27.5139 24.4211 89
## -------
## seqinfo: 711 sequences (1 circular) from hg38 genome
Create a GRangesList given a BEDset id with bb_to_grangeslist
.
bedset_id <- "lola_hg38_ucsc_features"
bb_to_grangeslist(bedbase, bedset_id)
## GRangesList object of length 11:
## [[1]]
## GRanges object with 864 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 690078-6272609 *
## [2] chr1 690078-2326424 *
## [3] chr1 771707-6806566 *
## [4] chr1 771707-3153758 *
## [5] chr1 805477-4942653 *
## ... ... ... ...
## [860] chrY 23762211-26011096 *
## [861] chrY 23762211-26011096 *
## [862] chrY 23774007-25910251 *
## [863] chrY 26011096-26174983 *
## [864] chrY 26312489-26653776 *
## -------
## seqinfo: 711 sequences (1 circular) from hg38 genome
##
## ...
## <10 more elements>
Save BED files or BEDsets with bb_save
:
bb_save(bedbase, ex_bed$id, tempdir())
Because bedbaser uses the AnVIL Service class, it’s possible to access any endpoint of the BEDbase API.
show(bedbase)
## service: bedbase
## host: api.bedbase.org
## tags(); use bedbase$<tab completion>:
## # A tibble: 44 × 3
## tag operation summary
## <chr> <chr> <chr>
## 1 base get_bedbase_db_stats_v1_genomes_get Get av…
## 2 base get_bedbase_db_stats_v1_stats_get Get su…
## 3 base get_detailed_stats_v1_detailed_stats_get Get de…
## 4 base get_detailed_usage_v1_detailed_usage_get Get de…
## 5 base redirect_to_download_v1_files__file_path__get Redire…
## 6 base service_info_v1_service_info_get GA4GH …
## 7 bed bed_to_bed_search_v1_bed_search_bed_post Search…
## 8 bed embed_bed_file_v1_bed_embed_post Get em…
## 9 bed get_bed_classification_v1_bed__bed_id__metadata_classification… Get cl…
## 10 bed get_bed_embedding_v1_bed__bed_id__embedding_get Get em…
## # ℹ 34 more rows
## tag values:
## base, bed, bedset, home, objects, search, NA
## schemas():
## AccessMethod, AccessURL, BaseListResponse, BedClassification,
## BedEmbeddingResult
## # ... with 42 more elements
For example, to access a BED file’s stats, access the endpoint with $
and use
httr to get the result. show
will display information about the
endpoint.
library(httr)
##
## Attaching package: 'httr'
## The following object is masked from 'package:Biobase':
##
## content
show(bedbase$get_bed_stats_v1_bed__bed_id__metadata_stats_get)
## get_bed_stats_v1_bed__bed_id__metadata_stats_get
## Get stats for a single BED record
## Description:
## Example bed_id: bbad85f21962bb8d972444f7f9a3a932
##
## Parameters:
## bed_id (string)
## BED digest
id <- "bbad85f21962bb8d972444f7f9a3a932"
rsp <- bedbase$get_bed_stats_v1_bed__bed_id__metadata_stats_get(id)
content(rsp)
## $number_of_regions
## [1] 26210
##
## $gc_content
## [1] 0.5
##
## $median_tss_dist
## [1] 31480
##
## $mean_region_width
## [1] 276.3
##
## $exon_frequency
## [1] 1358
##
## $exon_percentage
## [1] 0.0518
##
## $intron_frequency
## [1] 9390
##
## $intron_percentage
## [1] 0.3583
##
## $intergenic_percentage
## [1] 0.4441
##
## $intergenic_frequency
## [1] 11639
##
## $promotercore_frequency
## [1] 985
##
## $promotercore_percentage
## [1] 0.0376
##
## $fiveutr_frequency
## [1] 720
##
## $fiveutr_percentage
## [1] 0.0275
##
## $threeutr_frequency
## [1] 1074
##
## $threeutr_percentage
## [1] 0.041
##
## $promoterprox_frequency
## [1] 1044
##
## $promoterprox_percentage
## [1] 0.0398
Given a BED id, we can use liftOver to convert one genomic coordinate system to another.
Install liftOver and rtracklayer then load the packages.
if (!"BiocManager" %in% rownames(installed.packages())) {
install.packages("BiocManager")
}
BiocManager::install(c("liftOver", "rtracklayer"))
library(liftOver)
library(rtracklayer)
Create a GRanges object from a
mouse genome.
Create a BEDbase Service instance. Use the instance to create a GRanges
object from the BEDbase id
.
id <- "f2a5b06011706376560514c3f39648ea"
bedbase <- BEDbase()
gro <- bb_to_granges(bedbase, id)
gro
## GRanges object with 132610 ranges and 2 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <numeric>
## [1] chr1 3132268-3132768 + | chr1-21633 1
## [2] chr1 3185464-3185964 + | chr1-21634 1
## [3] chr1 3221560-3222060 + | chr1-6085 1
## [4] chr1 3476307-3476807 + | chr1-21635 1
## [5] chr1 3560226-3561000 + | chr1-4747 1
## ... ... ... ... . ... ...
## [132606] chrY 90737580-90739215 + | chrY-23 1
## [132607] chrY 90742758-90744732 + | chrY-35 1
## [132608] chrY 90810972-90814119 + | chrY-47 1
## [132609] chrY 90819248-90819748 + | chrY-131 1
## [132610] chrY 90828312-90828949 + | chrY-103 1
## -------
## seqinfo: 239 sequences (1 circular) from mm10 genome
Download the chain file from UCSC.
chain_url <- paste0(
"https://hgdownload.cse.ucsc.edu/goldenPath/mm10/liftOver/",
"mm10ToMm39.over.chain.gz"
)
tmpdir <- tempdir()
gz <- file.path(tmpdir, "mm10ToMm39.over.chain.gz")
download.file(chain_url, gz)
gunzip(gz, remove = FALSE)
Import the chain, set the sequence levels style, and set the genome for the GRanges object.
ch <- import.chain(file.path(tmpdir, "mm10ToMm39.over.chain"))
seqlevelsStyle(gro) <- "UCSC"
gro39 <- liftOver(gro, ch)
gro39 <- unlist(gro39)
genome(gro39) <- "mm39"
gro39
## GRanges object with 132675 ranges and 2 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <numeric>
## [1] chr1 3202491-3202991 + | chr1-21633 1
## [2] chr1 3255687-3256187 + | chr1-21634 1
## [3] chr1 3291783-3292283 + | chr1-6085 1
## [4] chr1 3546530-3547030 + | chr1-21635 1
## [5] chr1 3630449-3631223 + | chr1-4747 1
## ... ... ... ... . ... ...
## [132671] chrY 90748849-90750484 + | chrY-23 1
## [132672] chrY 90754027-90756001 + | chrY-35 1
## [132673] chrY 90822241-90825388 + | chrY-47 1
## [132674] chrY 90830517-90831017 + | chrY-131 1
## [132675] chrY 90839581-90840218 + | chrY-103 1
## -------
## seqinfo: 21 sequences from mm39 genome; no seqlengths
sessionInfo()
## R version 4.5.1 (2025-06-13)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB 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: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] BSgenome.Mmusculus.UCSC.mm10_1.4.3
## [2] httr_1.4.7
## [3] BSgenome.Hsapiens.UCSC.hg38_1.4.5
## [4] BSgenome_1.76.0
## [5] BiocIO_1.18.0
## [6] Biostrings_2.76.0
## [7] XVector_0.48.0
## [8] bedbaser_1.0.3
## [9] liftOver_1.32.0
## [10] Homo.sapiens_1.3.1
## [11] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [12] org.Hs.eg.db_3.21.0
## [13] GO.db_3.21.0
## [14] OrganismDbi_1.50.0
## [15] GenomicFeatures_1.60.0
## [16] AnnotationDbi_1.70.0
## [17] Biobase_2.68.0
## [18] gwascat_2.40.0
## [19] R.utils_2.13.0
## [20] R.oo_1.27.1
## [21] R.methodsS3_1.8.2
## [22] rtracklayer_1.68.0
## [23] GenomicRanges_1.60.0
## [24] GenomeInfoDb_1.44.0
## [25] IRanges_2.42.0
## [26] S4Vectors_0.46.0
## [27] BiocGenerics_0.54.0
## [28] generics_0.1.4
## [29] BiocStyle_2.36.0
##
## loaded via a namespace (and not attached):
## [1] jsonlite_2.0.0 magrittr_2.0.3
## [3] rmarkdown_2.29 zlibbioc_1.54.0
## [5] vctrs_0.6.5 memoise_2.0.1
## [7] Rsamtools_2.24.0 RCurl_1.98-1.17
## [9] htmltools_0.5.8.1 S4Arrays_1.8.1
## [11] BiocBaseUtils_1.10.0 progress_1.2.3
## [13] lambda.r_1.2.4 curl_6.4.0
## [15] SparseArray_1.8.1 sass_0.4.10
## [17] bslib_0.9.0 htmlwidgets_1.6.4
## [19] httr2_1.2.0 futile.options_1.0.1
## [21] cachem_1.1.0 GenomicAlignments_1.44.0
## [23] mime_0.13 lifecycle_1.0.4
## [25] pkgconfig_2.0.3 Matrix_1.7-3
## [27] R6_2.6.1 fastmap_1.2.0
## [29] GenomeInfoDbData_1.2.14 MatrixGenerics_1.20.0
## [31] shiny_1.11.1 digest_0.6.37
## [33] RSQLite_2.4.2 filelock_1.0.3
## [35] abind_1.4-8 compiler_4.5.1
## [37] withr_3.0.2 bit64_4.6.0-1
## [39] BiocParallel_1.42.1 DBI_1.2.3
## [41] biomaRt_2.64.0 rappdirs_0.3.3
## [43] DelayedArray_0.34.1 rjson_0.2.23
## [45] tools_4.5.1 httpuv_1.6.16
## [47] glue_1.8.0 restfulr_0.0.16
## [49] promises_1.3.3 grid_4.5.1
## [51] tidyr_1.3.1 data.table_1.17.8
## [53] hms_1.1.3 utf8_1.2.6
## [55] xml2_1.3.8 pillar_1.11.0
## [57] stringr_1.5.1 later_1.4.2
## [59] splines_4.5.1 dplyr_1.1.4
## [61] BiocFileCache_2.16.1 lattice_0.22-7
## [63] survival_3.8-3 bit_4.6.0
## [65] tidyselect_1.2.1 RBGL_1.84.0
## [67] miniUI_0.1.2 knitr_1.50
## [69] bookdown_0.43 SummarizedExperiment_1.38.1
## [71] snpStats_1.58.0 futile.logger_1.4.3
## [73] xfun_0.52 matrixStats_1.5.0
## [75] DT_0.33 stringi_1.8.7
## [77] UCSC.utils_1.4.0 yaml_2.3.10
## [79] evaluate_1.0.4 codetools_0.2-20
## [81] tibble_3.3.0 AnVILBase_1.2.0
## [83] BiocManager_1.30.26 graph_1.86.0
## [85] cli_3.6.5 AnVIL_1.20.1
## [87] xtable_1.8-4 jquerylib_0.1.4
## [89] Rcpp_1.1.0 dbplyr_2.5.0
## [91] png_0.1-8 XML_3.99-0.18
## [93] rapiclient_0.1.8 parallel_4.5.1
## [95] blob_1.2.4 prettyunits_1.2.0
## [97] bitops_1.0-9 txdbmaker_1.4.2
## [99] VariantAnnotation_1.54.1 purrr_1.1.0
## [101] crayon_1.5.3 rlang_1.1.6
## [103] KEGGREST_1.48.1 formatR_1.14