Title: | Integral Analysis of Diversity Based on Hill Numbers |
---|---|
Description: | Tools for analysing, comparing, visualising and partitioning diversity based on Hill numbers. 'hilldiv' is an R package that provides a set of functions to assist analysis of diversity for diet reconstruction, microbial community profiling or more general ecosystem characterisation analyses based on Hill numbers, using OTU/ASV tables and associated phylogenetic trees as inputs. The package includes functions for (phylo)diversity measurement, (phylo)diversity profile plotting, (phylo)diversity comparison between samples and groups, (phylo)diversity partitioning and (dis)similarity measurement. All of these grounded in abundance-based and incidence-based Hill numbers. The statistical framework developed around Hill numbers encompasses many of the most broadly employed diversity (e.g. richness, Shannon index, Simpson index), phylogenetic diversity (e.g. Faith's PD, Allen's H, Rao's quadratic entropy) and dissimilarity (e.g. Sorensen index, Unifrac distances) metrics. This enables the most common analyses of diversity to be performed while grounded in a single statistical framework. The methods are described in Jost et al. (2007) <DOI:10.1890/06-1736.1>, Chao et al. (2010) <DOI:10.1098/rstb.2010.0272> and Chiu et al. (2014) <DOI:10.1890/12-0960.1>; and reviewed in the framework of molecularly characterised biological systems in Alberdi & Gilbert (2019) <DOI:10.1111/1755-0998.13014>. |
Authors: | Antton Alberdi [aut, cre] |
Maintainer: | Antton Alberdi <[email protected]> |
License: | GPL-3 |
Version: | 1.5.3 |
Built: | 2024-11-16 06:22:43 UTC |
Source: | https://github.com/anttonalberdi/hilldiv |
Compute alpha diversity of a system comprised of multiple samples from a count (OTU/ASV) table. If a tree object is provided, the computed alpha diversity accounts for the phylogenetic relations across OTUs/ASVs.
alpha_div(countable,qvalue,tree,weight)
alpha_div(countable,qvalue,tree,weight)
countable |
A count table (matrix/data.frame) indicating the absolute or relative OTU/ASV abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
qvalue |
A positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the count table. Use the function match_data() if the count table and tree names do not match. |
weight |
A vector indicating the relative weight of each sample. The order needs to be identical to the order of the samples in the OTU table. The values need to sum up to 1. If empty, all samples are weighed the same. |
Alpha diversity computation (based on Hill numbers)
An alpha diversity value.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
div_part
, gamma_div
, match_data
data(bat.diet.otutable) data(bat.diet.tree) alpha_div(countable=bat.diet.otutable,qvalue=1) alpha_div(countable=bat.diet.otutable,qvalue=1,tree=bat.diet.tree) weight.vector = rep(1/ncol(bat.diet.otutable),ncol(bat.diet.otutable)) alpha_div(bat.diet.otutable,1,bat.diet.tree,weight.vector)
data(bat.diet.otutable) data(bat.diet.tree) alpha_div(countable=bat.diet.otutable,qvalue=1) alpha_div(countable=bat.diet.otutable,qvalue=1,tree=bat.diet.tree) weight.vector = rep(1/ncol(bat.diet.otutable),ncol(bat.diet.otutable)) alpha_div(bat.diet.otutable,1,bat.diet.tree,weight.vector)
Hierarchy table indicating the relationship between samples and their respective parent groups.
bat.diet.hierarchy
bat.diet.hierarchy
A data frame with 40 rows and 2 columns.
An OTU table containing the absolute read abundances of 363 OTUs in 40 faecal samples from 8 different bat species.
bat.diet.otutable
bat.diet.otutable
A data frame with 363 rows and 40 species.
Phylogenetic tree built from the representative sequences of the 363 OTUs included in the 'bat.diet.otutable' data set.
bat.diet.tree
bat.diet.tree
A phylo object with 363 tips and 362 internal nodes.
Compute dissimilarity or similarity values based on beta diversities (neutral or phylogenetic) and sample size.
beta_dis(beta, qvalue, N, metric, type)
beta_dis(beta, qvalue, N, metric, type)
beta |
A numeric beta diversity value or an object outputted by function div_part() (which contains all the information to compute (dis)similarities). |
qvalue |
A positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
N |
An integer indicating sample size, the number of sampling units to be used to compute the (dis)similarity measures. The argument is ovewritten if a 'div_part' object is used. |
metric |
A vector containing "C", "U", "V" or "S". C: Sørensen-type overlap or complement. U: Jaccard-type overlap or complement. V: Sørensen-type turnover or complement. S: Jaccard-type turnover or complement. See hilldiv wiki for further information. |
type |
A character object containing either "similarity" or "dissimilarity". If 'similarity' is used, similarity metrics (0: completely different composition - 1: identical composition) are returned. If 'dissimilarity' is used, dissimilarity metrics (0: identical composition - 1:completely different composition) are returned. |
(Dis)similarity computation from beta diversities based on Hill numbers
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
data(bat.diet.otutable) data(bat.diet.tree) #Manually indicating beta diversity, order of diversity and sample size beta_dis(beta=4.5,qvalue=1,N=8) beta_dis(beta=4.5,qvalue=1,N=8,metric="C",type="similarity") #Using an object created with the function div_part() divpartobject <- div_part(bat.diet.otutable,qvalue=0,tree=bat.diet.tree) beta_dis(divpartobject) beta_dis(divpartobject,metric="S",type="similarity")
data(bat.diet.otutable) data(bat.diet.tree) #Manually indicating beta diversity, order of diversity and sample size beta_dis(beta=4.5,qvalue=1,N=8) beta_dis(beta=4.5,qvalue=1,N=8,metric="C",type="similarity") #Using an object created with the function div_part() divpartobject <- div_part(bat.diet.otutable,qvalue=0,tree=bat.diet.tree) beta_dis(divpartobject) beta_dis(divpartobject,metric="S",type="similarity")
As DNA sequencing data include PCR and sequencing errors, copy number thresholds are commonly applied to discard the OTUs with low number of sequence copies. This threshold can be absolute or (ideally) relative to the sequencing depth of each sample.
copy_filt(abund, threshold, filter)
copy_filt(abund, threshold, filter)
abund |
A vector or a matrix/data.frame indicating the relative abundances of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
threshold |
An integer or a decimal number indicating the cut-off threshold. If an integer is provided, an absolute threshold is used (same threshold for all samples). If a decimal number is provided a relative copy number threshold is applied (dependent on the sequencing depth of each sample). |
filter |
Whether to remove the OTUs/ASVs with no read counts. Default=TRUE. |
OTU/ASV copy number filtering
Antton Alberdi, [email protected]
Alberdi A, Aizpurua O, Bohmann K, Gopalakrishnan S, Lynggaard C, Nielsen M, Gilbert MTP. 2019. Promises and pitfalls of using high-throughput sequencing for diet analysis. Molecular Ecology Resources, 19(2), 327-348.
data(bat.diet.otutable) #Remove singletons from all samples copy_filt(bat.diet.otutable,2) #Remove OTUs represented by less than 0.01% of the total reads per sample. copy_filt(bat.diet.otutable,0.0001)
data(bat.diet.otutable) #Remove singletons from all samples copy_filt(bat.diet.otutable,2) #Remove OTUs represented by less than 0.01% of the total reads per sample. copy_filt(bat.diet.otutable,0.0001)
The Sørensen-type overlap quantifies the effective average proportion of a sub-systems OTUs (or lineages in the case of phylodiversities) that is shared across all subsystems. This is thus a metric that quantifies overlap from the subsystems perspective. Its corresponding dissimilarity measure (1 - CqN) quantifies the effective average proportion of nonshared OTUs or lineages in a system. CqN is integrated in the functions beta_dis() and pair_dis().
CqN(beta, qvalue, N)
CqN(beta, qvalue, N)
beta |
A beta diversity value based on Hill numbers. |
qvalue |
The q value used to compute the beta diversity. It needs to be a positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
N |
An integer indicating sample size, the number of sampling units to be used to compute the similarity measure. |
Sørensen-type overlap
A Sørensen-type overlap value
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
CqN(beta=1.24,qvalue=1,N=3) CqN(1.24,1,3)
CqN(beta=1.24,qvalue=1,N=3) CqN(1.24,1,3)
Coverage of the estimated Hill numbers at different orders of diversity.
depth_cov(abund, qvalue)
depth_cov(abund, qvalue)
abund |
A vector or a matrix/data.frame indicating the relative abundances of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
qvalue |
A positive integer or decimal number (>=0), usually between 0 and 3. |
Depth coverage assessment
A matrix with observed diversity, estimated diversities and coverage
Antton Alberdi, [email protected]
Chao, A. & Jost, L. (2015) Estimating diversity and entropy profiles via discovery rates of new species. Methods in Ecology and Evolution, 6, 873-882.
Jost, L. (2006). Entropy and diversity. Oikos, 113, 363-375.
Hill, M. O. (1973). Diversity and evenness: a unifying notation and its consequences. Ecology, 54, 427-432.
data(bat.diet.otutable) depth_cov(bat.diet.otutable,0) depth_cov(bat.diet.otutable,qvalue=1)
data(bat.diet.otutable) depth_cov(bat.diet.otutable,0) depth_cov(bat.diet.otutable,qvalue=1)
Filter samples based on a minimum sequencing depth.
depth_filt(countable, threshold)
depth_filt(countable, threshold)
countable |
An OTU table (matrix/data.frame) indicating the absolute OTU abundances of multiple samples. Columns must refer to samples and rows to OTUs. |
threshold |
A number indicating the minimum sequencing depth required to keep the sample. |
Sequencing depth filtering
Antton Alberdi, [email protected]
Alberdi A, Aizpurua O, Bohmann K, Gopalakrishnan S, Lynggaard C, Nielsen M, Gilbert MTP. 2019. Promises and pitfalls of using high-throughput sequencing for diet analysis. Molecular Ecology Resources, 19(2), 327-348.
data(bat.diet.otutable) depth_filt(bat.diet.otutable,5000) depth_filt(bat.diet.otutable,threshold=20000)
data(bat.diet.otutable) depth_filt(bat.diet.otutable,5000) depth_filt(bat.diet.otutable,threshold=20000)
Visualisation of pairwise dissimilarities
dis_nmds(distance, hierarchy, colour, plot, centroids, labels, legend, runs)
dis_nmds(distance, hierarchy, colour, plot, centroids, labels, legend, runs)
distance |
Matrix of pairwise dissimilarities, usually one of the matrices listed in the output object of the pair_dis() function. |
hierarchy |
The first column containing the sample names while the second containing the groups names. If provided, dots are coloured according to groups and group centroids can also be visualised. |
colour |
The number of vector items (colours, e.g. '#34k235'), must be of length one if no hierarchy table is added, or must equal the number of groups if the hierarchy table is provided. |
plot |
Whether to plot a NMDS or a Shepard plot. Default: "NMDS". |
centroids |
Whether to link sample dots with group centroids or not. A hierarchy table is necessary to draw centroids. Default: FALSE |
labels |
Whether to print sample or group labels or both. A hierarchy table is necessary to plot grpup names. Default: "none". |
legend |
Whether to print the legend or not. Default: TRUE. |
runs |
Number of iterations for the NMDS function. Default: 100. |
Dissimilarity NMDS plot
An NMDS or Shepard plot.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pairdisres <- pair_dis(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) dis_nmds(pairdisres$L1_CqN) dis_nmds(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy, centroids=TRUE) dis_nmds(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy, centroids=TRUE, labels="group")
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pairdisres <- pair_dis(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) dis_nmds(pairdisres$L1_CqN) dis_nmds(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy, centroids=TRUE) dis_nmds(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy, centroids=TRUE, labels="group")
Multi-level diversity partitioning following the multiplicative definition based on Hill numbers. Hierarchical levels are defined from L1 (minimum, sample) to Ln (maximum, whole system), and as many intermediate levels as wanted can be defined in between. The hierarchical structure of the system is defined with the hierarchy table. If no hierarchy table is inputed, the function yields a simple two-level partitioning between alpha (L1), beta and gamma (L2).
div_part(countable, qvalue, tree, hierarchy)
div_part(countable, qvalue, tree, hierarchy)
countable |
A count table (matrix/data.frame) indicating the absolute or relative OTU/ASV abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
qvalue |
A positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the OTU table. Use the function match_data() if the OTU names do not match. |
hierarchy |
A matrix indicating the relation between samples (first column) and parent group(s). |
Multi-level diversity partitioning (based on Hill numbers)
A list object containing details of hierarchical diversity partitioning.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427–2439.
div_part
, gamma_div
, match_data
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #Two level examples (L1=sample (alpha diversity), L2=whole system (gamma diversity)) div_part(bat.diet.otutable,qvalue=1) div_part(bat.diet.otutable,qvalue=0,tree=bat.diet.tree) #Three-level example (L1=sample, L2=species, L3=whole system) div_part(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy)
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #Two level examples (L1=sample (alpha diversity), L2=whole system (gamma diversity)) div_part(bat.diet.otutable,qvalue=1) div_part(bat.diet.otutable,qvalue=0,tree=bat.diet.tree) #Three-level example (L1=sample, L2=species, L3=whole system) div_part(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy)
Create diversity profile vectors (single sample or system) or tables (multiple samples or groups) from count tables.
div_profile(count, qvalues, tree, hierarchy, level)
div_profile(count, qvalues, tree, hierarchy, level)
count |
A vector or a matrix indicating the (relative) OTU/ASV counts of one or multiple samples. If a matrix is provided, columns must refer to samples and rows to OTUs. |
qvalues |
A vector of sequential orders of diversity (default from 0 to 5). qvalue=seq(from = 0, to = 5, by = (0.1)) |
tree |
A tree of class 'phylo'. The tip labels must match the names of the vector values (if one sample) or matrix rows (if multiple samples). |
hierarchy |
A two-column matrix indicating the relation between samples (first column) and groups (second column). |
level |
Whether to compute alpha or gamma diversities of the system or the groups specified in the hierarchy table. |
Diversity profile
A vector or matrix containing diversity values at different orders of diversity (as specified in qvalues).
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Jost, L. (2014). Unifying species diversity, phylogenetic diversity, functional diversity, and related similarity and differentiation measures through hill numbers. Annual Review of Ecology Evolution and Systematics, 45, 297-324.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample example bat.diet.sample <- bat.diet.otutable[,1] div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1))) #One sample example (phylogenetic Hill numbers) names(bat.diet.sample) <- rownames(bat.diet.otutable) div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1)),tree=bat.diet.tree) #Multiple samples div_profile(bat.diet.otutable) #Multiple groups (gamma diversity) div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="gamma") #Multiple groups (alpha diversity) div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="alpha")
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample example bat.diet.sample <- bat.diet.otutable[,1] div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1))) #One sample example (phylogenetic Hill numbers) names(bat.diet.sample) <- rownames(bat.diet.otutable) div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1)),tree=bat.diet.tree) #Multiple samples div_profile(bat.diet.otutable) #Multiple groups (gamma diversity) div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="gamma") #Multiple groups (alpha diversity) div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="alpha")
Plot diversity profiles from objects generated with the function div_profile().
div_profile_plot(profile, colour, log, legend)
div_profile_plot(profile, colour, log, legend)
profile |
A div_profile() object or a vector/matrix containg diversity profile(s), with columns indicating samples/groups and rows indicating orders of diversity (q-values). |
colour |
A vector of RGB colours e.g. c("#34k235","#99cc00"). The number of vector items, must equal the number of samples or groups that are intended to plot. |
log |
Whether to transform Hill numbers to logarithmic scale (TRUE) or not (FALSE). This is useful when there are large differences between q values (e.g. sharp drop from q=0 to q=1), which might complicate visualization. Default: log=FALSE |
legend |
Whether to display the legend (TRUE) or not (FALSE) in diversity profiles containing multiple samples/groups. Default TRUE in multi-sample charts. |
Diversity profile plot
A diversity profile plot.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Jost, L. (2014). Unifying species diversity, phylogenetic diversity, functional diversity, and related similarity and differentiation measures through hill numbers. Annual Review of Ecology Evolution and Systematics, 45, 297-324.
data(bat.diet.otutable) data(bat.diet.hierarchy) #One sample example bat.diet.sample <- bat.diet.otutable[,1] profile.onesample <- div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1))) div_profile_plot(profile.onesample) #Multiple samples profile.multiplesamples <- div_profile(bat.diet.otutable) div_profile_plot(profile.multiplesamples) #Multiple groups (gamma diversity) profile.multiplegroups <- div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="gamma") div_profile_plot(profile.multiplegroups)
data(bat.diet.otutable) data(bat.diet.hierarchy) #One sample example bat.diet.sample <- bat.diet.otutable[,1] profile.onesample <- div_profile(count=bat.diet.sample,qvalues=seq(from = 0, to = 5, by = (0.1))) div_profile_plot(profile.onesample) #Multiple samples profile.multiplesamples <- div_profile(bat.diet.otutable) div_profile_plot(profile.multiplesamples) #Multiple groups (gamma diversity) profile.multiplegroups <- div_profile(bat.diet.otutable,hierarchy=bat.diet.hierarchy,level="gamma") div_profile_plot(profile.multiplegroups)
Diversity comparison test between groups of samples. The function automatically assesses whether the data meets the properties for parametric statistics and performs the appropriate test accordingly: Students' T, ANOVA, Wilcoxon or Kruskal-Wallis. If the posthoc argument is set as TRUE, multiple group comparisons are complemented with post hoc pairwise tests, either Tukey test (parametric) or Dunn test with Benjamini-Hochberg correction (non-parametric).
div_test(countable, qvalue, hierarchy, tree, posthoc)
div_test(countable, qvalue, hierarchy, tree, posthoc)
countable |
A matrix indicating the relative abundances of multiple samples. Columns should be samples and rows OTUs. |
qvalue |
A positive integer or decimal number (>=0), usually between 0 and 3. |
hierarchy |
A two-column matrix indicating the relation between samples (first column) and groups (second column). |
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the OTU table. Use the function match_data() if the OTU names do not match. |
posthoc |
Whether to run post hoc pairwise analyses or not. If TRUE, an ANOVA will be complemented with a Tukey test and a Kruskal-Wallis test will be complemented with a Dunn test. |
Diversity test
A statistical test output.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Jost, L. (2014). Unifying species diversity, phylogenetic diversity, functional diversity, and related similarity and differentiation measures through hill numbers. Annual Review of Ecology Evolution and Systematics, 45, 297-324.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) div_test(bat.diet.otutable,qvalue=1,hierarchy=bat.diet.hierarchy,tree=bat.diet.tree) div_test(bat.diet.otutable,2,bat.diet.hierarchy,bat.diet.tree) div_test(bat.diet.otutable,qvalue=1,hierarchy=bat.diet.hierarchy,posthoc=TRUE)
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) div_test(bat.diet.otutable,qvalue=1,hierarchy=bat.diet.hierarchy,tree=bat.diet.tree) div_test(bat.diet.otutable,2,bat.diet.hierarchy,bat.diet.tree) div_test(bat.diet.otutable,qvalue=1,hierarchy=bat.diet.hierarchy,posthoc=TRUE)
Plot of diversity comparison between groups of samples
div_test_plot(divtest, chart, colour, posthoc, threshold)
div_test_plot(divtest, chart, colour, posthoc, threshold)
divtest |
Object outputed by the div_test() function |
chart |
Chart type, either 'box' for boxplot, 'jitter' for jitter plot or 'violin' for violin plot. chart="box" |
colour |
The number of vector items (colours, e.g. '#34k235'), must equal the number of groups that are intended to plot. |
posthoc |
If 'TRUE' pairwise p-values of the posthoc analyses will be ploted. It requires the div_test() object to contain posthoc results. |
threshold |
Maximum p-value to show in pairwise posthoc results (usually 0.05, but could be any other number between 0 an 1). P-values above the threshold will not be showed. |
Diversity test plotting
Chart of (mean) diversities of contrasting groups with optional posthoc results.
Antton Alberdi, [email protected]
data(bat.diet.otutable) data(bat.diet.hierarchy) divtestres <- div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) div_test_plot(divtestres,chart="box") div_test_plot(divtestres,chart="violin") divtest.res.ph <- div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy,posthoc=TRUE) div_test_plot(divtest.res.ph,chart="jitter",posthoc=TRUE,threshold=0.5)
data(bat.diet.otutable) data(bat.diet.hierarchy) divtestres <- div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) div_test_plot(divtestres,chart="box") div_test_plot(divtestres,chart="violin") divtest.res.ph <- div_test(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy,posthoc=TRUE) div_test_plot(divtest.res.ph,chart="jitter",posthoc=TRUE,threshold=0.5)
Compute gamma diversity of a system from a matrix (OTU table) containing multiple samples. If a tree is provided, the computed gamma diversity accounts for the phylogenetic relations across OTUs.
gamma_div(countable,qvalue,tree,weight)
gamma_div(countable,qvalue,tree,weight)
countable |
A count table (matrix/data.frame) indicating the absolute or relative OTU/ASV abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
qvalue |
A positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the count table. Use the function match_data() if the count table and tree names do not match. |
weight |
A vector indicating the relative weight of each sample. The order needs to be identical to the order of the samples in the OTU table. The values need to sum up to 1. If empty, all samples are weighed the same. |
Gamma diversity computation (based on Hill numbers)
A gamma diversity value.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
div_part
, alpha_div
, match_data
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) gamma_div(countable=bat.diet.otutable,qvalue=1) gamma_div(countable=bat.diet.otutable,qvalue=1,tree=bat.diet.tree) weight.vector = rep(1/ncol(bat.diet.otutable),ncol(bat.diet.otutable)) gamma_div(bat.diet.otutable,1,bat.diet.tree,weight.vector)
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) gamma_div(countable=bat.diet.otutable,qvalue=1) gamma_div(countable=bat.diet.otutable,qvalue=1,tree=bat.diet.tree) weight.vector = rep(1/ncol(bat.diet.otutable),ncol(bat.diet.otutable)) gamma_div(bat.diet.otutable,1,bat.diet.tree,weight.vector)
Compute neutral or phylogenetic Hill numbers from a single sample (vector) or count table (matrix). Hill numbers or numbers equivalents of diversity indices are diversity measures that compute diversity in effective number of OTUs, i.e. the number of equally abundant OTUs that would be needed to give the same value of diversity.
hill_div(count, qvalue, tree, dist)
hill_div(count, qvalue, tree, dist)
count |
A vector or a matrix/data.frame indicating the (relative) counts of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
qvalue |
A positive integer or decimal number (>=0), usually between 0 and 3. |
tree |
An ultrametic tree of class 'phylo'. The tip labels must match the names of the vector values (if one sample) or matrix rows (if multiple samples). Use the function match_data() if the OTU names do not match. |
dist |
A dist object indicating the pairwise distances between samples. NOT implemented yet |
Hill numbers computation
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Jost, L. (2006). Entropy and diversity. Oikos, 113, 363-375.
Hill, M. O. (1973). Diversity and evenness: a unifying notation and its consequences. Ecology, 54, 427-432.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample bat.diet.sample <- bat.diet.otutable[,1] hill_div(bat.diet.sample,0) hill_div(bat.diet.sample,qvalue=1) #One sample (phylogenetic) names(bat.diet.sample) <- rownames(bat.diet.otutable) hill_div(bat.diet.sample,1,bat.diet.tree) #Multiple samples hill_div(bat.diet.otutable,0) #Incidence-based bat.diet.otutable.incidence <- to.incidence(bat.diet.otutable,bat.diet.hierarchy) hill_div(bat.diet.otutable.incidence,qvalue=1) hill_div(to.incidence(bat.diet.otutable,bat.diet.hierarchy),1)
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample bat.diet.sample <- bat.diet.otutable[,1] hill_div(bat.diet.sample,0) hill_div(bat.diet.sample,qvalue=1) #One sample (phylogenetic) names(bat.diet.sample) <- rownames(bat.diet.otutable) hill_div(bat.diet.sample,1,bat.diet.tree) #Multiple samples hill_div(bat.diet.otutable,0) #Incidence-based bat.diet.otutable.incidence <- to.incidence(bat.diet.otutable,bat.diet.hierarchy) hill_div(bat.diet.otutable.incidence,qvalue=1) hill_div(to.incidence(bat.diet.otutable,bat.diet.hierarchy),1)
Computes common diversity indices related to Hill numbers. If the input is a vector, the function computes the indices of a single sample, while if the input is a matrix (OTU table), the function computes individual diversity indices for each sample (column). An ultrametic OTU tree is required for computing phylogenetic diversity indices (Faith's PD, Allen's H and Rao's Q). If the relative abundances of each sample (vector or each column of the matrix) do not sum to 1, TSS normalisation is applied.
index_div(abund, tree, index)
index_div(abund, tree, index)
abund |
A vector or a matrix/data.frame indicating the relative abundances of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
tree |
An ultrametic tree of class 'phylo'. The tip labels must match the names of the vector values (if one sample) or matrix rows (if multiple samples). Use the function match_data() if the OTU names do not match. |
index |
Diversity index to be computed ("richness", "shannon", "simpson", "faith", "allen", "rao"). Default without tree argument: index="richness". Default with tree argument: index="faith". |
Diversity index computation
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Jost, L. (2006). Entropy and diversity. Oikos, 113, 363-375.
Rao, C. R. (1982). Diversity and dissimilarity coefficients: A unified approach. Theoretical Population Biology, 21, 24-43.
Shannon, C. E. (1948). A mathematical theory of communication. The Bell System Technical Journal, 27, 379-423.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample bat.diet.sample <- bat.diet.otutable[,1] index_div(bat.diet.sample) index_div(bat.diet.sample,index="shannon") #Multiple samples index_div(bat.diet.otutable) index_div(bat.diet.otutable,tree=bat.diet.tree,index="faith") #Incidence-based bat.diet.otutable.incidence <- to.incidence(bat.diet.otutable,bat.diet.hierarchy) index_div(bat.diet.otutable.incidence) index_div(bat.diet.otutable.incidence,index="simpson") index_div(to.incidence(bat.diet.otutable,bat.diet.hierarchy),tree=bat.diet.tree)
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) #One sample bat.diet.sample <- bat.diet.otutable[,1] index_div(bat.diet.sample) index_div(bat.diet.sample,index="shannon") #Multiple samples index_div(bat.diet.otutable) index_div(bat.diet.otutable,tree=bat.diet.tree,index="faith") #Incidence-based bat.diet.otutable.incidence <- to.incidence(bat.diet.otutable,bat.diet.hierarchy) index_div(bat.diet.otutable.incidence) index_div(bat.diet.otutable.incidence,index="simpson") index_div(to.incidence(bat.diet.otutable,bat.diet.hierarchy),tree=bat.diet.tree)
Multi-level diversity partitioning requires the groups at different hierarchical levels to be nested. i.e. two samples that belong to a common parent group cannot have different grandparent groups. The best example of nested hierarchy is taxonomy: e.g. two species that belong to the same genus cannot belong to different families. This function checks whether the groups specified in a hierarchy table have a nested structure.
is.nested(hierarchy)
is.nested(hierarchy)
hierarchy |
A matrix indicating the relation between samples (first column) and parent groups. |
Check if hierachy is nested
A logical value (TRUE/FALSE).
Antton Alberdi, [email protected]
data(bat.diet.hierarchy) is.nested(bat.diet.hierarchy)
data(bat.diet.hierarchy) is.nested(bat.diet.hierarchy)
Filter count tables and OTU/ASV phylogenetic trees to match OTUs/ASVs present in both data files..
match_data(countable, tree, output, silent)
match_data(countable, tree, output, silent)
countable |
A count table (matrix/data.frame) indicating the absolute or relative OTU/ASV abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
tree |
An ultrametic tree of class 'phylo'. |
output |
Whether to output a filtered count table ('countable') or a filtered OTU tree ('tree'). Default is empty, which only yields a message. |
silent |
Whether to stop printing text on screen. Default=FALSE. |
Match data
Antton Alberdi, [email protected]
data(bat.diet.otutable) data(bat.diet.tree) match_data(bat.diet.otutable,bat.diet.tree,output="countable") match_data(bat.diet.otutable,bat.diet.tree,output="tree")
data(bat.diet.otutable) data(bat.diet.tree) match_data(bat.diet.otutable,bat.diet.tree,output="countable") match_data(bat.diet.otutable,bat.diet.tree,output="tree")
Combines samples into groups defined by the hierarchy table, with the possibility to convert abundances into incidence data.
merge_samples(countable,hierarchy,incidence)
merge_samples(countable,hierarchy,incidence)
countable |
A count table (matrix/data.frame) indicating the absolute or relative OTU/ASV abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
hierarchy |
A two-column matrix indicating the relation between samples (first column) and groups (second column). |
relative |
Whether to output relative values or not. Default=TRUE. |
incidence |
Whether to transform abundance into incidence data when merging. Default=FALSE. |
Merge samples
A count table
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
data(bat.diet.otutable) data(bat.diet.hierarchy) merge_samples(countable=bat.diet.otutable,hierarchy=bat.diet.hierarchy) merge_samples(bat.diet.otutable,bat.diet.hierarchy) merge_samples(countable=bat.diet.otutable,hierarchy=bat.diet.hierarchy, incidence=TRUE)
data(bat.diet.otutable) data(bat.diet.hierarchy) merge_samples(countable=bat.diet.otutable,hierarchy=bat.diet.hierarchy) merge_samples(bat.diet.otutable,bat.diet.hierarchy) merge_samples(countable=bat.diet.otutable,hierarchy=bat.diet.hierarchy, incidence=TRUE)
Computation of pairwise dissimilarities based on Hill numbers diversity partitioning
pair_dis(countable, qvalue, tree, hierarchy, metric)
pair_dis(countable, qvalue, tree, hierarchy, metric)
countable |
A matrix indicating the relative abundances of multiple samples. Columns should be samples and rows OTUs. |
qvalue |
A positive integer or decimal number (>=0), usually between 0 and 3. |
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the OTU table. Use the function match_data() if the OTU names do not match. |
hierarchy |
A matrix indicating the relation between samples (first column) and groups. |
metric |
A vector containing any combination of "C", "U", "V" or "S". If not provided, all metrics will be computed. metric="U", metric=c("U","S"). |
Pairwise dissimilarity
A list of matrices containing pairwise beta diversities and dissimilarity metrics.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pair_dis(bat.diet.otutable,qvalue=1) pair_dis(bat.diet.otutable,qvalue=1,tree=bat.diet.tree,metric="V")
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pair_dis(bat.diet.otutable,qvalue=1) pair_dis(bat.diet.otutable,qvalue=1,tree=bat.diet.tree,metric="V")
Visualisation of pairwise dissimilarities. Not this function is deprecated. Use dis_nmds() and dis_network() instead.
pair_dis_plot(distance, hierarchy, type, colour, magnify)
pair_dis_plot(distance, hierarchy, type, colour, magnify)
distance |
Matrix of pairwise dissimilarities, usually one of the matrices listed in the output object of the pair_dis() function. |
hierarchy |
The first column lists the sample names while the second lists the groups. If provided, group profiles are plotted instead of individual profiles. |
type |
Whether to plot a NMDS, Shepard or qgraph chart. type="NMDS". |
colour |
he number of vector items (colours, e.g. '#34k235'), must equal the number of samples or groups that are intended to plot with different colours. |
magnify |
Only relevant for qgraph. Whether the pairwise dissimilarity values are transformed to 0-1 scale, 0 corresponding to the minimum dissimilarity and 1 to the maximum dissimilarity value. magnify=FALSE. |
Pairwise dissimilarity plot (deprecated)
An NMDS or network plot.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pairdisres <- pair_dis(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) pair_dis_plot(pairdisres$L2_CqN,hierarchy=bat.diet.hierarchy,type="NMDS") pair_dis_plot(pairdisres$L2_CqN,type="qgraph") pair_dis_plot(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy,type="qgraph")
data(bat.diet.otutable) data(bat.diet.tree) data(bat.diet.hierarchy) pairdisres <- pair_dis(bat.diet.otutable,qvalue=0,hierarchy=bat.diet.hierarchy) pair_dis_plot(pairdisres$L2_CqN,hierarchy=bat.diet.hierarchy,type="NMDS") pair_dis_plot(pairdisres$L2_CqN,type="qgraph") pair_dis_plot(pairdisres$L1_CqN,hierarchy=bat.diet.hierarchy,type="qgraph")
The Jaccard-type turnover-complement is thecomplement of the Jaccard-type turnover, which quantifies the normalized OTU turnover rate with respect to the whole system (i.e. gamma). SqN is integrated in the functions beta_dis() and pair_dis().
SqN(beta, N)
SqN(beta, N)
beta |
A beta diversity value based on Hill numbers. |
N |
An integer indicating sample size, the number of sampling units to be used to compute the similarity measure. |
Jaccard-type turnover-complement
A Jaccard-type turnover-complement value
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
SqN(beta=1.24,N=2) SqN(1.24,2)
SqN(beta=1.24,N=2) SqN(1.24,2)
Transform a count (OTU/ASV) table from abundance to incidence.
to.incidence(otutable, hierarchy, relative)
to.incidence(otutable, hierarchy, relative)
otutable |
A matrix/data.frame indicating the (relative) abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
hierarchy |
A two-column matrix indicating the relation between samples (first column) and groups (second column). |
relative |
Whether to transform the incidence vector or matrix to relative (0-1) values. Default: relative=FALSE. |
To incidence
A vector of incidence data of a single system if no hierarchy table is specified and a matrix of incidence data of multiple systems if a hierarchy table is specified.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
data(bat.diet.otutable) data(bat.diet.hierarchy) to.incidence(bat.diet.otutable) to.incidence(bat.diet.otutable,bat.diet.hierarchy) to.incidence(bat.diet.otutable,bat.diet.hierarchy,relative=TRUE) to.incidence(otutable=bat.diet.otutable,hierarchy=bat.diet.hierarchy,relative=TRUE)
data(bat.diet.otutable) data(bat.diet.hierarchy) to.incidence(bat.diet.otutable) to.incidence(bat.diet.otutable,bat.diet.hierarchy) to.incidence(bat.diet.otutable,bat.diet.hierarchy,relative=TRUE) to.incidence(otutable=bat.diet.otutable,hierarchy=bat.diet.hierarchy,relative=TRUE)
Transform an absolute or relative abundance vector (one sample) or matrix (multipla samples) into an occurrence vector or matrix.
to.occurrences(abund)
to.occurrences(abund)
abund |
A vector or a matrix/data.frame indicating the absolute or relative abundances of multiple samples. Columns must refer to samples and rows to OTUs/ASVs. |
To occurrences
A vector (one sample) or matrix (multiple samples) of occurrence data.
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
data(bat.diet.otutable) to.occurrences(bat.diet.otutable) to.occurrences(bat.diet.otutable[,1])
data(bat.diet.otutable) to.occurrences(bat.diet.otutable) to.occurrences(bat.diet.otutable[,1])
Computes phylogenetic tree depth based from a phylogenetic tree and a vector of (relative) abundances.
tree_depth(tree, abund)
tree_depth(tree, abund)
tree |
A phylogenetic tree of class 'phylo'. The tip labels must match the row names in the OTU table. Use the function match_data() if the OTU names do not match. |
abund |
A vector or a matrix/data.frame indicating the relative abundances of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
Tree depth
A tree depth value
Antton Alberdi, [email protected]
div_part
, gamma_div
, match_data
data(bat.diet.otutable) data(bat.diet.tree) tree_depth(tree=bat.diet.tree,abund=bat.diet.otutable) tree_depth(bat.diet.tree,bat.diet.otutable)
data(bat.diet.otutable) data(bat.diet.tree) tree_depth(tree=bat.diet.tree,abund=bat.diet.otutable) tree_depth(bat.diet.tree,bat.diet.otutable)
Normalise a vector or count matrix to the range of 0-1.
tss(abund)
tss(abund)
abund |
A vector or a matrix/data.frame indicating the relative abundances of one or multiple samples, respectively. If a matrix/data.frame is provided, columns must refer to samples and rows to OTUs. |
Total Sum Scaling normalisation
Normalised vector or matrix.
Antton Alberdi, [email protected]
data(bat.diet.otutable) tss(bat.diet.otutable) bat.diet.sample <- bat.diet.otutable[,1] tss(bat.diet.sample)
data(bat.diet.otutable) tss(bat.diet.otutable) bat.diet.sample <- bat.diet.otutable[,1] tss(bat.diet.sample)
The Jaccard-type overlap quantifies the effective proportion of OTUs or lineages in a system that are shared across all subsystems. Hence, this metric quantifies overlap from the perspective of the overall system. Its corresponding dissimilarity (1 - UqN) quantifies the effective proportion of nonshared OTUs or lineages in the overall system. UqN is integrated in the functions beta_dis() and pair_dis().
UqN(beta, qvalue, N)
UqN(beta, qvalue, N)
beta |
A beta diversity value based on Hill numbers. |
qvalue |
The q value used to compute the beta diversity. It needs to be a positive number, usually between 0 and 5, but most commonly 0, 1 or 2. It can be an integer or contain decimals. |
N |
An integer indicating sample size, the number of sampling units to be used to compute the similarity measure. |
Jaccard-type overlap
A Jaccard-type overlap value
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
UqN(beta=1.24,qvalue=1,N=2) UqN(1.24,1,2)
UqN(beta=1.24,qvalue=1,N=2) UqN(1.24,1,2)
The Sørensen-type turnover-complement is the complement of the Sørensen-type turnover, which quantifies the normalized OTU turnover rate with respect to the average subsystem (i.e., alpha), thus provides the proportion of a typical subsystem that changes across subsystems. VqN is integrated in the functions beta_dis() and pair_dis().
VqN(beta, N)
VqN(beta, N)
beta |
A beta diversity value based on Hill numbers. |
N |
An integer indicating sample size, the number of sampling units to be used to compute the similarity measure. |
Sørensen-type turnover-complement
A Sørensen-type turnover-complement value
Antton Alberdi, [email protected]
Alberdi, A., Gilbert, M.T.P. (2019). A guide to the application of Hill numbers to DNA-based diversity analyses. Molecular Ecology Resources, 19, 804-817.
Chao, A., Chiu, C.H., & Hsieh, T. C. (2012). Proposing a resolution to debates on diversity partitioning. Ecology, 93, 2037-2051.
Jost, L. (2007). Partitioning diversity into independent alpha and beta components. Ecology, 88, 2427-2439.
VqN(beta=1.24,N=2) VqN(1.24,2)
VqN(beta=1.24,N=2) VqN(1.24,2)