calculateUniFrac {mia} | R Documentation |
This function calculates the (Fast) UniFrac distance for all sample-pairs
in a TreeSummarizedExperiment
object.
## S4 method for signature 'ANY,phylo' calculateUniFrac( x, tree, weighted = FALSE, normalized = TRUE, BPPARAM = SerialParam() ) ## S4 method for signature 'TreeSummarizedExperiment,missing' calculateUniFrac(x, exprs_values = "counts", transposed = FALSE, ...) runUniFrac( x, tree, weighted = FALSE, normalized = TRUE, BPPARAM = SerialParam() )
x |
a numeric matrix or a
Please note that |
tree |
if |
weighted |
|
normalized |
|
BPPARAM |
A
|
exprs_values |
a single |
transposed |
Logical scalar, is x transposed with cells in rows? |
... |
optional arguments not used. |
Please note that if calculateUniFrac
is used as a FUN
for
runMDS
, the argument ntop
has to be set to nrow(x)
.
a sample-by-sample distance matrix, suitable for NMDS, etc.
Paul J. McMurdie. Adapted for mia by Felix G.M. Ernst
http://bmf.colorado.edu/unifrac/
The main implementation (Fast UniFrac) is adapted from the algorithm's description in:
Hamady, Lozupone, and Knight, “Fast UniFrac: facilitating high-throughput phylogenetic analyses of microbial communities including analysis of pyrosequencing and PhyloChip data.” The ISME Journal (2010) 4, 17–27.
See also additional descriptions of UniFrac in the following articles:
Lozupone, Hamady and Knight, “UniFrac - An Online Tool for Comparing Microbial Community Diversity in a Phylogenetic Context.”, BMC Bioinformatics 2006, 7:371
Lozupone, Hamady, Kelley and Knight, “Quantitative and qualitative (beta) diversity measures lead to different insights into factors that structure microbial communities.” Appl Environ Microbiol. 2007
Lozupone C, Knight R. “UniFrac: a new phylogenetic method for comparing microbial communities.” Appl Environ Microbiol. 2005 71 (12):8228-35.
data(esophagus) library(scater) calculateUniFrac(esophagus, weighted = FALSE) calculateUniFrac(esophagus, weighted = TRUE) calculateUniFrac(esophagus, weighted = TRUE, normalized = FALSE) # for using calculateUniFrac in conjunction with runMDS the tree argument # has to be given separately. In addition, subsetting using ntop must # be disabled esophagus <- runMDS(esophagus, FUN = calculateUniFrac, name = "UniFrac", tree = rowTree(esophagus), exprs_values = "counts", ntop = nrow(esophagus)) reducedDim(esophagus)