Last data update: 2014.03.03

R: measuring spatial and phylogenetic structuring of diversity...
spacodi.calcR Documentation

measuring spatial and phylogenetic structuring of diversity in communities

Description

Considering species-, phylogenetic-, or trait-diversities, this function measures diversity structuring of community samples.

Usage

spacodi.calc(sp.plot, phy = NULL, sp.traits = NULL, all.together = TRUE, prune = TRUE, pairwise = FALSE, ...)

Arguments

sp.plot

a community dataset in spacodiR format (see as.spacodi)

phy

a phylogenetic tree of class phylo or evolutionary distance matrix between species (see cophenetic.phylo)

sp.traits

a species-by-trait(s) dataframe or a species traits distance matrix (see dist)

all.together

Boolean; whether to treat all traits together or separately

prune

Boolean; whether to dynamically prune datasets if mismatches occur

pairwise

Boolean; whether to return pairwise diversity measures amongst all plots

...

additional arguments to be passed to match.spacodi.data

Details

spacodi.calc requires a community dataset (species-by-plots matrix; sp.plot) of numerical abundance, relative abundance, or presence | absence for plots. spacodi.calc returns statistics of diversity partitioning of plots, considering species diversity and, if additional information is provided, either trait or phylogenetic diversities among plots. If phy=NULL and sp.traits=NULL, a measure of partitioning for species diversity will be returned.

In treating each pair of plots as a community unto its own, pairwise=TRUE will return estimates for diversity structuring for all pairwise combinations of plots.

If a phylogeny or trait dataset is supplied with species that are not present in the community dataset (i.e., sp.plot) or vice versa, the user has the option to dynamically prune these datasets to match (prune=TRUE). If prune=FALSE and dataset mismatches occur, the function will inevitably return NaN where plots have fewer than two distinct species sampled.

