variance_adjusted_unifrac | R Documentation |
Variance Adjusted UniFrac beta diversity metric.
variance_adjusted_unifrac(counts, tree = NULL, pairs = NULL, cpus = n_cpus())
counts |
An OTU abundance matrix where each column is a sample, and
each row is an OTU. Any object coercible with |
tree |
A |
pairs |
Which combinations of samples should distances be
calculated for? The default value ( |
cpus |
How many parallel processing threads should be used. The
default, |
A dist
object.
Given n
branches with lengths L
and a pair of samples'
abundances (A
and B
) on each of those branches:
D = \displaystyle \frac{\sum_{i = 1}^{n} L_i\displaystyle \frac{|\frac{A_i}{A_T} - \frac{B_i}{B_T}|}{\sqrt{(A_i + B_i)(A_T + B_T - A_i - B_i)}} }{\sum_{i = 1}^{n} L_i\displaystyle \frac{\frac{A_i}{A_T} + \frac{B_i}{B_T}}{\sqrt{(A_i + B_i)(A_T + B_T - A_i - B_i)}} }
See vignette('unifrac')
for details and a worked example.
Chang Q, Luan Y, Sun F 2011. Variance adjusted weighted UniFrac: a powerful beta diversity measure for comparing communities based on phylogeny. BMC Bioinformatics, 12. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.1186/1471-2105-12-118")}
Other beta_diversity:
bray_curtis()
,
canberra()
,
euclidean()
,
generalized_unifrac()
,
gower()
,
jaccard()
,
kulczynski()
,
manhattan()
,
unweighted_unifrac()
,
weighted_normalized_unifrac()
,
weighted_unifrac()
# Example counts matrix
ex_counts
# Variance Adjusted UniFrac distance matrix
variance_adjusted_unifrac(ex_counts, tree = ex_tree)
# Only calculate distances for A vs all.
variance_adjusted_unifrac(ex_counts, tree = ex_tree, pairs = 1:3)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.