Skip to contents

Introduction

When working on your own genome project or when using publicly available genomes for comparative analyses, it is critical to assess the quality of your data. Over the past years, several tools have been developed and several metrics have been proposed to assess the quality of a genome assembly and annotation. cogeqc helps users interpret their genome assembly statistics by comparing them with statistics on publicly available genomes on the NCBI. Additionally, cogeqc also provides an interface to BUSCO (Simão et al. 2015), a popular tool to assess gene space completeness. Graphical functions are available to make publication-ready plots that summarize the results of quality control.

Installation

You can install cogeqc from Bioconductor with the following code:

if(!requireNamespace('BiocManager', quietly = TRUE))
  install.packages('BiocManager')
BiocManager::install("cogeqc")
# Load package after installation
library(cogeqc)

Assessing genome assembly quality: statistics in a context

When analyzing and interpreting genome assembly statistics, it is often useful to place your stats in a context by comparing them with stats from genomes of closely-related or even the same species. cogeqc provides users with an interface to the NCBI Datasets API, which can be used to retrieve summary stats for genomes on NCBI. In this section, we will guide you on how to retrieve such information and use it as a reference to interpret your data.

Obtaining assembly statistics for NCBI genomes

To obtain a data frame of summary statistics for NCBI genomes of a particular taxon, you will use the function get_genome_stats(). In the taxon parameter, you must specify the taxon from which data will be extracted. This can be done either by passing a character scalar with taxon name or by passing a numeric scalar with NCBI Taxonomy ID. For example, the code below demonstrates two ways of extracting stats on maize (Zea mays) genomes on NCBI:

# Example 1: get stats for all maize genomes using taxon name
maize_stats <- get_genome_stats(taxon = "Zea mays")
head(maize_stats)
#>         accession  source species_taxid species_name species_common_name
#> 1 GCA_902167145.1 GENBANK          4577     Zea mays               maize
#> 2 GCF_902167145.1  REFSEQ          4577     Zea mays               maize
#> 3 GCA_022117705.1 GENBANK          4577     Zea mays               maize
#> 4 GCA_029775835.1 GENBANK          4577     Zea mays               maize
#> 5 GCA_905067065.1 GENBANK          4577     Zea mays               maize
#> 6 GCA_902714155.1 GENBANK          4577     Zea mays               maize
#>   species_ecotype species_strain species_isolate species_cultivar
#> 1            <NA>             NA            <NA>              B73
#> 2            <NA>             NA            <NA>              B73
#> 3            <NA>             NA            <NA>        Mo17-2021
#> 4            <NA>             NA            <NA>           LT2357
#> 5            <NA>             NA            <NA>             <NA>
#> 6            <NA>             NA            <NA>         B73 Ab10
#>   assembly_level assembly_status                      assembly_name
#> 1     Chromosome         current           Zm-B73-REFERENCE-NAM-5.0
#> 2     Chromosome         current           Zm-B73-REFERENCE-NAM-5.0
#> 3           <NA>         current Zm-Mo17-REFERENCE-CAU-T2T-assembly
#> 4     Chromosome         current                       ASM2977583v1
#> 5     Chromosome         current       Zm-LH244-REFERENCE-BAYER-1.0
#> 6     Chromosome         current     Zm-B73_AB10-REFERENCE-NAM-1.0b
#>   assembly_type submission_date                      submitter
#> 1       haploid              NA                       MaizeGDB
#> 2       haploid              NA                       MaizeGDB
#> 3       haploid              NA   China Agriculture University
#> 4       haploid              NA Beijing Lantron Seed Co., LTD.
#> 5       haploid              NA              BAYER CROPSCIENCE
#> 6       haploid              NA                       MaizeGDB
#>                sequencing_technology atypical       refseq_category
#> 1                               <NA>    FALSE                  <NA>
#> 2                               <NA>    FALSE representative genome
#> 3 Oxford Nanopore PromethION; PacBio    FALSE                  <NA>
#> 4                      PacBio Sequel    FALSE                  <NA>
#> 5                               <NA>    FALSE                  <NA>
#> 6                               <NA>    FALSE                  <NA>
#>   chromosome_count sequence_length ungapped_length contig_count contig_N50
#> 1               10      2182075994      2178268108         1393   47037903
#> 2               10      2182075994      2178268108         1393   47037903
#> 3               10      2178604320      2178604320           10  220303002
#> 4               10      2106865080      2106637080          460   15883073
#> 5               10      2147745480      2107651308        56173      84946
#> 6               10      2243621556      2241350720         1016  161994764
#>   contig_L50 scaffold_count scaffold_N50 scaffold_L50 GC_percent
#> 1         16            685    226353449            5       47.0
#> 2         16            685    226353449            5       47.0
#> 3          5             10    220303002            5       47.0
#> 4         41             10    222005600            5       47.0
#> 5       7498             10    225452224            5       47.0
#> 6          6            936    225306452            5       46.5
#>   annotation_provider annotation_release_date gene_count_total
#> 1                <NA>                    <NA>               NA
#> 2         NCBI RefSeq              2020-08-09            49897
#> 3                <NA>                    <NA>               NA
#> 4                <NA>                    <NA>               NA
#> 5                <NA>                    <NA>               NA
#> 6                <NA>                    <NA>               NA
#>   gene_count_coding gene_count_noncoding gene_count_pseudogene gene_count_other
#> 1                NA                   NA                    NA               NA
#> 2             34337                10366                  5194               NA
#> 3                NA                   NA                    NA               NA
#> 4                NA                   NA                    NA               NA
#> 5                NA                   NA                    NA               NA
#> 6                NA                   NA                    NA               NA
#>   CC_ratio
#> 1    139.3
#> 2    139.3
#> 3      1.0
#> 4     46.0
#> 5   5617.3
#> 6    101.6
str(maize_stats)
#> 'data.frame':    110 obs. of  36 variables:
#>  $ accession              : chr  "GCA_902167145.1" "GCF_902167145.1" "GCA_022117705.1" "GCA_029775835.1" ...
#>  $ source                 : chr  "GENBANK" "REFSEQ" "GENBANK" "GENBANK" ...
#>  $ species_taxid          : int  4577 4577 4577 4577 4577 4577 4577 4577 4577 4577 ...
#>  $ species_name           : chr  "Zea mays" "Zea mays" "Zea mays" "Zea mays" ...
#>  $ species_common_name    : chr  "maize" "maize" "maize" "maize" ...
#>  $ species_ecotype        : chr  NA NA NA NA ...
#>  $ species_strain         : logi  NA NA NA NA NA NA ...
#>  $ species_isolate        : chr  NA NA NA NA ...
#>  $ species_cultivar       : chr  "B73" "B73" "Mo17-2021" "LT2357" ...
#>  $ assembly_level         : Factor w/ 4 levels "Complete","Chromosome",..: 2 2 NA 2 2 2 2 2 2 2 ...
#>  $ assembly_status        : chr  "current" "current" "current" "current" ...
#>  $ assembly_name          : chr  "Zm-B73-REFERENCE-NAM-5.0" "Zm-B73-REFERENCE-NAM-5.0" "Zm-Mo17-REFERENCE-CAU-T2T-assembly" "ASM2977583v1" ...
#>  $ assembly_type          : chr  "haploid" "haploid" "haploid" "haploid" ...
#>  $ submission_date        : logi  NA NA NA NA NA NA ...
#>  $ submitter              : chr  "MaizeGDB" "MaizeGDB" "China Agriculture University" "Beijing Lantron Seed Co., LTD." ...
#>  $ sequencing_technology  : chr  NA NA "Oxford Nanopore PromethION; PacBio" "PacBio Sequel" ...
#>  $ atypical               : logi  FALSE FALSE FALSE FALSE FALSE FALSE ...
#>  $ refseq_category        : chr  NA "representative genome" NA NA ...
#>  $ chromosome_count       : int  10 10 10 10 10 10 10 10 10 10 ...
#>  $ sequence_length        : num  2.18e+09 2.18e+09 2.18e+09 2.11e+09 2.15e+09 ...
#>  $ ungapped_length        : num  2.18e+09 2.18e+09 2.18e+09 2.11e+09 2.11e+09 ...
#>  $ contig_count           : int  1393 1393 10 460 56173 1016 1191 1747 2634 863 ...
#>  $ contig_N50             : int  47037903 47037903 220303002 15883073 84946 161994764 49888411 49071148 38847693 35863069 ...
#>  $ contig_L50             : int  16 16 5 41 7498 6 12 13 21 21 ...
#>  $ scaffold_count         : int  685 685 10 10 10 936 800 1379 2223 413 ...
#>  $ scaffold_N50           : int  226353449 226353449 220303002 222005600 225452224 225306452 220287990 223168870 222765871 223950520 ...
#>  $ scaffold_L50           : int  5 5 5 5 5 5 5 5 5 5 ...
#>  $ GC_percent             : num  47 47 47 47 47 46.5 47 46.5 47 47 ...
#>  $ annotation_provider    : chr  NA "NCBI RefSeq" NA NA ...
#>  $ annotation_release_date: chr  NA "2020-08-09" NA NA ...
#>  $ gene_count_total       : int  NA 49897 NA NA NA NA NA NA NA NA ...
#>  $ gene_count_coding      : int  NA 34337 NA NA NA NA NA NA NA NA ...
#>  $ gene_count_noncoding   : int  NA 10366 NA NA NA NA NA NA NA NA ...
#>  $ gene_count_pseudogene  : int  NA 5194 NA NA NA NA NA NA NA NA ...
#>  $ gene_count_other       : int  NA NA NA NA NA NA NA NA NA NA ...
#>  $ CC_ratio               : num  139 139 1 46 5617 ...