For proper display, please view the package manual online (http://cran.r-project.org/web/packages/spacodiR/spacodiR.pdf)

GLOBAL MEASURES

  • N: number of local communities sampled

  • n_k: number of individuals sampled in local community k

  • f_{ik}: observed relative abundance of species i in the local community k (∑_i{f_{ik}=1})

  • p_{ik}: presence (1) or absence (0) of species i in the local community k

  • δ_{ij}: phyletic or functional (trait) distance between species i and j

INDICES

  • indices using species abundances: Mean phyletic or functional distance between two individuals in the same local community is expressed as follows:

    D_k = frac{n_k}{n_k-1} ∑_i ∑_j{δ_{ij}f_{ik}f_{jk} }

    This measure is sample-size corrected by the term (frac{n_k}{n_k-1}) and is only applied where data are counts of individuals (not where data represent a relative measure of species abundances). Mean functional or phyletic distance between individuals of different species in the same local community is computed as follows:

    D_k* = ∑_i∑_{j \neq i}{δ_{ij}f_{ik}f_{jk} } / ∑_i∑_{j \neq i}{f_{ik}f_{jk} }

    Rao's quadratic entropy within samples (average inter-individual distance among samples) is computed either with (D_S) or without (D_S*) sample-size correction:

    D_S = frac{1}{N} ∑_{k ≤q N} D_k

    D_S^* = frac{1}{N} ∑_{k ≤q N} D_k^*

    Rao's quadratic entropy among samples (average distance between individuals from different samples) is as follows, including or excluding intraspecific comparisons (respectively):

    D_T = frac{1}{N(N-1)} ∑_k ∑_{l \neq k} ∑_i ∑_j δ_{ij}f_{ik}f_{jl}

    D_T^* = frac{1}{N(N-1)} ∑_k ∑_{l \neq k} [ ∑_i ∑_{j \neq i} δ_{ij}f_{ik}f_{jl} / ∑_i ∑_{j \neq i} f_{ik}f_{jl} ]

    • P_{ST} integrates both species turnover and phylogenetic turnover, weighting species according to local abundances.

      P_{ST}=T_{ST}=1-frac{D_S}{D_T}

      for phylogenetic (P_{ST}) or functional (T_{ST}) distances between individuals.

    • I_{ST} is a measure of species turnover accounting for local abundances, but with different properties than measures of species turnover based on species sharing. It is the analogue of G_{ST} in population genetics.

      I_{ST}=1-frac{D_S}{D_T}

      in the special case of spatial partitioning of species diversity, where δ_{ij}=0 if i=j and where δ_{ij}=1 if i \neq j

    • B_{ST} and U_{ST} integrate only phylogenetic or functional (phenotypic) turnover, respectively, weighting species according to local abundances.

      B_{ST}=U_{ST}=1-frac{D_S^*}{D_T^*}

      for phylogenetic (B_{ST}) or functional (U_{ST}) distances between individuals, excluding intraspecific comparisons.

  • indices using species occurrences: If we let p_{ik}=1 where species i occurs in sample k and p_{ik}=0 otherwise, we compute the following local measure of phyletic or functional distinctiveness:

    Δ_k=∑_i ∑_{j \neq i} {δ_{ij}p_{ik}p_{jk} } / ∑_i ∑_{j \neq i} {p_{ik}p_{jk} }

    Δ_k is the mean interspecific distance in local community k and mean distance between distinct species within local communities is then

    Δ_S = frac{1}{N} ∑_{k ≤q N} δ_k

    Mean distance between distinct species among local communities is:

    Δ_T = frac{1}{N} ∑_k ∑_{l \neq k} [ ∑_i ∑_{j \neq i} δ_{ij}p_{ik}p_{jl} / ∑_i ∑_{j \neq i} p_{ik}p_{jl} ]

    This formula corrects a typographic error in eq. 16 of Hardy and Senterre (2007) where j \neq i were missing.

    Π_{ST} expresses phylogenetic turnover without accounting for species local abundances (it is a presence-absence version of B_{ST}):

    Π_{ST} = 1 - frac{Δ_S}{Δ_T}

  • assumptions:

    • N is much is much smaller than the total number of local communities constituting the whole community of which the sampled local communities are representative (this justifies that Δ_T and Δ_T are estimated excluding pairs of species from the same local community)

    • n_k is much smaller than the total number of individuals constituting local community k (this justifies the frac{n_k}{n_k-1} correction factors)

Value

A named list of at least one element (Ist) is returned. The size of the returned list is wholly dependent upon given arguments.

SPECIES DIVERSITY STRUCTURING

  • Ist: a measure of local species identity excess between individuals, expressing species turnover. It is a form of spatial partition of Gini-Simpson diversity (equivalent to Fst in population genetics). Ist considers only abundances (or presences) in the species-by-plots matrix.

PHYLOGENETIC DIVERSITY STRUCTURING

  • Pst: a measure of local phyletic proximity excess between individuals, expressing species + phylogenetic turnover. It is a form of spatial partition of Rao's quadratic entropy (equivalent to Nst in population genetics). Tst is the analogue for trait data, estimating the spatial partitioning of mean trait-divergence between individuals.

  • Bst: a measure of local phyletic proximity excess between individuals of distinct species, expressing phylogenetic turnover (independently of species turnover). Ust is the analogue for trait data, estimating the spatial partitioning of mean trait-divergence between individuals that belong to distinct species.

  • PIst: Bst analogue for presence/absence data, expressing phylogenetic turnover (independently of species turnover). TAUst is the analogue for trait data, estimating mean trait-divergence between distinct species.

TRAIT DIVERSITY STRUCTURING

  • Measures analogous to those under PHYLOGENETIC DIVERSITY STRUCTURING can be computed from trait data. For trait data, these analogues are Tst (see Pst), Ust (see Bst), and TAUst (see PIst). Note: elsewhere, Ust will be referred to as T*st but here has been renamed to avoid issues of indexing in R. Trait values are not assumed to follow any particular model of evolution; rather, distances between observed species traits are expected to be uniform in distribution.

  • If all.together=TRUE, all traits will be used to generate distance a distance matrix for sampled species. Where all.together=FALSE is used, output is generated for each trait independently.

INTERPRETATION

  • spatial clustering: species within plots are more phylogenetically related on average than species from distinct plots where Pst > Ist, Bst > 0, or PIst > 0. Species are functionally more similar locally than those from distinct plots where Tst > Ist, Ust > 0, or TAUst > 0

  • spatial overdispersion: species within plots are less phylogenetically related on average than species from distinct plots where Pst < Ist, Bst < 0, or PIst < 0. Species are functionally less similar locally than are species from distinct plots where Tst < Ist, Ust < 0, or TAUst < 0

Author(s)

Olivier Hardy, Timothy Paine, and Jonathan Eastman

References

HARDY OJ and B SENTERRE. 2007. Characterizing the phylogenetic structure of communities by an additive partitioning of phylogenetic diversity. Journal of Ecology 95:493-506.

HARDY OJ. 2008. Testing the spatial phylogenetic structure of local communities: statistical performances of different null models and test statistics on a locally neutral community. Journal of Ecology 96:914-926.

HARDY OJ and L JOST. 2008. Interpreting and estimating measures of community phylogenetic structuring. Journal of Ecology 96:849-852.

See Also

match.spacodi.data; as.spacodi

Examples

# load a species-by-plots matrix, along with a tree
data(sp.example)
attributes(sp.example)
attach(sp.example)
spl
phy

# community diversity statistics of Hardy and Senterre (2007): tree-based
spacodi.calc(sp.plot = spl, phy = phy)

# community diversity statistics: trait-based with pairwise comparisons
spacodi.calc(sp.plot = spl, phy = phy, pairwise=TRUE)

# community diversity for a pair of traits
spacodi.calc(sp.plot = spl, sp.traits = trt, all.together=TRUE)

# community diversity for a pair of traits, each singly
spacodi.calc(sp.plot = spl, sp.traits = trt, all.together=FALSE)

# Ist: using abundance data only				
spacodi.calc(sp.plot = spl)	

# calculations with missing taxa between tree and sp.plot
# excluding the last five species in sp.plot, 
spacodi.calc(sp.plot = spl[1:15,], phy = phy, prune=TRUE)

# as before but with 'manual' pruning of the datasets
match.spacodi.data(sp.plot=spl[1:15,],phy=phy) -> prn.data
spacodi.calc(sp.plot=prn.data$sp.plot, phy=prn.data$sp.tree)
prn.data$sp.plot
prn.data$sp.tree


							

Results