Title: | Read/Write, Transform, and Summarize 'BIOM' Data |
---|---|
Description: | A toolkit for working with Biological Observation Matrix ('BIOM') files. Features include reading/writing all 'BIOM' formats, rarefaction, alpha diversity, beta diversity (including 'UniFrac'), summarizing counts by taxonomic level, and sample subsetting. Standalone functions for reading, writing, and subsetting phylogenetic trees are also provided. All CPU intensive operations are encoded in C with multi-thread support. |
Authors: | Daniel P. Smith [aut, cre, cph] |
Maintainer: | Daniel P. Smith <[email protected]> |
License: | AGPL-3 |
Version: | 1.0.3 |
Built: | 2024-08-16 06:37:47 UTC |
Source: | CRAN |
Estimate the diversity of each sample.
alpha.div(biom, rarefy = FALSE)
alpha.div(biom, rarefy = FALSE)
biom |
A |
rarefy |
Control how/whether rarefactions are done prior to alpha diversity computations. Options are:
|
A data frame of four diversity values for each sample in
biom
. The column names are Sample, Depth and the
diversity metrics: OTUs, Shannon, Chao1, Simpson,
and InvSimpson. The row names are the sample names, except when
multiple rarefactions are done.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ad <- alpha.div(biom) head(ad) biom <- subset(biom, `Body Site` == "Saliva" & Age < 26) ad <- alpha.div(biom, "multi") boxplot(Shannon ~ Depth, data=ad, xlab="Reads", ylab="Diversity")
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ad <- alpha.div(biom) head(ad) biom <- subset(biom, `Body Site` == "Saliva" & Age < 26) ad <- alpha.div(biom, "multi") boxplot(Shannon ~ Depth, data=ad, xlab="Reads", ylab="Diversity")
Make a distance matrix of samples vs samples.
beta.div(biom, method, weighted = TRUE, tree = NULL)
beta.div(biom, method, weighted = TRUE, tree = NULL)
biom |
A |
method |
The distance algorithm to use. Options are:
“manhattan”, “euclidean”,
“bray-curtis”, “jaccard”, and
“unifrac”. Non-ambiguous abbrevations of the method
names are also accepted. A phylogentic tree must be present in
|
weighted |
Take relative abundances into account. When
|
tree |
A |
A distance matrix.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) biom <- select(biom, 1:10) dm <- beta.div(biom, 'unifrac') as.matrix(dm)[1:4,1:4] plot(hclust(dm))
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) biom <- select(biom, 1:10) dm <- beta.div(biom, 'unifrac') as.matrix(dm)[1:4,1:4] plot(hclust(dm))
Get the abundance counts.
counts(biom)
counts(biom)
biom |
A |
A numeric matrix of the sample abundance counts in biom
.
Other accessor functions:
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) counts(biom)[1:4,1:5]
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) counts(biom)[1:4,1:5]
Get biom's misc information.
info(biom)
info(biom)
biom |
A |
A data frame of the metadata in biom
.
Other accessor functions:
counts()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) info(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) info(biom)
Get the sample metadata.
metadata(biom)
metadata(biom)
biom |
A |
A data frame of the metadata in biom
.
Other accessor functions:
counts()
,
info()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) metadata(biom)[1:4,1:3]
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) metadata(biom)[1:4,1:3]
Number of samples in a BIOM.
nsamples(biom)
nsamples(biom)
biom |
A |
The number of samples present.
Other accessor functions:
counts()
,
info()
,
metadata()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) nsamples(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) nsamples(biom)
Number of taxa in a BIOM.
ntaxa(biom)
ntaxa(biom)
biom |
A |
The number of taxa present.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ntaxa(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ntaxa(biom)
Get the phylogenetic tree.
phylogeny(biom)
phylogeny(biom)
biom |
A |
A phylo
class object of the tree in biom
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) summary(phylogeny(biom))
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) summary(phylogeny(biom))
Summarize the contents of a BIOM object
## S3 method for class 'BIOM' print(x, ...)
## S3 method for class 'BIOM' print(x, ...)
x |
A BIOM object, as returned from read.biom. |
... |
Not used. |
NULL (invisibly)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) print(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) print(biom)
Subset counts so that all samples have the same number of observations.
rarefy(biom, depth = NULL, seed = 0)
rarefy(biom, depth = NULL, seed = 0)
biom |
A |
depth |
The number of observations to keep, per sample. If set to
|
seed |
An integer to use for seeding the random number generator. If
you need to create different random rarefactions of the same |
A matrix
, simple_triplet_matrix
, or BIOM
object, depending on the input object type. The type of object provided
is the same type that is returned. The retained observations are randomly
selected, based on a seed value derived from the BIOM
object.
Therefore, rarefying the same biom to the same depth will always produce
the same resultant rarification.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) range(slam::col_sums(biom$counts)) biom <- rarefy(biom, depth=1000) range(slam::col_sums(biom$counts))
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) range(slam::col_sums(biom$counts)) biom <- rarefy(biom, depth=1000) range(slam::col_sums(biom$counts))
A toolkit for working with Biological Observation Matrix (BIOM) files. Features include reading/writing all BIOM formats, rarefaction, alpha diversity, beta diversity (including UniFrac), summarizing counts by taxonomic level, and sample subsetting. Standalone functions for reading, writing, and subsetting phylogenetic trees are also provided. All CPU intensive operations are encoded in C with multi-thread support.
Many rbiom
functions support multithreading:
The default behavior of these function is to run on as many cores as are
available in the local compute environment. If you wish to limit the number
of simultaneous threads, set RcppParallel
's numThreads
option.
For instance:
RcppParallel::setThreadOptions(numThreads = 4)
Extracts counts, metadata, taxonomy, and phylogeny from a biom file.
read.biom(src, tree = "auto", prune = FALSE)
read.biom(src, tree = "auto", prune = FALSE)
src |
Input data as either a file path, URL, or JSON string.
|
tree |
The default value of |
prune |
Should samples and taxa with zero observations be discarded? |
A BIOM
class object containing the parsed data. This object
can be treated as a list with the following named elements:
A numeric slam
sparse matrix of observation
counts. Taxa (OTUs) as rows and samples as columns.
A data frame containing any embedded metadata. Row names are sample IDs.
Character matrix of taxonomic names, if given. Row names are taxa (OTU) IDs. Column rows are named Kingdom, Phylum, Class, Order, Family, Genus, Species, and Strain, or TaxLvl.1, TaxLvl.2, ... , TaxLvl.N when more than 8 levels of taxonomy are encoded in the biom file.
An object of class phylo
defining the
phylogenetic relationships between the taxa. Although the
official specification for BIOM only includes phylogenetic trees
in BIOM version 2.1, if a BIOM version 1.0 file includes a
phylogeny
entry with newick data, then it will be loaded
here as well. The ape package has additional functions
for working with phylo
objects.
A named character vector, where the names are taxonomic identifiers and the values are the sequences they represent. These values are not part of the official BIOM specification, but will be read and written when defined.
A list of other attributes defined in the BIOM file,
such as id
, type
, format
, format_url
,
generated_by
, date
, matrix_type
,
matrix_element_type
, Comment
, and shape
metadata
, taxonomy
, and phylogeny
are optional
components of the BIOM file specification and therefore will be empty
in the returned object when they are not provided by the BIOM file.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) summary(biom) # Taxa Abundances as.matrix(biom$counts[1:4,1:4]) top5 <- names(head(rev(sort(slam::row_sums(biom$counts))), 5)) biom$taxonomy[top5,c('Family', 'Genus')] as.matrix(biom$counts[top5, 1:6]) # Metadata table(biom$metadata$Sex, biom$metadata$`Body Site`) sprintf("Mean age: %.1f", mean(biom$metadata$Age)) # Phylogenetic tree tree <- biom$phylogeny top5.tree <- rbiom::subtree(tree, top5) ape::plot.phylo(top5.tree)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) summary(biom) # Taxa Abundances as.matrix(biom$counts[1:4,1:4]) top5 <- names(head(rev(sort(slam::row_sums(biom$counts))), 5)) biom$taxonomy[top5,c('Family', 'Genus')] as.matrix(biom$counts[top5, 1:6]) # Metadata table(biom$metadata$Sex, biom$metadata$`Body Site`) sprintf("Mean age: %.1f", mean(biom$metadata$Age)) # Phylogenetic tree tree <- biom$phylogeny top5.tree <- rbiom::subtree(tree, top5) ape::plot.phylo(top5.tree)
Parse a fasta file into a named character vector.
read.fasta(file, ids = NULL)
read.fasta(file, ids = NULL)
file |
A file with fasta-formatted sequences. Can optionally be compressed with gzip, bzip2, xz, or lzma. |
ids |
Character vector of IDs to retrieve. The default, |
A named character vector in which names are the fasta headers and values are the sequences.
A phylogenetic tree is required for computing UniFrac distance matrices. You can load a tree either from a file or by providing the tree string directly. This tree must be in Newick format, also known as parenthetic format and New Hampshire format.
read.tree(src)
read.tree(src)
src |
Input data as either a file path, URL, or Newick string. URLs
must begin with http://, https://, ftp://, or
ftps://. Newick strings must have |
A phylo
class object representing the tree.
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) tree <- read.tree(" (t9:0.99,((t5:0.87,t2:0.89):0.51,(((t10:0.16,(t7:0.83,t4:0.96) :0.94):0.69,(t6:0.92,(t3:0.62,t1:0.85):0.54):0.23):0.74,t8:0.1 2):0.43):0.67);")
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) tree <- read.tree(" (t9:0.99,((t5:0.87,t2:0.89):0.51,(((t10:0.16,(t7:0.83,t4:0.96) :0.94):0.69,(t6:0.92,(t3:0.62,t1:0.85):0.54):0.23):0.74,t8:0.1 2):0.43):0.67);")
Get the sample names.
sample.names(biom)
sample.names(biom)
biom |
A |
A character vector of the sample IDs / names in biom
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sequences()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) sample.names(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) sample.names(biom)
Reduce samples to a specific list
select(biom, samples = NULL, nTop = NULL, nRandom = NULL, seed = 0)
select(biom, samples = NULL, nTop = NULL, nRandom = NULL, seed = 0)
biom |
A BIOM object, as returned from read.biom. |
samples |
Sample names, indices, or a logical vector identifying
the samples to keep. The latter two should be based on the order of
sample names given by |
nTop |
Selects this number of samples, taking the sample with the most
observations first, then the sample with the second-most observations,
etc. If |
nRandom |
Randomly selects this number of samples. If higher than the number of samples in the dataset, the entire dataset will be returned. See note. |
seed |
Random seed, used when selecting Note: Generally, you will specify only one of the filters: |
A BIOM
object.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ex1 <- select(biom, c('HMP14', 'HMP22', 'HMP03')) ex2 <- select(biom, c(32, 11, 28, 16, 46, 5)) ex3 <- select(biom, 1:50 %% 6 == 0) ex4 <- select(biom, nRandom = 10) ex5 <- select(biom, nTop = 5) ex6 <- select(biom, samples = 10:40, nTop = 20, nRandom = 10)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ex1 <- select(biom, c('HMP14', 'HMP22', 'HMP03')) ex2 <- select(biom, c(32, 11, 28, 16, 46, 5)) ex3 <- select(biom, 1:50 %% 6 == 0) ex4 <- select(biom, nRandom = 10) ex5 <- select(biom, nTop = 5) ex6 <- select(biom, samples = 10:40, nTop = 20, nRandom = 10)
DNA sequence associated with each taxonomic identifier.
sequences(biom)
sequences(biom)
biom |
A |
A named character vector of sequences in biom
. If this data
is not present, then returns NULL
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
taxa.names()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) sequences(biom)[1:4] # Write to a compressed fasta file in the temporary directory: seqs <- sequences(biom) conn <- bzfile(file.path(tempdir(), "Sequences.fa.bz2"), "w") cat(sprintf(">%s\n%s", names(seqs), seqs), file=conn, sep="\n") close(conn) # You can also use the write.fasta function for this task: write.fasta(biom, file.path(tempdir(), "Sequences.fa.gz"))
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) sequences(biom)[1:4] # Write to a compressed fasta file in the temporary directory: seqs <- sequences(biom) conn <- bzfile(file.path(tempdir(), "Sequences.fa.bz2"), "w") cat(sprintf(">%s\n%s", names(seqs), seqs), file=conn, sep="\n") close(conn) # You can also use the write.fasta function for this task: write.fasta(biom, file.path(tempdir(), "Sequences.fa.gz"))
Subset samples using the BIOM object's metadata
## S3 method for class 'BIOM' subset(x, ...)
## S3 method for class 'BIOM' subset(x, ...)
x |
A BIOM object, as returned from read.biom. |
... |
Test to run on the metadata to identify samples to retain. |
A BIOM
object.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ex1 <- subset(biom, Age > 30) ex2 <- subset(biom, `Body Site` %in% c("Saliva", "Stool")) ex3 <- subset(biom, Age < 25 & BMI > 22)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) ex1 <- subset(biom, Age > 30) ex2 <- subset(biom, `Body Site` %in% c("Saliva", "Stool")) ex3 <- subset(biom, Age < 25 & BMI > 22)
Create a subtree by specifying tips to keep.
subtree(tree, tips)
subtree(tree, tips)
tree |
A phylo object, as returned from read.tree. |
tips |
A character, numeric, or logical vector of tips to keep. |
A phylo
object for the subtree.
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) leafs <- tips(tree) subtree <- subtree(tree, head(leafs))
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) leafs <- tips(tree) subtree <- subtree(tree, head(leafs))
Get the taxa names.
taxa.names(biom)
taxa.names(biom)
biom |
A |
A character vector of the taxa IDs / names in biom
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.ranks()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxa.names(biom) %>% head()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxa.names(biom) %>% head()
Get the taxa ranks.
taxa.ranks(biom)
taxa.ranks(biom)
biom |
A |
A character vector of the taxa ranks in biom
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxonomy()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxa.ranks(biom)
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxa.ranks(biom)
Generate a matrix of samples by taxa, at the specified taxonomic rank.
taxa.rollup(biom, rank = "OTU", map = NULL, lineage = FALSE, sparse = FALSE)
taxa.rollup(biom, rank = "OTU", map = NULL, lineage = FALSE, sparse = FALSE)
biom |
A |
rank |
The taxonomic rank. E.g. “OTU”, “Phylum”, etc. May also be given numerically: 0 for OTU, 1 for the highest level (i.e. Kingdom), and extending to the number of taxonomic ranks encoded in the original biom file. See example below to fetch the names of all available ranks. |
map |
A character matrix defining the value that each taxa IDs is
assigned for each taxonomic rank. If |
lineage |
Include all ranks in the name of the taxa. For instance,
setting to |
sparse |
If true, returns a sparse matrix as described by
|
A numeric matrix with samples as column names, and taxonomic identifiers as row names.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) colnames(biom$taxonomy) phyla <- taxa.rollup(biom, 'Phylum') phyla[1:4,1:6] # Custom matrices should be formatted like so: counts <- as.matrix(biom$counts) map <- biom$taxonomy counts[1:3,1:6] map[1:3,1:4] phyla <- taxa.rollup(counts, 'Phylum', map=map) phyla[1:3,1:6]
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) colnames(biom$taxonomy) phyla <- taxa.rollup(biom, 'Phylum') phyla[1:4,1:6] # Custom matrices should be formatted like so: counts <- as.matrix(biom$counts) map <- biom$taxonomy counts[1:3,1:6] map[1:3,1:4] phyla <- taxa.rollup(counts, 'Phylum', map=map) phyla[1:3,1:6]
Get the taxonomy table.
taxonomy(biom)
taxonomy(biom)
biom |
A |
A character matrix of the named taxonomies in biom
.
Other accessor functions:
counts()
,
info()
,
metadata()
,
nsamples()
,
ntaxa()
,
phylogeny()
,
sample.names()
,
sequences()
,
taxa.names()
,
taxa.ranks()
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxonomy(biom)[1:4,]
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) taxonomy(biom)[1:4,]
Names of a phylogenetic tree's tips/leafs.
tips(x)
tips(x)
x |
A phylo object, as returned from read.tree.. |
A character vector with the leaf names.
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) leafs <- tips(tree) subtree <- subtree(tree, head(leafs))
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) leafs <- tips(tree) subtree <- subtree(tree, head(leafs))
This is the function called internally by beta.div, but is made
visible here so you can use it with matrices and trees without having to
first convert them to BIOM
objects.
unifrac(biom, weighted = TRUE, tree = NULL)
unifrac(biom, weighted = TRUE, tree = NULL)
biom |
A |
weighted |
Use weighted UniFrac, which takes abundance into account rather than simply presence/absence. |
tree |
A |
A distance matrix of class dist
.
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) biom <- select(biom, 1:10) dm <- unifrac(biom) plot(hclust(dm), cex=.8) as.matrix(dm)[1:4,1:4] # Using a custom matrix and tree mtx <- matrix(sample.int(12*20), ncol=20) dimnames(mtx) <- list(LETTERS[1:12], letters[1:20]) tree <- ape::as.phylo(hclust(dist(mtx))) dm <- unifrac(mtx, tree=tree) as.matrix(dm)[1:4,1:4]
library(rbiom) infile <- system.file("extdata", "hmp50.bz2", package = "rbiom") biom <- read.biom(infile) biom <- select(biom, 1:10) dm <- unifrac(biom) plot(hclust(dm), cex=.8) as.matrix(dm)[1:4,1:4] # Using a custom matrix and tree mtx <- matrix(sample.int(12*20), ncol=20) dimnames(mtx) <- list(LETTERS[1:12], letters[1:20]) tree <- ape::as.phylo(hclust(dist(mtx))) dm <- unifrac(mtx, tree=tree) as.matrix(dm)[1:4,1:4]
Write counts, metadata, taxonomy, and phylogeny to a biom file.
write.biom(biom, file, format = "json")
write.biom(biom, file, format = "json")
biom |
The BIOM object to save to the file. |
file |
Path to the output file. |
format |
Options are “tab”,
“json”, and “hdf5”,
corresponding to classic tabular format, biom format version 1.0 and
biom version 2.1, respectively. Abbreviations are also accepted. See
http://biom-format.org/documentation/ for details. NOTE: to write
HDF5 formatted BIOM files, the BioConductor R package |
On success, returns NULL
invisibly.
Write sequences from a BIOM object to a file in fasta format.
write.fasta(seqs, outfile)
write.fasta(seqs, outfile)
seqs |
A named character vector where names are sequence names and
values are the sequences. Also accepts a |
outfile |
Path to the output fasta file. Files ending in |
On success, returns NULL
invisibly.
Write a newick formatted phylogenetic tree.
write.tree(tree = NULL, file = NULL)
write.tree(tree = NULL, file = NULL)
tree |
A |
file |
Filename or connection to write the newick file to (optional). |
If file is NULL, the newick string as a character vector. Otherwise,
the return value from writeChar
, typically invsible(NULL).
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) newick <- write.tree(tree)
library(rbiom) infile <- system.file("extdata", "newick.tre", package = "rbiom") tree <- read.tree(infile) newick <- write.tree(tree)
Write data and summary information to a Microsoft Excel-compatible workbook.
write.xlsx(biom, outfile, depth = NULL, seed = 0)
write.xlsx(biom, outfile, depth = NULL, seed = 0)
biom |
The BIOM object to save to the file. |
outfile |
Path to the output xlsx file. |
depth |
Depth to rarefy to. See |
seed |
Random seed to use in rarefying. See |
On success, returns NULL
invisibly.
Any data frame attributes on biom
will be included as separate
worksheets. An attribute named 'Reads Per Step' is treated specially and
merged with the usual 'Reads Per Sample' tab - if provided, its row names
should match those in biom
exactly.