Create a SpatialFeatureExperiment object
Lambda Moses
2024-05-17
Source:vignettes/create_sfe.Rmd
create_sfe.Rmd
library(Voyager)
library(SpatialFeatureExperiment)
library(rjson)
library(Matrix)
library(SFEData)
library(RBioFormats)
The following is reproduced from the SFE vignette.
Visium Space Ranger output
10x Genomics Space Ranger output from a Visium experiment can be read
in a similar manner as in SpatialExperiment
; the
SpatialFeatureExperiment
SFE object has the
spotPoly
column geometry for the spot polygons. If the
filtered matrix (i.e. only spots in the tissue) is read in, then a
column graph called visium
will also be present for the
spatial neighborhood graph of the Visium spots on the tissue. The graph
is not computed if all spots are read in regardless of whether they are
on tissue.
dir <- system.file("extdata", package = "SpatialFeatureExperiment")
sample_ids <- c("sample01", "sample02")
(samples <- file.path(dir, sample_ids, "outs"))
#> [1] "/Users/runner/work/_temp/Library/SpatialFeatureExperiment/extdata/sample01/outs"
#> [2] "/Users/runner/work/_temp/Library/SpatialFeatureExperiment/extdata/sample02/outs"
The results for each tissue capture should be in the
outs
directory. Inside the outs
directory
there are two directories: raw_reature_bc_matrix
has the
unfiltered gene count matrix, and spatial
has the spatial
information.
list.files(samples[1])
#> [1] "filtered_feature_bc_matrix" "spatial"
The DropletUtils
package has a function read10xCounts()
which reads the gene
count matrix. SPE reads in the spatial information, and SFE uses the
spatial information to construct Visium spot polygons and spatial
neighborhood graphs. Inside the spatial
directory:
list.files(file.path(samples[1], "spatial"))
#> [1] "aligned_fiducials.jpg" "barcode_fluorescence_intensity.csv"
#> [3] "detected_tissue_image.jpg" "scalefactors_json.json"
#> [5] "spatial_enrichment.csv" "tissue_hires_image.png"
#> [7] "tissue_lowres_image.png" "tissue_positions.csv"
tissue_lowres_image.png
is a low resolution image of the
tissue.
Inside the scalefactors_json.json
file:
fromJSON(file = file.path(samples[1], "spatial", "scalefactors_json.json"))
#> $tissue_hires_scalef
#> [1] 0.0751202
#>
#> $tissue_lowres_scalef
#> [1] 0.02253606
#>
#> $fiducial_diameter_fullres
#> [1] 304.6116
#>
#> $spot_diameter_fullres
#> [1] 188.5691
spot_diameter_fullres
is the diameter of each Visium
spot in the full resolution H&E image in pixels.
tissue_hires_scalef
and tissue_lowres_scalef
are the ratio of the size of the high resolution (but not full
resolution) and low resolution H&E image to the full resolution
image. fiducial_diameter_fullres
is the diameter of each
fiducial spot used to align the spots to the H&E image in pixels in
the full resolution image.
The tissue_positions_list.csv
file contains information
for the spatial coordinates of the spots and whether each spot is in
tissue as automatically detected by Space Ranger or manually annotated
in the Loupe browser. If the polygon of the tissue boundary is
available, whether from image processing or manual annotation, geometric
operations as supported by the SFE package, which is based on the
sf
package, can be used to find which spots intersect with
the tissue and which spots are contained in the tissue. Geometric
operations can also find the polygons of the intersections between spots
and the tissue, but the results can get messy since the intersections
can have not only polygons but also points and lines.
Now we read in the toy data that is in the Space Ranger output
format. Since Bioconductor version 3.17 (Voyager version 1.2.0), the
image is read as a SpatRaster
object with the terra
package, so it is not loaded into memory unless necessary. When plotting
a large image, it will be downsampled and thus not fully loaded into
memory. The unit can be set with the unit
argument, and can
be either pixels in full resolution image or microns. The latter is
calculated from the former based on spacing between spots, which is
known to be 100 microns.
(sfe3 <- read10xVisiumSFE(samples, dirs = samples, sample_id = sample_ids,
type = "sparse", data = "filtered", images = "lowres",
unit = "full_res_image_pixel"))
#> class: SpatialFeatureExperiment
#> dim: 5 25
#> metadata(0):
#> assays(1): counts
#> rownames(5): ENSG00000014257 ENSG00000142515 ENSG00000263639
#> ENSG00000163810 ENSG00000149591
#> rowData names(14): symbol Feature.Type ...
#> Median.Normalized.Average.Counts_sample02
#> Barcodes.Detected.per.Feature_sample02
#> colnames(25): GTGGCGTGCACCAGAG-1 GGTCCCATAACATAGA-1 ...
#> TGCAATTTGGGCACGG-1 ATGCCAATCGCTCTGC-1
#> colData names(10): in_tissue array_row ... channel3_mean channel3_stdev
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
#> imgData names(4): sample_id image_id data scaleFactor
#>
#> unit: full_res_image_pixel
#> Geometries:
#> colGeometries: spotPoly (POLYGON)
#>
#> Graphs:
#> sample01: col: visium
#> sample02: col: visium
Space Ranger output includes the gene count matrix, spot coordinates, and spot diameter. The Space Ranger output does NOT include nuclei segmentation or pathologist annotation of histological regions. Extra image processing, such as with ImageJ and QuPath, are required for those geometries.
Vizgen MERFISH output
The commercialized MERFISH from Vizgen has a standard output format,
that can be read into SFE with readVizgen()
. Because the
cell segmentation from each field of view (FOV) has a separate HDF5 file
and a MERFISH dataset can have hundreds of FOVs, we strongly recommend
reading the MERFISH output on a server with a large number of CPU cores.
Alternatively, some but not all MERFISH datasets store cell segmentation
in a parquet
file, which can be more easily read into R.
This requires the installation of arrow
. Here we read a toy
dataset which is the first FOV from a real dataset:
fp <- tempdir()
dir_use <- VizgenOutput(file_path = file.path(fp, "vizgen"))
#> see ?SFEData and browseVignettes('SFEData') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> The downloaded files are in /private/var/folders/yv/tw23g49x7kb1jh_s6mf3zvyh0000gn/T/RtmpBKRjJC/vizgen/vizgen_cellbound
list.files(dir_use)
#> [1] "cell_boundaries" "cell_boundaries.parquet"
#> [3] "cell_by_gene.csv" "cell_metadata.csv"
#> [5] "detected_transcripts.csv" "images"
The optional add_molecules
argument can be set to
TRUE
to read in the transcript spots
(sfe_mer <- readVizgen(dir_use, z = 3L, image = "PolyT", add_molecules = TRUE))
#> >>> 1 `.parquet` files exist:
#> /private/var/folders/yv/tw23g49x7kb1jh_s6mf3zvyh0000gn/T/RtmpBKRjJC/vizgen/vizgen_cellbound/cell_boundaries.parquet
#> >>> using -> /private/var/folders/yv/tw23g49x7kb1jh_s6mf3zvyh0000gn/T/RtmpBKRjJC/vizgen/vizgen_cellbound/cell_boundaries.parquet
#> >>> Cell segmentations are found in `.parquet` file
#> Removing 35 cells with area less than 15
#> >>> filtering geometries to match 1023 cells with counts > 0
#> >>> Reading transcript coordinates
#> >>> Converting transcript spots to geometry
#> >>> Writing reformatted transcript spots to disk
#> class: SpatialFeatureExperiment
#> dim: 88 1023
#> metadata(0):
#> assays(1): counts
#> rownames(88): CD4 TLL1 ... Blank-38 Blank-39
#> rowData names(0):
#> colnames(1023): 112824700230101267 112824700230101269 ...
#> 112824700330100848 112824700330100920
#> colData names(11): fov volume ... solidity sample_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : center_x center_y
#> imgData names(4): sample_id image_id data scaleFactor
#>
#> unit: micron
#> Geometries:
#> colGeometries: centroids (POINT), cellSeg (POLYGON)
#> rowGeometries: txSpots (MULTIPOINT)
#>
#> Graphs:
#> sample01:
The unit is always in microns. To make it easier and faster to read the data next time, the processed cell segmentation geometries and transcript spots are written to the same directory where the data resides:
list.files(dir_use)
#> [1] "cell_boundaries" "cell_boundaries.parquet"
#> [3] "cell_by_gene.csv" "cell_metadata.csv"
#> [5] "detected_transcripts.csv" "detected_transcripts.parquet"
#> [7] "images"
10X Xenium output
SFE supports reading the output from Xenium Onboarding Analysis (XOA)
v1 and v2 with the function readXenium()
. Especially for
XOA v2, arrow
is strongly recommended. The cell and nuclei
polygon vertices and transcript spot coordinates are in
parquet
files Similar to readVizgen()
,
readXenium()
makes sf
data frames from the
vertices and transcript spots and saves them as GeoParquet files.
dir_use <- XeniumOutput("v2", file_path = file.path(fp, "xenium"))
#> see ?SFEData and browseVignettes('SFEData') for documentation
#> loading from cache
#> The downloaded files are in /private/var/folders/yv/tw23g49x7kb1jh_s6mf3zvyh0000gn/T/RtmpBKRjJC/xenium/xenium2
list.files(dir_use)
#> [1] "cell_boundaries.csv.gz" "cell_boundaries.parquet"
#> [3] "cell_feature_matrix.h5" "cells.csv.gz"
#> [5] "cells.parquet" "experiment.xenium"
#> [7] "morphology_focus" "nucleus_boundaries.csv.gz"
#> [9] "nucleus_boundaries.parquet" "transcripts.csv.gz"
#> [11] "transcripts.parquet"
# RBioFormats issue: https://github.com/aoles/RBioFormats/issues/42
try(sfe_xen <- readXenium(dir_use, add_molecules = TRUE))
#> >>> Must use gene symbols as row names when adding transcript spots.
#> Error in .jcall(.jcall("RBioFormats", "Lloci/formats/meta/MetadataStore;", :
#> java.lang.NullPointerException: Cannot invoke "loci.formats.DimensionSwapper.setMetadataFiltered(boolean)" because "RBioFormats.reader" is null
(sfe_xen <- readXenium(dir_use, add_molecules = TRUE))
#> >>> Must use gene symbols as row names when adding transcript spots.
#> >>> Cell segmentations are found in `.parquet` file(s)
#> >>> Reading cell and nucleus segmentations
#> >>> Making MULTIPOLYGON nuclei geometries
#> >>> Making POLYGON cell geometries
#> >>> Checking polygon validity
#> >>> Saving geometries to parquet files
#> >>> Reading cell metadata -> `cells.csv`
#> >>> Reading h5 gene count matrix
#> >>> filtering cellSeg geometries to match 6272 cells with counts > 0
#> >>> filtering nucSeg geometries to match 6158 cells with counts > 0
#> >>> Reading transcript coordinates
#> >>> Converting transcript spots to geometry
#> >>> Writing reformatted transcript spots to disk
#> >>> Total of 116 features/genes with no transcript detected or `min_phred` < 20 are removed from SFE object
#> >>> To keep all features -> set `min_phred = NULL`
#> class: SpatialFeatureExperiment
#> dim: 398 6272
#> metadata(1): Samples
#> assays(1): counts
#> rownames(398): ABCC11 ACE2 ... UnassignedCodeword_0488
#> UnassignedCodeword_0497
#> rowData names(3): ID Symbol Type
#> colnames(6272): abclkehb-1 abcnopgp-1 ... odmgoega-1 odmgojlc-1
#> colData names(9): transcript_counts control_probe_counts ...
#> nucleus_area sample_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : x_centroid y_centroid
#> imgData names(4): sample_id image_id data scaleFactor
#>
#> unit: micron
#> Geometries:
#> colGeometries: centroids (POINT), cellSeg (POLYGON), nucSeg (MULTIPOLYGON)
#> rowGeometries: txSpots (MULTIPOINT)
#>
#> Graphs:
#> sample01:
list.files(dir_use)
#> [1] "cell_boundaries_sf.parquet" "cell_boundaries.csv.gz"
#> [3] "cell_boundaries.parquet" "cell_feature_matrix.h5"
#> [5] "cells.csv.gz" "cells.parquet"
#> [7] "experiment.xenium" "morphology_focus"
#> [9] "nucleus_boundaries_sf.parquet" "nucleus_boundaries.csv.gz"
#> [11] "nucleus_boundaries.parquet" "transcripts.csv.gz"
#> [13] "transcripts.parquet" "tx_spots.parquet"
Nanostring CosMX output
This is similar to readVizgen()
and
readXenium()
, except that the output doesn’t come with
images.
dir_use <- CosMXOutput(file_path = file.path(fp, "cosmx"))
#> see ?SFEData and browseVignettes('SFEData') for documentation
#> downloading 1 resources
#> retrieving 1 resource
#> loading from cache
#> The downloaded files are in /private/var/folders/yv/tw23g49x7kb1jh_s6mf3zvyh0000gn/T/RtmpBKRjJC/cosmx/cosmx
list.files(dir_use)
#> [1] "Run5642_S3_Quarter_exprMat_file.csv"
#> [2] "Run5642_S3_Quarter_metadata_file.csv"
#> [3] "Run5642_S3_Quarter_tx_file.csv"
#> [4] "Run5642_S3_Quarter-polygons.csv"
(sfe_cosmx <- readCosMX(dir_use, add_molecules = TRUE))
#> >>> Constructing cell polygons
#> >>> Reading transcript coordinates
#> >>> Converting transcript spots to geometry
#> >>> Writing reformatted transcript spots to disk
#> class: SpatialFeatureExperiment
#> dim: 960 27
#> metadata(0):
#> assays(1): counts
#> rownames(960): Chrna4 Slc6a1 ... NegPrb9 NegPrb10
#> rowData names(0):
#> colnames(27): 367_1 368_1 ... 581_1 583_1
#> colData names(19): fov cell_ID ... Max.DAPI sample_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : CenterX_global_px CenterY_global_px
#> imgData names(0):
#>
#> unit: full_res_image_pixel
#> Geometries:
#> colGeometries: centroids (POINT), cellSeg (POLYGON)
#> rowGeometries: txSpots (MULTIPOINT)
#>
#> Graphs:
#> sample01:
list.files(dir_use)
#> [1] "cell_boundaries_sf.parquet"
#> [2] "Run5642_S3_Quarter_exprMat_file.csv"
#> [3] "Run5642_S3_Quarter_metadata_file.csv"
#> [4] "Run5642_S3_Quarter_tx_file.csv"
#> [5] "Run5642_S3_Quarter-polygons.csv"
#> [6] "tx_spots.parquet"
Other technologies
A read function for Visium HD is in progress. Contribution for Akoya, Molecular Cartography, and Curio Seeker are welcome. See the issues.
Create SFE object from scratch
An SFE object can be constructed from scratch with the assay matrices
and metadata. In this toy example, dgCMatrix
is used, but
since SFE inherits from SingleCellExperiment (SCE), other types of
arrays supported by SCE such as delayed arrays should also work.
# Visium barcode location from Space Ranger
data("visium_row_col")
coords1 <- visium_row_col[visium_row_col$col < 6 & visium_row_col$row < 6,]
coords1$row <- coords1$row * sqrt(3)
# Random toy sparse matrix
set.seed(29)
col_inds <- sample(1:13, 13)
row_inds <- sample(1:5, 13, replace = TRUE)
values <- sample(1:5, 13, replace = TRUE)
mat <- sparseMatrix(i = row_inds, j = col_inds, x = values)
colnames(mat) <- coords1$barcode
rownames(mat) <- sample(LETTERS, 5)
This should be sufficient to create an SPE object, and an SFE object,
even though no sf
data frame was constructed for the
geometries. The constructor behaves similarly to the SPE constructor.
The centroid coordinates of the Visium spots in the example can be
converted into spot polygons with the spotDiameter
argument, which can also be relevant to other technologies with round
spots or beads, such as Slide-seq. Spot diameter in pixels in full
resolution images can be found in the
scalefactors_json.json
file in Space Ranger output.
sfe3 <- SpatialFeatureExperiment(list(counts = mat), colData = coords1,
spatialCoordsNames = c("col", "row"),
spotDiameter = 0.7)
More geometries and spatial graphs can be added after calling the constructor.
Geometries can also be supplied in the constructor.
Coercion from SpatialExperiment
SPE objects can be coerced into SFE objects. If column geometries or spot diameter are not specified, then a column geometry called “centroids” will be created.
spe <- SpatialExperiment::read10xVisium(samples, sample_ids, type = "sparse",
data = "filtered",
images = "hires", load = FALSE)
For the coercion, column names must not be duplicate.
colnames(spe) <- make.unique(colnames(spe), sep = "-")
rownames(spatialCoords(spe)) <- colnames(spe)
(sfe3 <- toSpatialFeatureExperiment(spe))
#> class: SpatialFeatureExperiment
#> dim: 5 25
#> metadata(0):
#> assays(1): counts
#> rownames(5): ENSG00000014257 ENSG00000142515 ENSG00000263639
#> ENSG00000163810 ENSG00000149591
#> rowData names(1): symbol
#> colnames(25): GTGGCGTGCACCAGAG-1 GGTCCCATAACATAGA-1 ...
#> TGCAATTTGGGCACGG-1 ATGCCAATCGCTCTGC-1
#> colData names(4): in_tissue array_row array_col sample_id
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> spatialCoords names(2) : pxl_col_in_fullres pxl_row_in_fullres
#> imgData names(4): sample_id image_id data scaleFactor
#>
#> unit:
#> Geometries:
#> colGeometries: centroids (POINT)
#>
#> Graphs:
#> sample01:
#> sample02:
If images are present in the SPE object, they will be converted into
SpatRaster
when the SPE object is converted into SFE.
Plotting functions in the Voyager
package relies on
SpatRaster
to plot the image behind the geometries.
Coercion from Seurat
Seurat objects can be coerced into SFE objects though coercion from SFE to Seurat is not yet implemented.
dir_extdata <- system.file("extdata", package = "SpatialFeatureExperiment")
obj_vis <- readRDS(file.path(dir_extdata, "seu_vis_toy.rds"))
sfe_conv_vis <-
toSpatialFeatureExperiment(x = obj_vis,
image_scalefactors = "lowres",
unit = "micron",
BPPARAM = BPPARAM)
sfe_conv_vis
Session info
# Clean up
unlink(file.path(fp, "vizgen"), recursive = TRUE)
unlink(file.path(fp, "xenium"), recursive = TRUE)
unlink(file.path(fp, "cosmx"), recursive = TRUE)
sessionInfo()
#> R version 4.4.0 (2024-04-24)
#> Platform: x86_64-apple-darwin20
#> Running under: macOS Ventura 13.6.6
#>
#> Matrix products: default
#> BLAS: /Library/Frameworks/R.framework/Versions/4.4-x86_64/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.4-x86_64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.0
#>
#> locale:
#> [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#>
#> time zone: UTC
#> tzcode source: internal
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] RBioFormats_1.4.0 SFEData_1.6.0
#> [3] Matrix_1.7-0 rjson_0.2.21
#> [5] Voyager_1.6.0 SpatialFeatureExperiment_1.6.1
#>
#> loaded via a namespace (and not attached):
#> [1] jsonlite_1.8.8 wk_0.9.1
#> [3] magrittr_2.0.3 magick_2.8.3
#> [5] rmarkdown_2.27 fs_1.6.4
#> [7] zlibbioc_1.50.0 ragg_1.3.2
#> [9] vctrs_0.6.5 spdep_1.3-3
#> [11] memoise_2.0.1 DelayedMatrixStats_1.26.0
#> [13] RCurl_1.98-1.14 terra_1.7-71
#> [15] htmltools_0.5.8.1 S4Arrays_1.4.0
#> [17] curl_5.2.1 AnnotationHub_3.12.0
#> [19] BiocNeighbors_1.22.0 Rhdf5lib_1.26.0
#> [21] s2_1.1.6 SparseArray_1.4.3
#> [23] rhdf5_2.48.0 sass_0.4.9
#> [25] spData_2.3.0 KernSmooth_2.23-24
#> [27] bslib_0.7.0 htmlwidgets_1.6.4
#> [29] desc_1.4.3 cachem_1.1.0
#> [31] igraph_2.0.3 mime_0.12
#> [33] lifecycle_1.0.4 pkgconfig_2.0.3
#> [35] R6_2.5.1 fastmap_1.2.0
#> [37] GenomeInfoDbData_1.2.12 MatrixGenerics_1.16.0
#> [39] digest_0.6.35 colorspace_2.1-0
#> [41] ggnewscale_0.4.10 AnnotationDbi_1.66.0
#> [43] patchwork_1.2.0 S4Vectors_0.42.0
#> [45] dqrng_0.4.0 RSpectra_0.16-1
#> [47] ExperimentHub_2.12.0 textshaping_0.3.7
#> [49] GenomicRanges_1.56.0 RSQLite_2.3.6
#> [51] beachmat_2.20.0 filelock_1.0.3
#> [53] fansi_1.0.6 httr_1.4.7
#> [55] abind_1.4-5 compiler_4.4.0
#> [57] proxy_0.4-27 withr_3.0.0
#> [59] bit64_4.0.5 tiff_0.1-12
#> [61] BiocParallel_1.38.0 sfarrow_0.4.1
#> [63] DBI_1.2.2 HDF5Array_1.32.0
#> [65] R.utils_2.12.3 rappdirs_0.3.3
#> [67] DelayedArray_0.30.1 classInt_0.4-10
#> [69] bluster_1.14.0 tools_4.4.0
#> [71] units_0.8-5 R.oo_1.26.0
#> [73] glue_1.7.0 dbscan_1.1-12
#> [75] EBImage_4.46.0 rhdf5filters_1.16.0
#> [77] grid_4.4.0 sf_1.0-16
#> [79] cluster_2.1.6 generics_0.1.3
#> [81] memuse_4.2-3 gtable_0.3.5
#> [83] R.methodsS3_1.8.2 class_7.3-22
#> [85] data.table_1.15.4 xml2_1.3.6
#> [87] sp_2.1-4 utf8_1.2.4
#> [89] XVector_0.44.0 BiocGenerics_0.50.0
#> [91] BiocVersion_3.19.1 pillar_1.9.0
#> [93] limma_3.60.0 rJava_1.0-11
#> [95] dplyr_1.1.4 BiocFileCache_2.12.0
#> [97] lattice_0.22-6 bit_4.0.5
#> [99] deldir_2.0-4 tidyselect_1.2.1
#> [101] SingleCellExperiment_1.26.0 locfit_1.5-9.9
#> [103] Biostrings_2.72.0 scuttle_1.14.0
#> [105] sfheaders_0.4.4 knitr_1.46
#> [107] IRanges_2.38.0 edgeR_4.2.0
#> [109] SummarizedExperiment_1.34.0 stats4_4.4.0
#> [111] xfun_0.44 Biobase_2.64.0
#> [113] statmod_1.5.0 DropletUtils_1.24.0
#> [115] matrixStats_1.3.0 UCSC.utils_1.0.0
#> [117] fftwtools_0.9-11 yaml_2.3.8
#> [119] boot_1.3-30 evaluate_0.23
#> [121] codetools_0.2-20 tibble_3.2.1
#> [123] BiocManager_1.30.23 cli_3.6.2
#> [125] arrow_14.0.0.2 systemfonts_1.1.0
#> [127] munsell_0.5.1 jquerylib_0.1.4
#> [129] Rcpp_1.0.12 GenomeInfoDb_1.40.0
#> [131] dbplyr_2.5.0 zeallot_0.1.0
#> [133] png_0.1-8 parallel_4.4.0
#> [135] assertthat_0.2.1 blob_1.2.4
#> [137] pkgdown_2.0.9 ggplot2_3.5.1
#> [139] jpeg_0.1-10 sparseMatrixStats_1.16.0
#> [141] bitops_1.0-7 SpatialExperiment_1.14.0
#> [143] scales_1.3.0 e1071_1.7-14
#> [145] purrr_1.0.2 crayon_1.5.2
#> [147] scico_1.5.0 rlang_1.1.3
#> [149] KEGGREST_1.44.0