# Example 2: get stats for all maize genomes using NCBI Taxonomy ID
maize_stats2 <- get_genome_stats(taxon = 4577)

# Checking if objects are the same
identical(maize_stats, maize_stats2)
#> [1] TRUE

As you can see, there are 110 maize genomes on the NCBI. You can also include filters in your searches by passing a list of key-value pairs with keys in list names and values in elements. For instance, to obtain only chromosome-scale and annotated maize genomes, you would run:

# Get chromosome-scale maize genomes with annotation
## Create list of filters
filt <- list(
    filters.has_annotation = "true",
    filters.assembly_level = "chromosome"
)
filt
#> $filters.has_annotation
#> [1] "true"
#> 
#> $filters.assembly_level
#> [1] "chromosome"

## Obtain data
filtered_maize_genomes <- get_genome_stats(taxon = "Zea mays", filters = filt)
dim(filtered_maize_genomes)
#> [1]  4 36

For a full list of filtering parameters and possible arguments, see the API documentation.

Comparing custom stats with NCBI stats

Now, suppose you sequenced a genome, obtained assembly and annotation stats, and want to compare them to NCBI genomes to identify potential issues. Examples of situations you may encounter include:

  • The genome you assembled is huge and you think there might be a problem with your assembly.

  • Your gene annotation pipeline predicted n genes, but you are not sure if this number is reasonable compared to other assemblies of the same species or closely-related species.

To compare user-defined summary stats with NCBI stats, you will use the function compare_genome_stats(). This function will include the values you observed for each statistic into a distribution (based on NCBI stats) and return the percentile and rank of your observed values in each distribution.

As an example, let’s go back to our maize stats we obtained in the previous section. Suppose you sequenced a new maize genome and observed the following values:

  1. Genome size = 2.4 Gb
  2. Number of genes = 50,000
  3. CC ratio = 2 1

