Requires the relevant Bioconductor R package to be installed:
convert_to_SE
-convert_to_TSE
-
Arguments
- biom
An rbiom object, such as from
as_rbiom()
. Any value accepted byas_rbiom()
can also be given here.- ...
Not Used.
Details
A SummarizedExperiment object includes counts, metadata, and taxonomy.
TreeSummarizedExperiment additionally includes the tree and sequences.
Examples
library(rbiom)
print(hmp50)
#>
#> ══ Human Microbiome Project - 50 Sample Demo ═══════════════
#>
#> Oral, nasal, vaginal, and fecal samples from a diverse set
#> of healthy volunteers. Source: Human Microbiome Project
#> (<https://hmpdacc.org>).
#>
#> 50 Samples: HMP01, HMP02, HMP03, ..., and HMP50
#> 490 OTUs: Unc01yki, Unc53100, LtbAci52, ...
#> 7 Ranks: .otu, Kingdom, Phylum, ..., and Genus
#> 5 Fields: .sample, Age, BMI, Body Site, and Sex
#> Tree: <present>
#>
#> ── 182 - 22k reads/sample ──────────────────── 2023-09-22 ──
#>
# Requires 'SummarizedExperiment' Bioconductor R package
if (nzchar(system.file(package = "SummarizedExperiment"))) {
se <- convert_to_SE(hmp50)
print(se)
}
#> class: SummarizedExperiment
#> dim: 490 50
#> metadata(0):
#> assays(1): OTU table
#> rownames(490): Unc01yki Unc53100 ... UncOr431 UncTr598
#> rowData names(6): Kingdom Phylum ... Family Genus
#> colnames(50): HMP01 HMP02 ... HMP49 HMP50
#> colData names(4): Age BMI Body Site Sex
# Requires 'TreeSummarizedExperiment' Bioconductor R package
if (nzchar(system.file(package = "TreeSummarizedExperiment"))) {
tse <- convert_to_TSE(hmp50)
print(tse)
}
#> class: TreeSummarizedExperiment
#> dim: 490 50
#> metadata(0):
#> assays(1): OTU table
#> rownames(490): Unc01yki Unc53100 ... UncOr431 UncTr598
#> rowData names(6): Kingdom Phylum ... Family Genus
#> colnames(50): HMP01 HMP02 ... HMP49 HMP50
#> colData names(4): Age BMI Body Site Sex
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
#> rowLinks: a LinkDataFrame (490 rows)
#> rowTree: 1 phylo tree(s) (490 leaves)
#> colLinks: NULL
#> colTree: NULL
#> referenceSeq: a DNAStringSet (490 sequences)