anansi wrapper for the MultiAssayExperiment class
Source:R/AllGenerics.R
, R/getAnansi-methods.R
getAnansi.Rd
getAnansi
provides a wrapper to execute the anansi pipeline on a
MultiAssayExperiment::MultiAssayExperiment()
object. It applies the functions weaveWeb()
and
anansi()
in the given order.
Usage
getAnansi(x, ...)
# S4 method for class 'MultiAssayExperiment'
getAnansi(
x,
formula,
link = NULL,
tableY = NULL,
tableX = NULL,
...,
force_new = FALSE
)
Arguments
- x
- ...
additional parameters that can be passed to
AnansiWeb()
oranansi()
.- formula
A formula object. Used to assess differential associations.
- link
One of the following:
Character scalar with value "none"
data.frame with two columns
list with two such data.frames
- tableY, tableX
Character scalar
ornumeric scalar
. Selects experiment corresponding totableY
andtableX
fromexperiments(x)
ofMultiAssayExperiment
object by name or index, name is recommended. (Default slots:Y = 1L
,X = 2L
).- force_new
boolean
If x already has a dictionaryMatrix
in metadata, ignore it and generate a new object anyway? (Default: FALSE).
Value
If return.format
is "table"
(default), a wide format data.frame
intended to be compatible with ggplot2
, or specialized plotting functions
(See plotAnansi()
). If return.format
is
"list"
, a list with aforementioned table, as well as input and
additional information. If return.format
is "raw"
, a list of
raw output (used for testing purposes).
Details
This wrapper of anansi()
allows to perform a complete anansi
analysis directly on objects of class
MultiAssayExperiment::MultiAssayExperiment()
. First, the assays from
experiments specified by tableX
and tableY
are passed to AnansiWeb()
to
build an AnansiWeb object. If there are more than one assay in that
experiment, the specific assay can by specified using assay.type1
and
assay.type2
. Next, this object is fed to the main anansi()
function to
compute interactions between the two assays.
Examples
# Import libraries
library(mia)
#> Loading required package: MultiAssayExperiment
#> Loading required package: SummarizedExperiment
#> Loading required package: MatrixGenerics
#> Loading required package: matrixStats
#>
#> Attaching package: ‘MatrixGenerics’
#> The following objects are masked from ‘package:matrixStats’:
#>
#> colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
#> colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
#> colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
#> colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
#> colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
#> colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
#> colWeightedMeans, colWeightedMedians, colWeightedSds,
#> colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
#> rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
#> rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
#> rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
#> rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
#> rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
#> rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
#> rowWeightedSds, rowWeightedVars
#> Loading required package: GenomicRanges
#> Loading required package: stats4
#> Loading required package: BiocGenerics
#> Loading required package: generics
#>
#> Attaching package: ‘generics’
#> The following objects are masked from ‘package:base’:
#>
#> as.difftime, as.factor, as.ordered, intersect, is.element, setdiff,
#> setequal, union
#>
#> Attaching package: ‘BiocGenerics’
#> The following objects are masked from ‘package:stats’:
#>
#> IQR, mad, sd, var, xtabs
#> The following objects are masked from ‘package:base’:
#>
#> Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
#> as.data.frame, basename, cbind, colnames, dirname, do.call,
#> duplicated, eval, evalq, get, grep, grepl, is.unsorted, lapply,
#> mapply, match, mget, order, paste, pmax, pmax.int, pmin, pmin.int,
#> rank, rbind, rownames, sapply, saveRDS, table, tapply, unique,
#> unsplit, which.max, which.min
#> Loading required package: S4Vectors
#>
#> Attaching package: ‘S4Vectors’
#> The following object is masked from ‘package:tidyr’:
#>
#> expand
#> The following object is masked from ‘package:utils’:
#>
#> findMatches
#> The following objects are masked from ‘package:base’:
#>
#> I, expand.grid, unname
#> Loading required package: IRanges
#> Loading required package: GenomeInfoDb
#> Loading required package: Biobase
#> Welcome to Bioconductor
#>
#> Vignettes contain introductory material; view with
#> 'browseVignettes()'. To cite Bioconductor, see
#> 'citation("Biobase")', and for packages 'citation("pkgname")'.
#>
#> Attaching package: ‘Biobase’
#> The following object is masked from ‘package:MatrixGenerics’:
#>
#> rowMedians
#> The following objects are masked from ‘package:matrixStats’:
#>
#> anyMissing, rowMedians
#> Loading required package: SingleCellExperiment
#> Loading required package: TreeSummarizedExperiment
#> Loading required package: Biostrings
#> Loading required package: XVector
#>
#> Attaching package: ‘Biostrings’
#> The following object is masked from ‘package:base’:
#>
#> strsplit
#> This is mia version 1.17.0
#> - Online documentation and vignettes: https://microbiome.github.io/mia/
#> - Online book 'Orchestrating Microbiome Analysis (OMA)': https://microbiome.github.io/OMA/docs/devel/
library(TreeSummarizedExperiment)
library(MultiAssayExperiment)
web <- randomWeb(n_samples = 100)
mae <- as(web, "MultiAssayExperiment")
# Perform anansi analysis
out <- getAnansi(mae,
tableY = "y", tableX = "x",
formula = ~group_ab
)
#> Fitting least-squares for following model:
#> ~ x + group_ab + x:group_ab
#> Running correlations for the following groups:
#> b, a
# View subset of results
head(out, 5)
#> feature_Y feature_X All_r.values All_t.values All_p.values b_r.values
#> 1 y_1 x_1 0.00258535 0.02559375 0.9796334 -0.05643566
#> 2 y_2 x_1 -0.04556572 0.45154659 0.6525934 0.05769589
#> 3 y_6 x_1 0.11093642 1.10503532 0.2718499 0.04420351
#> 4 y_7 x_1 -0.12850659 1.28278636 0.2025926 0.02209890
#> 5 y_8 x_1 -0.05910079 0.58609243 0.5591609 0.09887132
#> b_t.values b_p.values a_r.values a_t.values a_p.values full_r.squared
#> 1 0.3875210 0.7001199 0.03648797 0.2555860 0.79934039 0.002127277
#> 2 0.3962031 0.6937483 -0.10485754 0.7380716 0.46399127 0.007609320
#> 3 0.3033405 0.7629691 0.15972600 1.1326233 0.26288445 0.014444228
#> 4 0.1515394 0.8801988 -0.27146220 1.9743750 0.05398934 0.034501000
#> 5 0.6811652 0.4991095 -0.20359035 1.4556187 0.15187548 0.022783722
#> full_f.values full_p.values disjointed_group_ab_r.squared
#> 1 0.1023270 0.9028327 0.002115178
#> 2 0.3680479 0.6930563 0.006464888
#> 3 0.7034842 0.4973909 0.002894405
#> 4 1.7152250 0.1853906 0.018168026
#> 5 1.1191163 0.3307907 0.020469724
#> disjointed_group_ab_f.values disjointed_group_ab_p.values
#> 1 0.2034875 0.6529378
#> 2 0.6246676 0.4312661
#> 3 0.2786695 0.5987927
#> 4 1.7764043 0.1857461
#> 5 2.0061590 0.1598982
#> emergent_group_ab_r.squared emergent_group_ab_f.values
#> 1 0.01245822 1.198461
#> 2 0.01245822 1.198461
#> 3 0.01245822 1.198461
#> 4 0.01245822 1.198461
#> 5 0.01245822 1.198461
#> emergent_group_ab_p.values All_q.values b_q.values a_q.values full_q.values
#> 1 0.2763973 0.9807904 0.8401439 0.9841348 0.9708096
#> 2 0.2763973 0.9764323 0.8401439 0.9493973 0.9708096
#> 3 0.2763973 0.9145004 0.8422043 0.9493973 0.9066463
#> 4 0.2763973 0.9145004 0.8801988 0.8638295 0.9066463
#> 5 0.2763973 0.9764323 0.8239828 0.9460500 0.9066463
#> disjointed_group_ab_q.values emergent_group_ab_q.values
#> 1 0.9155455 0.5527945
#> 2 0.9155455 0.5527945
#> 3 0.9155455 0.5527945
#> 4 0.9155455 0.5527945
#> 5 0.9155455 0.5527945