The “dar” package is a versatile and user-friendly tool designed to accept inputs in a variety of formats. It primarily utilizes the phyloseq
format but also supports the TreeSummarizedExperiment
format. This flexibility allows users to conduct differential abundance analysis smoothly, irrespective of their initial data format. To facilitate this, a detailed guide is available to aid users in converting other prevalent data formats, such as biome
, mothur
, metaphlan
, and more, into the necessary phyloseq
or TreeSummarizedExperiment
formats.
biome
FormatThe biome
format is a commonly used format in bioinformatics to represent microbiome sequencing data. Here’s how you can import data in biome
format to both phyloseq
and TreeSummarizedExperiment.
To convert data from the biome
format to the phyloseq
format, you can use the phyloseq::import_biom()
function. Here’s a step-by-step example of how to perform this conversion:
# Example of a rich dense biom file
rich_dense_biom <-
system.file("extdata", "rich_dense_otu_table.biom", package = "phyloseq")
# Import biom as a phyloseq-class object
phy <- phyloseq::import_biom(
rich_dense_biom,
parseFunction = parse_taxonomy_greengenes
)
phy
#> phyloseq-class experiment-level object
#> otu_table() OTU Table: [ 5 taxa and 6 samples ]
#> sample_data() Sample Data: [ 6 samples by 4 sample variables ]
#> tax_table() Taxonomy Table: [ 5 taxa by 7 taxonomic ranks ]
# Print sample_data
phyloseq::sample_data(phy)
#> BarcodeSequence LinkerPrimerSequence BODY_SITE Description
#> Sample1 CGCTTATCGAGA CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample2 CATACCAGTAGC CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample3 CTCTCTACCTGT CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample4 CTCTCGGCCTGT CATGCTGCCTCCCGTAGGAGT skin human skin
#> Sample5 CTCTCTACCAAT CATGCTGCCTCCCGTAGGAGT skin human skin
#> Sample6 CTAACTACCAAT CATGCTGCCTCCCGTAGGAGT skin human skin
# Print tax_table
phyloseq::tax_table(phy)
#> Taxonomy Table: [5 taxa by 7 taxonomic ranks]:
#> Kingdom Phylum Class Order
#> GG_OTU_1 "Bacteria" "Proteobacteria" "Gammaproteobacteria" "Enterobacteriales"
#> GG_OTU_2 "Bacteria" "Cyanobacteria" "Nostocophycideae" "Nostocales"
#> GG_OTU_3 "Archaea" "Euryarchaeota" "Methanomicrobia" "Methanosarcinales"
#> GG_OTU_4 "Bacteria" "Firmicutes" "Clostridia" "Halanaerobiales"
#> GG_OTU_5 "Bacteria" "Proteobacteria" "Gammaproteobacteria" "Enterobacteriales"
#> Family Genus Species
#> GG_OTU_1 "Enterobacteriaceae" "Escherichia" NA
#> GG_OTU_2 "Nostocaceae" "Dolichospermum" NA
#> GG_OTU_3 "Methanosarcinaceae" "Methanosarcina" NA
#> GG_OTU_4 "Halanaerobiaceae" "Halanaerobium" "Halanaerobiumsaccharolyticum"
#> GG_OTU_5 "Enterobacteriaceae" "Escherichia" NA
# Recipe init
rec <- dar::recipe(phy, var_info = "BODY_SITE", tax_info = "Genus")
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 5 taxa and 6 samples
#> ℹ variable of interes BODY_SITE (class: character, levels: gut, skin)
#> ℹ taxonomic level Genus
To convert data from the biome
format to the TreeSummarizedExperiment
format, you can use the mia::loadFromBiom()
function. Here’s a step-by-step example of how to perform this conversion:
# Example of a rich dense biom file
rich_dense_biom <-
system.file("extdata", "rich_dense_otu_table.biom", package = "phyloseq")
# Import biom as a phyloseq-class object
tse <- mia::loadFromBiom(rich_dense_biom)
#> Warning in mia::loadFromBiom(rich_dense_biom): 'loadFromBiom' is deprecated.
#> Use 'importBIOM' instead.
tse
#> class: TreeSummarizedExperiment
#> dim: 5 6
#> metadata(0):
#> assays(1): counts
#> rownames(5): GG_OTU_1 GG_OTU_2 GG_OTU_3 GG_OTU_4 GG_OTU_5
#> rowData names(7): taxonomy1 taxonomy2 ... taxonomy6 taxonomy7
#> colnames(6): Sample1 Sample2 ... Sample5 Sample6
#> colData names(4): BarcodeSequence LinkerPrimerSequence BODY_SITE
#> Description
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> rowLinks: NULL
#> rowTree: NULL
#> colLinks: NULL
#> colTree: NULL
# Print sample_data
colData(tse)
#> DataFrame with 6 rows and 4 columns
#> BarcodeSequence LinkerPrimerSequence BODY_SITE Description
#> <character> <character> <character> <character>
#> Sample1 CGCTTATCGAGA CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample2 CATACCAGTAGC CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample3 CTCTCTACCTGT CATGCTGCCTCCCGTAGGAGT gut human gut
#> Sample4 CTCTCGGCCTGT CATGCTGCCTCCCGTAGGAGT skin human skin
#> Sample5 CTCTCTACCAAT CATGCTGCCTCCCGTAGGAGT skin human skin
#> Sample6 CTAACTACCAAT CATGCTGCCTCCCGTAGGAGT skin human skin
# Print tax_table
rowData(tse)
#> DataFrame with 5 rows and 7 columns
#> taxonomy1 taxonomy2 taxonomy3
#> <character> <character> <character>
#> GG_OTU_1 k__Bacteria p__Proteobacteria c__Gammaproteobacteria
#> GG_OTU_2 k__Bacteria p__Cyanobacteria c__Nostocophycideae
#> GG_OTU_3 k__Archaea p__Euryarchaeota c__Methanomicrobia
#> GG_OTU_4 k__Bacteria p__Firmicutes c__Clostridia
#> GG_OTU_5 k__Bacteria p__Proteobacteria c__Gammaproteobacteria
#> taxonomy4 taxonomy5 taxonomy6
#> <character> <character> <character>
#> GG_OTU_1 o__Enterobacteriales f__Enterobacteriaceae g__Escherichia
#> GG_OTU_2 o__Nostocales f__Nostocaceae g__Dolichospermum
#> GG_OTU_3 o__Methanosarcinales f__Methanosarcinaceae g__Methanosarcina
#> GG_OTU_4 o__Halanaerobiales f__Halanaerobiaceae g__Halanaerobium
#> GG_OTU_5 o__Enterobacteriales f__Enterobacteriaceae g__Escherichia
#> taxonomy7
#> <character>
#> GG_OTU_1 s__
#> GG_OTU_2 s__
#> GG_OTU_3 s__
#> GG_OTU_4 s__Halanaerobiumsacc..
#> GG_OTU_5 s__
# Change the column names of the tax_table
colnames(rowData(tse)) <-
c("Kingdom", "Phylum", "Class", "Order", "Family", "Genus", "Species")
rowData(tse)
#> DataFrame with 5 rows and 7 columns
#> Kingdom Phylum Class
#> <character> <character> <character>
#> GG_OTU_1 k__Bacteria p__Proteobacteria c__Gammaproteobacteria
#> GG_OTU_2 k__Bacteria p__Cyanobacteria c__Nostocophycideae
#> GG_OTU_3 k__Archaea p__Euryarchaeota c__Methanomicrobia
#> GG_OTU_4 k__Bacteria p__Firmicutes c__Clostridia
#> GG_OTU_5 k__Bacteria p__Proteobacteria c__Gammaproteobacteria
#> Order Family Genus
#> <character> <character> <character>
#> GG_OTU_1 o__Enterobacteriales f__Enterobacteriaceae g__Escherichia
#> GG_OTU_2 o__Nostocales f__Nostocaceae g__Dolichospermum
#> GG_OTU_3 o__Methanosarcinales f__Methanosarcinaceae g__Methanosarcina
#> GG_OTU_4 o__Halanaerobiales f__Halanaerobiaceae g__Halanaerobium
#> GG_OTU_5 o__Enterobacteriales f__Enterobacteriaceae g__Escherichia
#> Species
#> <character>
#> GG_OTU_1 s__
#> GG_OTU_2 s__
#> GG_OTU_3 s__
#> GG_OTU_4 s__Halanaerobiumsacc..
#> GG_OTU_5 s__
# Recipe init
rec <- dar::recipe(tse, var_info = "BODY_SITE", tax_info = "Genus")
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 5 taxa and 6 samples
#> ℹ variable of interes BODY_SITE (class: character, levels: gut, skin)
#> ℹ taxonomic level Genus
qiime
FormatThe qiime
format is another commonly used format in bioinformatics for microbiome sequencing data. Here’s how you can import data in qiime
format to both Phyloseq
and TreeSummarizedExperiment.
To convert data from the qiime
format to the Phyloseq
format, you can use the phyloseq::import_qiime()
function. Here’s a step-by-step example of how to perform this conversion:
# Import QIIME data
phy_qiime <- phyloseq::import_qiime(
otufilename = system.file("extdata", "GP_otu_table_rand_short.txt.gz", package = "phyloseq"),
mapfilename = system.file("extdata", "master_map.txt", package = "phyloseq"),
treefilename = system.file("extdata", "GP_tree_rand_short.newick.gz", package = "phyloseq")
)
#> Processing map file...
#> Processing otu/tax file...
#> Reading file into memory prior to parsing...
#> Detecting first header line...
#> Header is on line 2
#> Converting input file to a table...
#> Defining OTU table...
#> Parsing taxonomy table...
#> Processing phylogenetic tree...
#> /home/biocbuild/bbs-3.19-bioc/R/site-library/phyloseq/extdata/GP_tree_rand_short.newick.gz ...
phy_qiime
#> phyloseq-class experiment-level object
#> otu_table() OTU Table: [ 500 taxa and 26 samples ]
#> sample_data() Sample Data: [ 26 samples by 7 sample variables ]
#> tax_table() Taxonomy Table: [ 500 taxa by 7 taxonomic ranks ]
#> phy_tree() Phylogenetic Tree: [ 500 tips and 499 internal nodes ]
# Recipe init
rec <- dar::recipe(phy_qiime, var_info = "SampleType", tax_info = "Genus")
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 500 taxa and 26 samples
#> ℹ variable of interes SampleType (class: character, levels: Feces, Freshwater, Freshwater (creek), Mock, Ocean, Sediment (estuary), Skin, Soil, Tongue)
#> ℹ taxonomic level Genus
To convert data from the qiime
format to the TreeSummarizedExperiment
format, you can use the mia::loadFromQIIME2()
function. Here’s a step-by-step example of how to perform this conversion:
# Import QIIME data to tse
tse_qiime <- mia::loadFromQIIME2(
featureTableFile = system.file("extdata", "table.qza", package = "mia"),
taxonomyTableFile = system.file("extdata", "taxonomy.qza", package = "mia"),
sampleMetaFile = system.file("extdata", "sample-metadata.tsv", package = "mia"),
refSeqFile = system.file("extdata", "refseq.qza", package = "mia"),
phyTreeFile = system.file("extdata", "tree.qza", package = "mia")
)
#> Warning in mia::loadFromQIIME2(featureTableFile = system.file("extdata", :
#> 'loadFromQIIME2' is deprecated. Use 'importQIIME2' instead.
tse_qiime
#> class: TreeSummarizedExperiment
#> dim: 770 34
#> metadata(0):
#> assays(1): counts
#> rownames(770): 4b5eeb300368260019c1fbc7a3c718fc
#> fe30ff0f71a38a39cf1717ec2be3a2fc ... 98d250a339a635f20e26397dafc6ced3
#> 1830c14ead81ad012f1db0e12f8ab6a4
#> rowData names(8): Kingdom Phylum ... Species Confidence
#> colnames(34): L1S8 L1S57 ... L6S68 L6S93
#> colData names(9): sample.id barcode.sequence ...
#> reported.antibiotic.usage days.since.experiment.start
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> rowLinks: a LinkDataFrame (770 rows)
#> rowTree: 1 phylo tree(s) (770 leaves)
#> colLinks: NULL
#> colTree: NULL
#> referenceSeq: a DNAStringSet (770 sequences)
# Recipe init
rec <- dar::recipe(tse_qiime, var_info = "body.site", tax_info = "Genus")
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 770 taxa and 34 samples
#> ℹ variable of interes body.site (class: character, levels: gut, left palm, right palm, tongue)
#> ℹ taxonomic level Genus
mothur
FormatThe mothur
format is another commonly used format in bioinformatics for microbiome sequencing data. Here’s how you can import data in mothur
format to both Phyloseq
and TreeSummarizedExperiment.
To convert data from the mothur
format to the Phyloseq
format, you can use the phyloseq::import_mothur()
function. Here’s a step-by-step example of how to perform this conversion:
# Import Mothur data
phy_mothur <- phyloseq::import_mothur(
mothur_list_file = system.file("extdata", "esophagus.fn.list.gz", package = "phyloseq"),
mothur_group_file = system.file("extdata", "esophagus.good.groups.gz", package = "phyloseq"),
mothur_tree_file = system.file("extdata", "esophagus.tree.gz", package = "phyloseq")
)
phy_mothur
#> phyloseq-class experiment-level object
#> otu_table() OTU Table: [ 591 taxa and 3 samples ]
#> phy_tree() Phylogenetic Tree: [ 591 tips and 590 internal nodes ]
# Recipe init
rec <- dar::recipe(phy_mothur)
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 591 taxa and 3 samples
#> ✖ undefined variable of interest. Use add_var() to add it to Recipe!
#> ✖ undefined taxonomic level. Use add_tax() to add it to Recipe!
To convert data from the mothur
format to the TreeSummarizedExperiment
format, you can use the mia::loadFromMothur()
function. Here’s a step-by-step example of how to perform this conversion:
# Import Mothur data to TreeSummarizedExperiment
tse_mothur <- mia::loadFromMothur(
sharedFile = system.file("extdata", "mothur_example.shared", package = "mia"),
taxonomyFile = system.file("extdata", "mothur_example.cons.taxonomy", package = "mia"),
designFile = system.file("extdata", "mothur_example.design", package = "mia")
) |> methods::as("TreeSummarizedExperiment")
#> Warning in mia::loadFromMothur(sharedFile = system.file("extdata",
#> "mothur_example.shared", : 'loadFromMothur' is deprecated. Use 'importMothur'
#> instead.
tse_mothur
#> class: TreeSummarizedExperiment
#> dim: 100 100
#> metadata(0):
#> assays(1): counts
#> rownames(100): Otu001 Otu002 ... Otu099 Otu100
#> rowData names(8): OTU Size ... Family Genus
#> colnames(100): Sample1 Sample2 ... Sample99 Sample100
#> colData names(7): group sex ... numOtus Group
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> rowLinks: NULL
#> rowTree: NULL
#> colLinks: NULL
#> colTree: NULL
# Recipe init
rec <- dar::recipe(tse_mothur, var_info = "drug", tax_info = "Genus")
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 100 taxa and 100 samples
#> ℹ variable of interes drug (class: character, levels: A, B)
#> ℹ taxonomic level Genus
metaphlan
FormatThe metaphlan
format is another commonly used format in bioinformatics for microbiome sequencing data. Here’s how you can import data in metaphlan
format to TreeSummarizedExperiment.
To convert data from the metaphlan
format to the TreeSummarizedExperiment
format, you can use the mia::loadFromMetaphlan()
function. Here’s a step-by-step example of how to perform this conversion:
# Importing data from Metaphlan
tse_metaphlan <- mia::loadFromMetaphlan(
file = system.file("extdata", "merged_abundance_table.txt", package = "mia")
)
#> Warning in mia::loadFromMetaphlan(file = system.file("extdata",
#> "merged_abundance_table.txt", : 'loadFromMetaphlan' is deprecated. Use
#> 'importMetaPhlAn' instead.
# Recipe init
rec <- dar::recipe(tse_metaphlan)
rec
#> ── DAR Recipe ──────────────────────────────────────────────────────────────────
#> Inputs:
#>
#> ℹ phyloseq object with 16 taxa and 6 samples
#> ✖ undefined variable of interest. Use add_var() to add it to Recipe!
#> ✖ undefined taxonomic level. Use add_tax() to add it to Recipe!
In this guide, we have explored various methods for importing microbiome sequencing data from different formats into Phyloseq
and TreeSummarizedExperiment.
We’ve covered the biome
, qiime
, mothur
, metaphlan
, and humann
formats, providing step-by-step examples for each.
The flexibility of these tools allows for a smooth transition between different data formats, making it easier to conduct your analysis irrespective of the initial data format. By following the steps outlined in this guide, you should be able to successfully convert your data and carry out further differential abundance analysis.
Remember, the specific details of your data may require you to adjust the parameters in the import functions. Always inspect your data after conversion to ensure it has been imported correctly.
devtools::session_info()
#> ─ Session info ───────────────────────────────────────────────────────────────
#> setting value
#> version R version 4.4.0 beta (2024-04-15 r86425)
#> os Ubuntu 22.04.4 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate en_US.UTF-8
#> ctype en_US.UTF-8
#> tz America/New_York
#> date 2024-04-30
#> pandoc 2.7.3 @ /usr/bin/ (via rmarkdown)
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-5 2016-07-21 [2] CRAN (R 4.4.0)
#> ade4 1.7-22 2023-02-06 [2] CRAN (R 4.4.0)
#> ape 5.8 2024-04-11 [2] CRAN (R 4.4.0)
#> assertthat 0.2.1 2019-03-21 [2] CRAN (R 4.4.0)
#> beachmat 2.20.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> beeswarm 0.4.0 2021-06-01 [2] CRAN (R 4.4.0)
#> Biobase * 2.64.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> BiocGenerics * 0.50.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> BiocNeighbors 1.22.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> BiocParallel 1.38.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> BiocSingular 1.20.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> biomformat 1.32.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> Biostrings * 2.72.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> bluster 1.14.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> brio 1.1.5 2024-04-24 [2] CRAN (R 4.4.0)
#> bslib 0.7.0 2024-03-29 [2] CRAN (R 4.4.0)
#> ca 0.71.1 2020-01-24 [2] CRAN (R 4.4.0)
#> cachem 1.0.8 2023-05-01 [2] CRAN (R 4.4.0)
#> cli 3.6.2 2023-12-11 [2] CRAN (R 4.4.0)
#> cluster 2.1.6 2023-12-01 [3] CRAN (R 4.4.0)
#> codetools 0.2-20 2024-03-31 [3] CRAN (R 4.4.0)
#> colorspace 2.1-0 2023-01-23 [2] CRAN (R 4.4.0)
#> crayon 1.5.2 2022-09-29 [2] CRAN (R 4.4.0)
#> crosstalk 1.2.1 2023-11-23 [2] CRAN (R 4.4.0)
#> dar * 1.0.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> data.table 1.15.4 2024-03-30 [2] CRAN (R 4.4.0)
#> DBI 1.2.2 2024-02-16 [2] CRAN (R 4.4.0)
#> DECIPHER 3.0.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> decontam 1.24.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> DelayedArray 0.30.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> DelayedMatrixStats 1.26.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> dendextend 1.17.1 2023-03-25 [2] CRAN (R 4.4.0)
#> devtools 2.4.5 2022-10-11 [2] CRAN (R 4.4.0)
#> digest 0.6.35 2024-03-11 [2] CRAN (R 4.4.0)
#> DirichletMultinomial 1.46.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.4.0)
#> ellipsis 0.3.2 2021-04-29 [2] CRAN (R 4.4.0)
#> evaluate 0.23 2023-11-01 [2] CRAN (R 4.4.0)
#> fansi 1.0.6 2023-12-08 [2] CRAN (R 4.4.0)
#> farver 2.1.1 2022-07-06 [2] CRAN (R 4.4.0)
#> fastmap 1.1.1 2023-02-24 [2] CRAN (R 4.4.0)
#> foreach 1.5.2 2022-02-02 [2] CRAN (R 4.4.0)
#> fs 1.6.4 2024-04-25 [2] CRAN (R 4.4.0)
#> furrr 0.3.1 2022-08-15 [2] CRAN (R 4.4.0)
#> future 1.33.2 2024-03-26 [2] CRAN (R 4.4.0)
#> generics 0.1.3 2022-07-05 [2] CRAN (R 4.4.0)
#> GenomeInfoDb * 1.40.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> GenomeInfoDbData 1.2.12 2024-04-16 [2] Bioconductor
#> GenomicRanges * 1.56.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> ggbeeswarm 0.7.2 2023-04-29 [2] CRAN (R 4.4.0)
#> ggplot2 3.5.1 2024-04-23 [2] CRAN (R 4.4.0)
#> ggrepel 0.9.5 2024-01-10 [2] CRAN (R 4.4.0)
#> globals 0.16.3 2024-03-08 [2] CRAN (R 4.4.0)
#> glue 1.7.0 2024-01-09 [2] CRAN (R 4.4.0)
#> gridExtra 2.3 2017-09-09 [2] CRAN (R 4.4.0)
#> gtable 0.3.5 2024-04-22 [2] CRAN (R 4.4.0)
#> heatmaply 1.5.0 2023-10-06 [2] CRAN (R 4.4.0)
#> highr 0.10 2022-12-22 [2] CRAN (R 4.4.0)
#> htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.4.0)
#> htmlwidgets 1.6.4 2023-12-06 [2] CRAN (R 4.4.0)
#> httpuv 1.6.15 2024-03-26 [2] CRAN (R 4.4.0)
#> httr 1.4.7 2023-08-15 [2] CRAN (R 4.4.0)
#> igraph 2.0.3 2024-03-13 [2] CRAN (R 4.4.0)
#> IRanges * 2.38.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> irlba 2.3.5.1 2022-10-03 [2] CRAN (R 4.4.0)
#> iterators 1.0.14 2022-02-05 [2] CRAN (R 4.4.0)
#> jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.4.0)
#> jsonlite 1.8.8 2023-12-04 [2] CRAN (R 4.4.0)
#> knitr 1.46 2024-04-06 [2] CRAN (R 4.4.0)
#> labeling 0.4.3 2023-08-29 [2] CRAN (R 4.4.0)
#> later 1.3.2 2023-12-06 [2] CRAN (R 4.4.0)
#> lattice 0.22-6 2024-03-20 [3] CRAN (R 4.4.0)
#> lazyeval 0.2.2 2019-03-15 [2] CRAN (R 4.4.0)
#> lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.4.0)
#> listenv 0.9.1 2024-01-29 [2] CRAN (R 4.4.0)
#> magrittr 2.0.3 2022-03-30 [2] CRAN (R 4.4.0)
#> MASS 7.3-60.2 2024-04-16 [3] local
#> Matrix 1.7-0 2024-03-22 [3] CRAN (R 4.4.0)
#> MatrixGenerics * 1.16.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> matrixStats * 1.3.0 2024-04-11 [2] CRAN (R 4.4.0)
#> memoise 2.0.1 2021-11-26 [2] CRAN (R 4.4.0)
#> mgcv 1.9-1 2023-12-21 [3] CRAN (R 4.4.0)
#> mia * 1.12.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> microbiome 1.26.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> mime 0.12 2021-09-28 [2] CRAN (R 4.4.0)
#> miniUI 0.1.1.1 2018-05-18 [2] CRAN (R 4.4.0)
#> MultiAssayExperiment * 1.30.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> multtest 2.60.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> munsell 0.5.1 2024-04-01 [2] CRAN (R 4.4.0)
#> nlme 3.1-164 2023-11-27 [3] CRAN (R 4.4.0)
#> parallelly 1.37.1 2024-02-29 [2] CRAN (R 4.4.0)
#> permute 0.9-7 2022-01-27 [2] CRAN (R 4.4.0)
#> phyloseq * 1.48.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> pillar 1.9.0 2023-03-22 [2] CRAN (R 4.4.0)
#> pkgbuild 1.4.4 2024-03-17 [2] CRAN (R 4.4.0)
#> pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.4.0)
#> pkgload 1.3.4 2024-01-16 [2] CRAN (R 4.4.0)
#> plotly 4.10.4 2024-01-13 [2] CRAN (R 4.4.0)
#> plyr 1.8.9 2023-10-02 [2] CRAN (R 4.4.0)
#> profvis 0.3.8 2023-05-02 [2] CRAN (R 4.4.0)
#> promises 1.3.0 2024-04-05 [2] CRAN (R 4.4.0)
#> purrr 1.0.2 2023-08-10 [2] CRAN (R 4.4.0)
#> R6 2.5.1 2021-08-19 [2] CRAN (R 4.4.0)
#> RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.4.0)
#> Rcpp 1.0.12 2024-01-09 [2] CRAN (R 4.4.0)
#> registry 0.5-1 2019-03-05 [2] CRAN (R 4.4.0)
#> remotes 2.5.0 2024-03-17 [2] CRAN (R 4.4.0)
#> reshape2 1.4.4 2020-04-09 [2] CRAN (R 4.4.0)
#> rhdf5 2.48.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> rhdf5filters 1.16.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> Rhdf5lib 1.26.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> rlang 1.1.3 2024-01-10 [2] CRAN (R 4.4.0)
#> rmarkdown 2.26 2024-03-05 [2] CRAN (R 4.4.0)
#> rsvd 1.0.5 2021-04-16 [2] CRAN (R 4.4.0)
#> Rtsne 0.17 2023-12-07 [2] CRAN (R 4.4.0)
#> S4Arrays 1.4.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> S4Vectors * 0.42.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> sass 0.4.9 2024-03-15 [2] CRAN (R 4.4.0)
#> ScaledMatrix 1.12.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> scales 1.3.0 2023-11-28 [2] CRAN (R 4.4.0)
#> scater 1.32.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> scuttle 1.14.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> seriation 1.5.5 2024-04-17 [2] CRAN (R 4.4.0)
#> sessioninfo 1.2.2 2021-12-06 [2] CRAN (R 4.4.0)
#> shiny 1.8.1.1 2024-04-02 [2] CRAN (R 4.4.0)
#> SingleCellExperiment * 1.26.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> SparseArray 1.4.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> sparseMatrixStats 1.16.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> stringi 1.8.3 2023-12-11 [2] CRAN (R 4.4.0)
#> stringr 1.5.1 2023-11-14 [2] CRAN (R 4.4.0)
#> SummarizedExperiment * 1.34.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> survival 3.6-4 2024-04-24 [3] CRAN (R 4.4.0)
#> testthat 3.2.1.1 2024-04-14 [2] CRAN (R 4.4.0)
#> tibble 3.2.1 2023-03-20 [2] CRAN (R 4.4.0)
#> tidyr 1.3.1 2024-01-24 [2] CRAN (R 4.4.0)
#> tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.4.0)
#> tidytree 0.4.6 2023-12-12 [2] CRAN (R 4.4.0)
#> treeio 1.28.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> TreeSummarizedExperiment * 2.12.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> TSP 1.2-4 2023-04-04 [2] CRAN (R 4.4.0)
#> UCSC.utils 1.0.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> UpSetR 1.4.0 2019-05-22 [2] CRAN (R 4.4.0)
#> urlchecker 1.0.1 2021-11-30 [2] CRAN (R 4.4.0)
#> usethis 2.2.3 2024-02-19 [2] CRAN (R 4.4.0)
#> utf8 1.2.4 2023-10-22 [2] CRAN (R 4.4.0)
#> vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.4.0)
#> vegan 2.6-4 2022-10-11 [2] CRAN (R 4.4.0)
#> vipor 0.4.7 2023-12-18 [2] CRAN (R 4.4.0)
#> viridis 0.6.5 2024-01-29 [2] CRAN (R 4.4.0)
#> viridisLite 0.4.2 2023-05-02 [2] CRAN (R 4.4.0)
#> webshot 0.5.5 2023-06-26 [2] CRAN (R 4.4.0)
#> withr 3.0.0 2024-01-16 [2] CRAN (R 4.4.0)
#> xfun 0.43 2024-03-25 [2] CRAN (R 4.4.0)
#> xtable 1.8-4 2019-04-21 [2] CRAN (R 4.4.0)
#> XVector * 0.44.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#> yaml 2.3.8 2023-12-11 [2] CRAN (R 4.4.0)
#> yulab.utils 0.1.4 2024-01-28 [2] CRAN (R 4.4.0)
#> zlibbioc 1.50.0 2024-04-30 [2] Bioconductor 3.19 (R 4.4.0)
#>
#> [1] /tmp/Rtmpa31KbR/Rinst2bd145554a745f
#> [2] /home/biocbuild/bbs-3.19-bioc/R/site-library
#> [3] /home/biocbuild/bbs-3.19-bioc/R/library
#>
#> ──────────────────────────────────────────────────────────────────────────────