| Title: | Annotation-Based Analysis of Specific Interactions |
|---|---|
| Description: | Studies including both microbiome and metabolomics data are becoming more common. Often, it would be helpful to integrate both datasets in order to see if they corroborate each others patterns. All vs all association is imprecise and likely to yield spurious associations. This package takes a knowledge-based approach to constrain association search space, only considering metabolite-function pairs that have been recorded in a pathway database. This package also provides a framework to assess differential association. |
| Authors: | Thomaz Bastiaanssen [aut, cre] (ORCID: <https://orcid.org/0000-0001-6891-734X>), Thomas Quinn [aut] (ORCID: <https://orcid.org/0000-0003-0286-6329>), Giulio Benedetti [aut] (ORCID: <https://orcid.org/0000-0002-8732-7692>), Tuomas Borman [aut] (ORCID: <https://orcid.org/0000-0002-8563-8884>), Leo Lahti [aut] (ORCID: <https://orcid.org/0000-0001-5537-637X>) |
| Maintainer: | Thomaz Bastiaanssen <[email protected]> |
| License: | GPL-3 |
| Version: | 1.3.0 |
| Built: | 2026-05-27 06:17:02 UTC |
| Source: | https://github.com/bioc/anansi |
Coercion functions for anansi
## S3 method for class ''anansi::AnansiWeb'' as.list(x, ...) ## S3 method for class ''anansi::MultiFactor'' as.list(x, ..., use.names = TRUE) ## S3 method for class ''anansi::AnansiWeb'' as.data.frame(x, row.names, optional, ...) asMAE(x) asTSE(x)## S3 method for class ''anansi::AnansiWeb'' as.list(x, ...) ## S3 method for class ''anansi::MultiFactor'' as.list(x, ..., use.names = TRUE) ## S3 method for class ''anansi::AnansiWeb'' as.data.frame(x, row.names, optional, ...) asMAE(x) asTSE(x)
x |
input object |
... |
additional arguments (currently not used). |
use.names |
|
row.names, optional
|
Not used. See ?base::data.frame |
An object of the desired class.
# AnansiWeb x <- randomWeb( n_samples = 5, n_features_x = 4, n_features_y = 6 ) as.list(x) as.data.frame(x) # AnansiWeb to MultiAssayExperiment asMAE(x) asTSE(x) # MultiFactor x <- randomMultiFactor(n_types = 3, n_features = 3) as.list(x, use.names = TRUE)# AnansiWeb x <- randomWeb( n_samples = 5, n_features_x = 4, n_features_y = 6 ) as.list(x) as.data.frame(x) # AnansiWeb to MultiAssayExperiment asMAE(x) asTSE(x) # MultiFactor x <- randomMultiFactor(n_types = 3, n_features = 3) as.list(x, use.names = TRUE)
AnansiTale is the main container that will hold your
stats output data coming out of the anansi pipeline.
AnansiTale( subject = character(0), type = character(0), df = integer(0), estimates = integer(0), f.values = integer(0), t.values = integer(0), p.values = integer(0) )AnansiTale( subject = character(0), type = character(0), df = integer(0), estimates = integer(0), f.values = integer(0), t.values = integer(0), p.values = integer(0) )
subject |
A character that describes the data that was queried. |
type |
A character that describes type of parameter contained in the
|
df |
a vector of length 2, containing df1 and df2 corresponding to the F-ratio considered. |
estimates |
A matrix containing the estimates for the parameters named
in the |
f.values |
A matrix containing the f-values, for least-squares. |
t.values |
A matrix containing the t-values, for correlations. |
p.values |
A matrix containing the p.values for the parameters named in
the |
an AnansiTale
subjectA character that describes the data that was queried.
typeA character that describes type of parameter contained in the
estimates slot. For example r.values for correlations or r.squared
for models.
dfa vector of length 2, containing df1 and df2 corresponding to the F-ratio considered.
estimatesA matrix containing the estimates for the parameters named in
the type slot.
f.valuesA matrix containing the f-values, for least-squares.
t.valuesA matrix containing the t-values, for correlations.
p.valuesA matrix containing the p.values for the parameters named in
the type slot.
AnansiTaleAnansiTale
AnansiWeb is an S7 class containing two feature tables as well as a
dictionary to link them. AnansiWeb is the main container that will
hold your input data throughout the anansi pipeline.
Typical use of the anansi package will involve generating an AnansiWeb
object using the weaveWeb() function.
The function AnansiWeb() constructs an AnansiWeb object from two
feature tables and an adjacency matrix.
## Constructor for `AnansiWeb` objects AnansiWeb(tableX, tableY, dictionary, metadata = data.frame())## Constructor for `AnansiWeb` objects AnansiWeb(tableX, tableY, dictionary, metadata = data.frame())
tableY, tableX
|
A table containing features of interest. Rows should be samples and columns should be features. Y and X refer to the position of the features in a formula: Y ~ X. |
dictionary |
A binary adjacency matrix of class |
metadata |
|
an AnansiWeb object, with sparse binary biadjacency matrix
with features from y as rows and features from x as columns in
dictionary slot.
tableY,tableXTwo matrix objects of measurements, data. Rows are
samples and columns are features. Access with @tableY and @tableX.
dictionaryMatrix, binary adjacency matrix. Optionally sparse.
Typically generated using the weaveWeb() function. Access with
@dictionary.
metadataOptional data.frame of sample metadata. Access with
@metadata.
randomAnansi: For generation of random AnansiWeb objects.
kegg_link(): For examples of input for link argument.
# Use AnansiWeb() to construct an AnansiWeb object from components: tX <- `dimnames<-`(replicate(5, (rnorm(36))), value = list( as.character(seq_len(36)), letters[1:5] ) ) tY <- `dimnames<-`(replicate(3, (rnorm(36))), value = list( as.character(seq_len(36)), LETTERS[1:3] ) ) d <- matrix(TRUE, nrow = NCOL(tY), ncol = NCOL(tX), # Note: Dictionary should have named dimensions dimnames = list( y_names = colnames(tY), x_names = colnames(tX) ) ) web <- AnansiWeb(tableX = tX, tableY = tY, dictionary = d)# Use AnansiWeb() to construct an AnansiWeb object from components: tX <- `dimnames<-`(replicate(5, (rnorm(36))), value = list( as.character(seq_len(36)), letters[1:5] ) ) tY <- `dimnames<-`(replicate(3, (rnorm(36))), value = list( as.character(seq_len(36)), LETTERS[1:3] ) ) d <- matrix(TRUE, nrow = NCOL(tY), ncol = NCOL(tX), # Note: Dictionary should have named dimensions dimnames = list( y_names = colnames(tY), x_names = colnames(tX) ) ) web <- AnansiWeb(tableX = tX, tableY = tY, dictionary = d)
Methods for AnansiWeb S7 container class
x |
input, AnansiWeb object |
The desired information from an AnansiWeb object
weaveWeb(): for general use.
AnansiWeb-pairwise: for methods for pairwise operations
# Setup web <- randomWeb(n_samp = 36) # Accessors dimnames(web) dim(web) names(web) # Getters and setters: tableX(web)[1:5, 1:5] tableY(web)[1:5, 1:5] dictionary(web) head(metadata(web)) # Assign some random metadata metadata(web) <- data.frame( id = row.names(tableY(web)), a = rnorm(36), b = sample(c("a", "b"), 36, TRUE), row.names = "id" ) # Coerce to list weblist <- as.list(web) # Coerce to Data.frame webdf <- as.data.frame(web) # Coerce to MultiAssayExperiment mae <- asMAE(web) # Coerce to TreeSummarizedExperiment tse <- asTSE(web)# Setup web <- randomWeb(n_samp = 36) # Accessors dimnames(web) dim(web) names(web) # Getters and setters: tableX(web)[1:5, 1:5] tableY(web)[1:5, 1:5] dictionary(web) head(metadata(web)) # Assign some random metadata metadata(web) <- data.frame( id = row.names(tableY(web)), a = rnorm(36), b = sample(c("a", "b"), 36, TRUE), row.names = "id" ) # Coerce to list weblist <- as.list(web) # Coerce to Data.frame webdf <- as.data.frame(web) # Coerce to MultiAssayExperiment mae <- asMAE(web) # Coerce to TreeSummarizedExperiment tse <- asTSE(web)
Methods to run pairwise analysis on multi-modal data.
X, x
|
input, AnansiWeb object |
... |
additional arguments |
FUN |
a function with at least two arguments. The variables |
MoreArgs, SIMPLIFY, USE.NAMES
|
see ?base::mapply |
which |
|
with.metadata |
|
pairs: A two-column array index, corresponding to i,j coordinates in matrix
notation.
getFeaturePairs: A list of data.frames with the paired data
web <- randomWeb() # Extract data.frames in pairs (only show first) getFeaturePairs(web)[1L] pairs(x = web) getFeaturePairs(x = web, which = NULL, with.metadata = FALSE) pairwiseApply( X = web, FUN = function(x, y) cor(x, y), MoreArgs = NULL, SIMPLIFY = TRUE, USE.NAMES = TRUE )web <- randomWeb() # Extract data.frames in pairs (only show first) getFeaturePairs(web)[1L] pairs(x = web) getFeaturePairs(x = web, which = NULL, with.metadata = FALSE) pairwiseApply( X = web, FUN = function(x, y) cor(x, y), MoreArgs = NULL, SIMPLIFY = TRUE, USE.NAMES = TRUE )
Get dictionary
dictionary(x, ...)dictionary(x, ...)
x |
input object |
... |
additional arguments |
dictionary slot of x
x <- randomWeb(10) dictionary(x)x <- randomWeb(10) dictionary(x)
kegg_link() is a convenience function to return a list containing two
data.frames; ec2cpd and ec2ko. This will be their most likely use.
ec2cpd and ec2ko are two data.frames, used to link ko, ecs and cpd
identifiers in the KEGG database.
data("ec2ko", package = "anansi") data("ec2cpd", package = "anansi") kegg_link()data("ec2ko", package = "anansi") data("ec2cpd", package = "anansi") kegg_link()
ec2ko: a data.frame of two columns, named "ec"
and "ko". The IDs refer to KEGG orthologues. Enzyme commission numbers,
ecs, typically describe reactions captured by them.
ec2cpd: a data.frame of two columns, named "ec" and "cpd".
The IDs refer to compounds in the KEGG database. Enzyme commission
numbers, ecs, typically describe reactions either producing or requiring
them.
kegg_link() returns a list containing the two aforementioned
data.frames, ec2cpd and ec2ko.
ec2ko: Adapted from https://www.genome.jp/kegg/, using
KEGGREST. Script to generate available in example.
ec2cpd: Adapted from https://www.genome.jp/kegg/ using
KEGGREST. Script to generate available in example.
kegg_link() # Generate ec2ko and ec2cpd: # Don't download during tests. set to `TRUE` to download. dry_run <- TRUE if (!dry_run) { ec2ko <- KEGGREST::keggLink("ec", "ko") ec2ko <- data.frame( ec = gsub("ec:", "", x = ec2ko, fixed = TRUE), ko = gsub("ko:", "", x = names(ec2ko), fixed = TRUE), row.names = NULL ) ec2cpd <- KEGGREST::keggLink("ec", "cpd") ec2cpd <- data.frame( ec = gsub("ec:", "", x = ec2cpd, fixed = TRUE), cpd = gsub("cpd:", "", x = names(ec2cpd), fixed = TRUE), row.names = NULL ) }kegg_link() # Generate ec2ko and ec2cpd: # Don't download during tests. set to `TRUE` to download. dry_run <- TRUE if (!dry_run) { ec2ko <- KEGGREST::keggLink("ec", "ko") ec2ko <- data.frame( ec = gsub("ec:", "", x = ec2ko, fixed = TRUE), ko = gsub("ko:", "", x = names(ec2ko), fixed = TRUE), row.names = NULL ) ec2cpd <- KEGGREST::keggLink("ec", "cpd") ec2cpd <- data.frame( ec = gsub("ec:", "", x = ec2cpd, fixed = TRUE), cpd = gsub("cpd:", "", x = names(ec2cpd), fixed = TRUE), row.names = NULL ) }
Piphillin was used to infer functions from the 16S sequencing data in terms of KOs. Unfortunately, the Piphillin algorithm is proprietary and has since been taken down.
data("FMT_data", package = "anansi")data("FMT_data", package = "anansi")
A marix object with 6468 rows, KOs, and 36 columns, samples.
doi:10.1038/s43587-021-00093-9
Snippet of the CLR-transformed hippocampal metabolomics data from the FMT Aging study.
data("FMT_data", package = "anansi")data("FMT_data", package = "anansi")
A matrix object with three rows, compounds, and 36 columns, samples.
doi:10.1038/s43587-021-00093-9
There were three treatment groups in the study that all received faecal microbiota transplantation (FMT). Young mice that received FMT from young mice (Young yFMT), aged mice that received FMT from aged mice (Aged oFMT) and aged mice that received FMT from young mice (Aged yFMT).
data("FMT_data", package = "anansi")data("FMT_data", package = "anansi")
A data.frame object with 36 rows, samples, and two columns, denoting sample ID and treatment group, respectively.
doi:10.1038/s43587-021-00093-9
Get a listof edges
getEdgeList(x, ...)getEdgeList(x, ...)
x |
input object |
... |
additional arguments |
a two-column data.frame that lists the content of each entry in the input MultiFactor
x <- randomMultiFactor(n_features = 10) getEdgeList(x)x <- randomMultiFactor(n_features = 10) getEdgeList(x)
Get a list of all pairs of features
getFeaturePairs(x, ...)getFeaturePairs(x, ...)
x |
input object |
... |
additional arguments for specific methods |
an list of two-column data.frames that represent all feature pairs.
x <- randomWeb(10) head(getFeaturePairs(x), 3)x <- randomWeb(10) head(getFeaturePairs(x), 3)
Get a graph object.
getGraph(x, ...)getGraph(x, ...)
x |
input |
... |
additional arguments (currently not used). |
a specified graph object.
# Show methods getGraph# Show methods getGraph
Enzymes and metabolites that make up the Krebs cycle. Features share a row when they play a role in the same reaction. Used in vignettes, exported to facilitate user exploration.
data("krebs", package = "anansi")data("krebs", package = "anansi")
A data.frame of two factors, named "Enzyme" and "Metabolite",
which contain these respective components of the krebs cycle.
Generated by hand based on the krebs cycle.
make a link data.frame for biobakery mapping files input
linkBiobakeryMap(map)linkBiobakeryMap(map)
map |
|
a two-column data.frame that can be converted into an adjacency
matrix used as input for weaveWeb().
# some dummy input, as a character vector of IDs separated by tabs. x <- c("x_1\ty_1\ty_2\ty_4", "x_2\ty_1\ty_3", "x_3\ty_2\ty_y") linkBiobakeryMap(x)# some dummy input, as a character vector of IDs separated by tabs. x <- c("x_1\ty_1\ty_2\ty_4", "x_2\ty_1\ty_3", "x_3\ty_2\ty_y") linkBiobakeryMap(x)
Get metadata.
x |
input object |
... |
additional arguments |
Compatible with S4Vectors generic.
metadata slot of x
x <- randomWeb(10) metadata(x)x <- randomWeb(10) metadata(x)
Set metadata.
x |
input object |
... |
additional arguments |
value |
replacement value, coerced to data.frame |
Compatible with S4Vectors generic.
x with modified metadata slot.
x <- randomWeb(10) metadata(x) <- cbind( metadata(x), new_groups = c("A", "B") )x <- randomWeb(10) metadata(x) <- cbind( metadata(x), new_groups = c("A", "B") )
MultiFactor is an S7 class to organize and manage multiple sets of factors,
for instance when tracing or converting feature IDs across databases. Methods
for MultiFactor aim to follow factor behaviour.
MultiFactor(x, levels = NULL, drop.unmatched = FALSE) ## Constructor for `MultiFactor` objects MultiFactor(x, levels = NULL, drop.unmatched = FALSE)MultiFactor(x, levels = NULL, drop.unmatched = FALSE) ## Constructor for `MultiFactor` objects MultiFactor(x, levels = NULL, drop.unmatched = FALSE)
x |
|
levels |
an optional named list of vectors of the unique values (as character strings) that x might have taken. The default is the unique set of values taken by lapply(x, as.character), sorted into increasing order of x. |
drop.unmatched |
|
The most straightforward way to construct a MultiFactor object is as a
named list of named data.frames. The columns of the data.frames indicate the
category of factor in that column.
A MultiFactor object presents itself similar to a data.frame, in the
sense that level types can be called as columns and individual data.frame
components can be called as rows.
a MultiFactor object.
indexNamed list of named integer data frames of at least two columns
each. The column names correspond to names in the levels slot. Similar
to factors, the integers in those columns correspond to the characters
in that level. Accessed through regular list methods (e.g., [, [[).
levelsNamed list of character vectors. Accessed through levels(x)
map(sparse) Matrix specifying which elements contain which levels.
Accesses through x@dictionary.
kegg_link(): for an example of valid input.
# Generate some random linkage input x <- data.frame( a = sample(letters[seq(3)], 10, replace = TRUE), A = sample(LETTERS[seq(3)], 10, replace = TRUE) ) MultiFactor(x)# Generate some random linkage input x <- data.frame( a = sample(letters[seq(3)], 10, replace = TRUE), A = sample(LETTERS[seq(3)], 10, replace = TRUE) ) MultiFactor(x)
droplevels() droplevels(MultiFactor) returns a MultiFactor
with unused levels removed, Analogous to the factor method.
... |
|
value |
Replacement value, typically of same type as that which is to be replaced. |
exclude |
|
select |
|
use.names, ignore.mcols
|
For compatibility, not used. |
x |
|
format |
|
Only one of select and exclude should be provided, as they are
each others complement.
A MultiFactor
a specified graph object.
igraph::graph_from_data_frame() and
igraph::as_graphnel(), which are used under the hood, from
igraph::igraph() package.
# Setup x <- MultiFactor(kegg_link()) x # Basic properties dim(x) dimnames(x) # Factor-like properties head(levels(x)$ko) droplevels(x) head(unfactor(x)$ec2ko) # Extract common output formats getEdgeList(x) droplevels(x, exclude = list(ko = "K00001")) droplevels(x, select = list(ko = "K00001")) # Generate an igraph object g <- getGraph(x = kegg_link(), format = "igraph") plot(g)# Setup x <- MultiFactor(kegg_link()) x # Basic properties dim(x) dimnames(x) # Factor-like properties head(levels(x)$ko) droplevels(x) head(unfactor(x)$ec2ko) # Extract common output formats getEdgeList(x) droplevels(x, exclude = list(ko = "K00001")) droplevels(x, select = list(ko = "K00001")) # Generate an igraph object g <- getGraph(x = kegg_link(), format = "igraph") plot(g)
Apply a function on each pair of features
pairwiseApply(X, ...)pairwiseApply(X, ...)
X |
input object |
... |
additional arguments |
a list containing the output of applying the function to each feature pair.
See ?base::mapply()
web <- randomWeb(10) # For each feature pair, was the value for x higher than the value for y? pairwise_gt <- pairwiseApply( X = web, FUN = function(x, y) x > y, MoreArgs = NULL, SIMPLIFY = FALSE, USE.NAMES = TRUE ) head(pairwise_gt) # Run cor.test() on each pair of features pairwise_cor <- pairwiseApply( X = web, FUN = function(x, y) cor.test(x, y), MoreArgs = NULL, SIMPLIFY = FALSE, USE.NAMES = TRUE ) pairwise_cor[1]web <- randomWeb(10) # For each feature pair, was the value for x higher than the value for y? pairwise_gt <- pairwiseApply( X = web, FUN = function(x, y) x > y, MoreArgs = NULL, SIMPLIFY = FALSE, USE.NAMES = TRUE ) head(pairwise_gt) # Run cor.test() on each pair of features pairwise_cor <- pairwiseApply( X = web, FUN = function(x, y) cor.test(x, y), MoreArgs = NULL, SIMPLIFY = FALSE, USE.NAMES = TRUE ) pairwise_cor[1]
plotAnansi generates an association plot from the output of
anansi() in the table format. It provides a convenient way to
visually assess relevant results from the anansi analysis, either in the
form of a dotplot or a graph.
x |
a |
... |
additional parameters |
layout |
|
association.type |
|
model.var |
|
group |
|
signif.threshold |
|
colour_by |
|
color_by |
|
fill_by |
|
size_by |
|
shape_by |
|
x_lab |
|
y_lab |
|
y_position |
|
show.cor |
|
plotAnansi provides a standardised method to visualise the results
of anansi by means of a differential association plot. The input for this
function should be generated from anansi() or
anansi(), with return.format = "table"
a figure that can be further modified using the ggplot2 suite
A ggplot2 object.
# Import libraries library(mia) library(TreeSummarizedExperiment) library(MultiAssayExperiment) library(ggraph) web <- randomWeb(n_samples = 100) mae <- asMAE(web) # Perform anansi analysis out <- weaveWeb(mae, tableY = "y", tableX = "x" ) |> anansi(formula = ~group_ab) # Select significant interactions out <- out[out$full_p.values < 0.05, ] # Visualise disjointed associations filled by group plotAnansi(out, association.type = "disjointed", model.var = "group_ab", signif.threshold = 0.05, fill_by = "group" ) # Visualise full associations filled by group plotAnansi(out, association.type = "full", signif.threshold = 0.05, fill_by = "group" ) # Visualise full associations as graph plotAnansi(out, layout = "graph", association.type = "full", signif.threshold = 0.05, show.cor = TRUE ) # Visualise disjointed associations as graph plotAnansi(out, layout = "graph", association.type = "disjointed", model.var = "group_ab", signif.threshold = 0.05 )# Import libraries library(mia) library(TreeSummarizedExperiment) library(MultiAssayExperiment) library(ggraph) web <- randomWeb(n_samples = 100) mae <- asMAE(web) # Perform anansi analysis out <- weaveWeb(mae, tableY = "y", tableX = "x" ) |> anansi(formula = ~group_ab) # Select significant interactions out <- out[out$full_p.values < 0.05, ] # Visualise disjointed associations filled by group plotAnansi(out, association.type = "disjointed", model.var = "group_ab", signif.threshold = 0.05, fill_by = "group" ) # Visualise full associations filled by group plotAnansi(out, association.type = "full", signif.threshold = 0.05, fill_by = "group" ) # Visualise full associations as graph plotAnansi(out, layout = "graph", association.type = "full", signif.threshold = 0.05, show.cor = TRUE ) # Visualise disjointed associations as graph plotAnansi(out, layout = "graph", association.type = "disjointed", model.var = "group_ab", signif.threshold = 0.05 )
Randomly generate a valid AnansiWeb or MultiFactor object.
randomWeb( n_samples = 10, n_reps = 1L, n_features_x = 8, n_features_y = 12, sparseness = 0.5, tableY = NULL, tableX = NULL, dictionary = NULL ) randomMultiFactor(n_types = 6, n_features = 100, sparseness = 0.5) krebsDemoWeb(n_samples = 100, n_reps = 4L)randomWeb( n_samples = 10, n_reps = 1L, n_features_x = 8, n_features_y = 12, sparseness = 0.5, tableY = NULL, tableX = NULL, dictionary = NULL ) randomMultiFactor(n_types = 6, n_features = 100, sparseness = 0.5) krebsDemoWeb(n_samples = 100, n_reps = 4L)
n_samples, n_reps
|
|
n_features_y, n_features_x
|
|
sparseness |
|
tableY, tableX
|
A table containing features of interest. Rows should be samples and columns should be features. Y and X refer to the position of the features in a formula: Y ~ X. |
dictionary |
A binary adjacency matrix of class |
n_types |
|
n_features |
|
a randomly generated object of the specified class.
# Make a random AnansiWeb object randomWeb() krebsDemoWeb() randomMultiFactor()# Make a random AnansiWeb object randomWeb() krebsDemoWeb() randomMultiFactor()
Get tableX
tableX(x, ...)tableX(x, ...)
x |
input object |
... |
additional arguments |
tableX slot of x
x <- randomWeb(10) tableX(x)x <- randomWeb(10) tableX(x)
Get tableY
tableY(x, ...)tableY(x, ...)
x |
input object |
... |
additional arguments |
tableY slot of x
x <- randomWeb(10) tableY(x)x <- randomWeb(10) tableY(x)
Weave an AnansiWeb object
weaveWeb(x, ...) weaveKEGG(x, ...)weaveWeb(x, ...) weaveKEGG(x, ...)
x |
input object |
... |
additional arguments |
an AnansiWeb object, with sparse binary biadjacency matrix
with features from y as rows and features from x as columns in
dictionary slot.
# Setup demo tables ec2ko <- kegg_link()[["ec2ko"]] ec2cpd <- kegg_link()[["ec2cpd"]] # Basic usage weaveWeb(cpd ~ ko, link = kegg_link()) weaveWeb(x = "ko", y = "ec", link = ec2ko) weaveWeb(ec ~ cpd, link = ec2cpd) # A wrapper is available for kegg ko, ec and cpd data generic <- weaveWeb(cpd ~ ko, link = kegg_link()) kegg_wrapper <- weaveKEGG(cpd ~ ko) identical(generic, kegg_wrapper) # The following are equivalent to transposition: a <- weaveWeb(ko ~ cpd, link = kegg_link()) |> dictionary() b <- weaveWeb(cpd ~ ko, link = kegg_link()) |> dictionary() identical(a, Matrix::t(b))# Setup demo tables ec2ko <- kegg_link()[["ec2ko"]] ec2cpd <- kegg_link()[["ec2cpd"]] # Basic usage weaveWeb(cpd ~ ko, link = kegg_link()) weaveWeb(x = "ko", y = "ec", link = ec2ko) weaveWeb(ec ~ cpd, link = ec2cpd) # A wrapper is available for kegg ko, ec and cpd data generic <- weaveWeb(cpd ~ ko, link = kegg_link()) kegg_wrapper <- weaveKEGG(cpd ~ ko) identical(generic, kegg_wrapper) # The following are equivalent to transposition: a <- weaveWeb(ko ~ cpd, link = kegg_link()) |> dictionary() b <- weaveWeb(cpd ~ ko, link = kegg_link()) |> dictionary() identical(a, Matrix::t(b))
Generate a biadjacency matrix, linking the features between two tables.
Return an AnansiWeb object which contains all three.
weaveWeb() is for general use and has flexible default settings.
weaveKEGG() is a wrapper that sets link to kegg_link().
All variants are special cases of weaveWeb().
x, y
|
|
link |
One of the following:
|
tableY, tableX
|
A table containing features of interest. Rows should be
samples and columns should be features. Y and X refer to the position of
the features in a formula: Y ~ X. |
metadata |
Optional |
verbose |
|
typeY, typeX
|
|
experiment1, experiment2
|
synonymous args to |
assay.type1, assay.type2
|
synonymous args to |
If the link argument is "none", all features will be considered
linked. If one or more data.frames, colnames should be as specified in
x and y.
an AnansiWeb object, with sparse binary biadjacency matrix
with features from y as rows and features from x as columns in
dictionary slot.
AnansiWeb: For general constructor and methods.
kegg_link(): For examples of input for link argument.
# Setup demo tables, see first vignette. data(FMT_data) t1 <- t(FMT_metab) t2 <- t(FMT_KOs) # Input objects and syntax: ## define `x` and `y` as characters web <- weaveWeb( x = "ko", y = "cpd", link = kegg_link(), tableX = t2, tableY = t1, metadata = NULL, verbose = TRUE ) ## define `x` and `y` with a formula web2 <- weaveWeb( x = cpd ~ ko, link = kegg_link(), tableX = t2, tableY = t1, metadata = NULL, verbose = TRUE ) identical(web, web2) # Method for MultiAssayExperiment S4 object mae <- asMAE(web) weaveWeb( x = mae, link = kegg_link(), tableY = "cpd", tableX = "ko", force_new = FALSE ) # Method for TreeSummarizedExperiment S4 object tse <- asTSE(web) weaveWeb( x = tse, link = kegg_link(), tableY = "cpd", tableX = "ko", force_new = FALSE )# Setup demo tables, see first vignette. data(FMT_data) t1 <- t(FMT_metab) t2 <- t(FMT_KOs) # Input objects and syntax: ## define `x` and `y` as characters web <- weaveWeb( x = "ko", y = "cpd", link = kegg_link(), tableX = t2, tableY = t1, metadata = NULL, verbose = TRUE ) ## define `x` and `y` with a formula web2 <- weaveWeb( x = cpd ~ ko, link = kegg_link(), tableX = t2, tableY = t1, metadata = NULL, verbose = TRUE ) identical(web, web2) # Method for MultiAssayExperiment S4 object mae <- asMAE(web) weaveWeb( x = mae, link = kegg_link(), tableY = "cpd", tableX = "ko", force_new = FALSE ) # Method for TreeSummarizedExperiment S4 object tse <- asTSE(web) weaveWeb( x = tse, link = kegg_link(), tableY = "cpd", tableX = "ko", force_new = FALSE )