To compare your observed values with those for publicly available maize genomes, you need to store them in a data frame. The column accession is mandatory, and any other column will be matched against columns in the data frame obtained with get_genome_stats(). Thus, make sure column names in your data frame match column names in the reference data frame. Then, you can compare both data frames as below:

# Check column names in the data frame of stats for maize genomes on the NCBI
names(maize_stats)
#>  [1] "accession"               "source"                 
#>  [3] "species_taxid"           "species_name"           
#>  [5] "species_common_name"     "species_ecotype"        
#>  [7] "species_strain"          "species_isolate"        
#>  [9] "species_cultivar"        "assembly_level"         
#> [11] "assembly_status"         "assembly_name"          
#> [13] "assembly_type"           "submission_date"        
#> [15] "submitter"               "sequencing_technology"  
#> [17] "atypical"                "refseq_category"        
#> [19] "chromosome_count"        "sequence_length"        
#> [21] "ungapped_length"         "contig_count"           
#> [23] "contig_N50"              "contig_L50"             
#> [25] "scaffold_count"          "scaffold_N50"           
#> [27] "scaffold_L50"            "GC_percent"             
#> [29] "annotation_provider"     "annotation_release_date"
#> [31] "gene_count_total"        "gene_count_coding"      
#> [33] "gene_count_noncoding"    "gene_count_pseudogene"  
#> [35] "gene_count_other"        "CC_ratio"

# Create a simulated data frame of stats for a maize genome
my_stats <- data.frame(
    accession = "my_lovely_maize",
    sequence_length = 2.4 * 1e9,
    gene_count_total = 50000,
    CC_ratio = 2
)

# Compare stats
compare_genome_stats(ncbi_stats = maize_stats, user_stats = my_stats)
#>         accession         variable percentile rank
#> 1 my_lovely_maize  sequence_length 0.98198198    3
#> 2 my_lovely_maize gene_count_total 1.00000000    1
#> 3 my_lovely_maize         CC_ratio 0.02857143    2

Visualizing summary assembly statistics

To have a visual representation of the summary stats obtained with get_genome_stats(), you will use the function plot_genome_stats().

# Summarize genome stats in a plot
plot_genome_stats(ncbi_stats = maize_stats)

Finally, you can pass your data frame of observed stats to highlight your values (as red points) in the distributions.

plot_genome_stats(ncbi_stats = maize_stats, user_stats = my_stats)

Assessing gene space completeness with BUSCO

One of the most common metrics to assess gene space completeness is BUSCO (best universal single-copy orthologs) (Simão et al. 2015). cogeqc allows users to run BUSCO from an R session and visualize results graphically. BUSCO summary statistics will help you assess which assemblies have high quality based on the percentage of complete BUSCOs.

Running BUSCO

To run BUSCO from R, you will use the function run_busco()2. Here, we will use an example FASTA file containing the first 1,000 lines of the Herbaspirilllum seropedicae SmR1 genome (GCA_000143225), which was downloaded from Ensembl Bacteria. We will run BUSCO using burkholderiales_odb10 as the lineage dataset. To view all available datasets, run list_busco_datasets().

# Path to FASTA file
sequence <- system.file("extdata", "Hse_subset.fa", package = "cogeqc")

# Path to directory where BUSCO datasets will be stored
download_path <- paste0(tempdir(), "/datasets")

# Run BUSCO if it is installed
if(busco_is_installed()) {
  run_busco(sequence, outlabel = "Hse", mode = "genome",
            lineage = "burkholderiales_odb10",
            outpath = tempdir(), download_path = download_path)
}

The output will be stored in the directory specified in outpath. You can read and parse BUSCO’s output with the function read_busco(). For example, let’s read the output of a BUSCO run using the genome of the green algae Ostreococcus tauri. The output directory is /extdata.

# Path to output directory
output_dir <- system.file("extdata", package = "cogeqc")

busco_summary <- read_busco(output_dir)
busco_summary
#>                Class Frequency           Lineage
#> 1        Complete_SC      1412 chlorophyta_odb10
#> 2 Complete_duplicate         4 chlorophyta_odb10
#> 3         Fragmented        35 chlorophyta_odb10
#> 4            Missing        68 chlorophyta_odb10

This is an example output for a BUSCO run with a single FASTA file. You can also specify a directory containing multiple FASTA files in the sequence argument of run_busco(). This way, BUSCO will be run in batch mode. Let’s see what the output of BUSCO in batch mode looks like:

data(batch_summary)
batch_summary
#>                Class Frequency               Lineage   File
#> 1        Complete_SC      98.5 burkholderiales_odb10 Hse.fa
#> 2        Complete_SC      98.8 burkholderiales_odb10 Hru.fa
#> 3 Complete_duplicate       0.7 burkholderiales_odb10 Hse.fa
#> 4 Complete_duplicate       0.7 burkholderiales_odb10 Hru.fa
#> 5         Fragmented       0.4 burkholderiales_odb10 Hse.fa
#> 6         Fragmented       0.3 burkholderiales_odb10 Hru.fa
#> 7            Missing       0.4 burkholderiales_odb10 Hse.fa
#> 8            Missing       0.2 burkholderiales_odb10 Hru.fa

The only difference between this data frame and the previous one is the column File, which contains information on the FASTA file. The example dataset batch_summary contains the output of run_busco() using a directory containing two genomes (Herbaspirillum seropedicae SmR1 and Herbaspirillum rubrisubalbicans M1) as parameter to the sequence argument.

Visualizing BUSCO summary statistics

After using run_busco() and parsing its output with read_busco(), users can visualize summary statistics with plot_busco().

# Single FASTA file - Ostreococcus tauri
plot_busco(busco_summary)


# Batch mode - Herbaspirillum seropedicae and H. rubrisubalbicans
plot_busco(batch_summary)

We usually consider genomes with >90% of complete BUSCOs as having high quality. Thus, we can conclude that the three genomes analyzed here are high-quality genomes.

Session information

This document was created under the following conditions:

sessioninfo::session_info()
#> ─ Session info ───────────────────────────────────────────────────────────────
#>  setting  value
#>  version  R version 4.3.2 (2023-10-31)
#>  os       Ubuntu 22.04.3 LTS
#>  system   x86_64, linux-gnu
#>  ui       X11
#>  language en
#>  collate  en_US.UTF-8
#>  ctype    en_US.UTF-8
#>  tz       UTC
#>  date     2024-01-12
#>  pandoc   3.1.1 @ /usr/local/bin/ (via rmarkdown)
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────
#>  package          * version   date (UTC) lib source
#>  ape                5.7-1     2023-03-13 [1] RSPM (R 4.3.0)
#>  aplot              0.2.2     2023-10-06 [1] RSPM (R 4.3.0)
#>  beeswarm           0.4.0     2021-06-01 [1] RSPM (R 4.3.0)
#>  BiocGenerics       0.48.1    2023-11-01 [1] Bioconductor
#>  BiocManager        1.30.22   2023-08-08 [2] CRAN (R 4.3.2)
#>  BiocStyle        * 2.30.0    2023-10-24 [1] Bioconductor
#>  Biostrings         2.70.1    2023-10-25 [1] Bioconductor
#>  bitops             1.0-7     2021-04-24 [1] RSPM (R 4.3.0)
#>  bookdown           0.37      2023-12-01 [1] RSPM (R 4.3.0)
#>  bslib              0.6.1     2023-11-28 [2] RSPM (R 4.3.0)
#>  cachem             1.0.8     2023-05-01 [2] RSPM (R 4.3.0)
#>  cli                3.6.2     2023-12-11 [2] RSPM (R 4.3.0)
#>  cogeqc           * 1.7.2     2024-01-12 [1] Bioconductor
#>  colorspace         2.1-0     2023-01-23 [1] RSPM (R 4.3.0)
#>  crayon             1.5.2     2022-09-29 [2] RSPM (R 4.3.0)
#>  desc               1.4.3     2023-12-10 [2] RSPM (R 4.3.0)
#>  digest             0.6.33    2023-07-07 [2] RSPM (R 4.3.0)
#>  dplyr              1.1.4     2023-11-17 [1] RSPM (R 4.3.0)
#>  evaluate           0.23      2023-11-01 [2] RSPM (R 4.3.0)
#>  fansi              1.0.6     2023-12-08 [2] RSPM (R 4.3.0)
#>  farver             2.1.1     2022-07-06 [1] RSPM (R 4.3.0)
#>  fastmap            1.1.1     2023-02-24 [2] RSPM (R 4.3.0)
#>  fs                 1.6.3     2023-07-20 [2] RSPM (R 4.3.0)
#>  generics           0.1.3     2022-07-05 [1] RSPM (R 4.3.0)
#>  GenomeInfoDb       1.38.5    2023-12-28 [1] Bioconductor 3.18 (R 4.3.2)
#>  GenomeInfoDbData   1.2.11    2024-01-03 [1] Bioconductor
#>  ggbeeswarm         0.7.2     2023-04-29 [1] RSPM (R 4.3.0)
#>  ggfun              0.1.3     2023-09-15 [1] RSPM (R 4.3.0)
#>  ggplot2            3.4.4     2023-10-12 [1] RSPM (R 4.3.0)
#>  ggplotify          0.1.2     2023-08-09 [1] RSPM (R 4.3.0)
#>  ggtree             3.10.0    2023-10-24 [1] Bioconductor
#>  glue               1.6.2     2022-02-24 [2] RSPM (R 4.3.0)
#>  gridGraphics       0.5-1     2020-12-13 [1] RSPM (R 4.3.0)
#>  gtable             0.3.4     2023-08-21 [1] RSPM (R 4.3.0)
#>  highr              0.10      2022-12-22 [2] RSPM (R 4.3.0)
#>  htmltools          0.5.7     2023-11-03 [2] RSPM (R 4.3.0)
#>  igraph             1.6.0     2023-12-11 [1] RSPM (R 4.3.0)
#>  IRanges            2.36.0    2023-10-24 [1] Bioconductor
#>  jquerylib          0.1.4     2021-04-26 [2] RSPM (R 4.3.0)
#>  jsonlite           1.8.8     2023-12-04 [2] RSPM (R 4.3.0)
#>  knitr              1.45      2023-10-30 [2] RSPM (R 4.3.0)
#>  labeling           0.4.3     2023-08-29 [1] RSPM (R 4.3.0)
#>  lattice            0.22-5    2023-10-24 [3] RSPM (R 4.3.0)
#>  lazyeval           0.2.2     2019-03-15 [1] RSPM (R 4.3.0)
#>  lifecycle          1.0.4     2023-11-07 [2] RSPM (R 4.3.0)
#>  magrittr           2.0.3     2022-03-30 [2] RSPM (R 4.3.0)
#>  memoise            2.0.1     2021-11-26 [2] RSPM (R 4.3.0)
#>  munsell            0.5.0     2018-06-12 [1] RSPM (R 4.3.0)
#>  nlme               3.1-164   2023-11-27 [3] RSPM (R 4.3.0)
#>  patchwork          1.2.0     2024-01-08 [1] RSPM (R 4.3.0)
#>  pillar             1.9.0     2023-03-22 [2] RSPM (R 4.3.0)
#>  pkgconfig          2.0.3     2019-09-22 [2] RSPM (R 4.3.0)
#>  pkgdown            2.0.7     2022-12-14 [2] RSPM (R 4.3.0)
#>  plyr               1.8.9     2023-10-02 [1] RSPM (R 4.3.0)
#>  purrr              1.0.2     2023-08-10 [2] RSPM (R 4.3.0)
#>  R6                 2.5.1     2021-08-19 [2] RSPM (R 4.3.0)
#>  ragg               1.2.7     2023-12-11 [2] RSPM (R 4.3.0)
#>  Rcpp               1.0.11    2023-07-06 [2] RSPM (R 4.3.0)
#>  RCurl              1.98-1.13 2023-11-02 [1] RSPM (R 4.3.0)
#>  reshape2           1.4.4     2020-04-09 [1] RSPM (R 4.3.0)
#>  rlang              1.1.2     2023-11-04 [2] RSPM (R 4.3.0)
#>  rmarkdown          2.25      2023-09-18 [2] RSPM (R 4.3.0)
#>  S4Vectors          0.40.2    2023-11-23 [1] Bioconductor 3.18 (R 4.3.2)
#>  sass               0.4.8     2023-12-06 [2] RSPM (R 4.3.0)
#>  scales             1.3.0     2023-11-28 [1] RSPM (R 4.3.0)
#>  sessioninfo        1.2.2     2021-12-06 [2] RSPM (R 4.3.0)
#>  stringi            1.8.3     2023-12-11 [2] RSPM (R 4.3.0)
#>  stringr            1.5.1     2023-11-14 [2] RSPM (R 4.3.0)
#>  systemfonts        1.0.5     2023-10-09 [2] RSPM (R 4.3.0)
#>  textshaping        0.3.7     2023-10-09 [2] RSPM (R 4.3.0)
#>  tibble             3.2.1     2023-03-20 [2] RSPM (R 4.3.0)
#>  tidyr              1.3.0     2023-01-24 [1] RSPM (R 4.3.0)
#>  tidyselect         1.2.0     2022-10-10 [1] RSPM (R 4.3.0)
#>  tidytree           0.4.6     2023-12-12 [1] RSPM (R 4.3.0)
#>  treeio             1.26.0    2023-10-24 [1] Bioconductor
#>  utf8               1.2.4     2023-10-22 [2] RSPM (R 4.3.0)
#>  vctrs              0.6.5     2023-12-01 [2] RSPM (R 4.3.0)
#>  vipor              0.4.7     2023-12-18 [1] RSPM (R 4.3.0)
#>  withr              2.5.2     2023-10-30 [2] RSPM (R 4.3.0)
#>  xfun               0.41      2023-11-01 [2] RSPM (R 4.3.0)
#>  XVector            0.42.0    2023-10-24 [1] Bioconductor
#>  yaml               2.3.8     2023-12-11 [2] RSPM (R 4.3.0)
#>  yulab.utils        0.1.3     2024-01-08 [1] RSPM (R 4.3.0)
#>  zlibbioc           1.48.0    2023-10-24 [1] Bioconductor
#> 
#>  [1] /__w/_temp/Library
#>  [2] /usr/local/lib/R/site-library
#>  [3] /usr/local/lib/R/library
#> 
#> ──────────────────────────────────────────────────────────────────────────────

References

Paul, Matt, Thomas Carroll, and Doug Barrows. 2021. Herper: The Herper Package Is a Simple Toolset to Install and Manage Conda Packages and Environments from r. https://github.com/RockefellerUniversity/Herper.
Simão, Felipe A, Robert M Waterhouse, Panagiotis Ioannidis, Evgenia V Kriventseva, and Evgeny M Zdobnov. 2015. “BUSCO: Assessing Genome Assembly and Annotation Completeness with Single-Copy Orthologs.” Bioinformatics 31 (19): 3210–12.
Wang, Peng, and Fei Wang. 2022. “A Proposed Metric Set for Evaluation of Genome Assembly Quality.” Trends in Genetics.