PD_ses: Phylogenetic diversity standardized for species richness

View source: R/PD_ses.R

PD_sesR Documentation

Phylogenetic diversity standardized for species richness

Description

This function computes the standard effect size of PD by correcting for changes in species richness. The novelty of this function is its ability to utilize sparse community matrix making it possible to efficiently randomize very large community matrices spanning thousands of taxa and sites.

Usage

PD_ses(
  x,
  phy,
  model = c("tipshuffle", "rowwise", "colwise"),
  reps = 10,
  metric = "pd",
  ...
)

Arguments

x

a (sparse) community matrix, i.e. an object of class matrix or Matrix.

phy

a phylogenetic tree (object of class phylo).

model

The null model for separating patterns from processes and for contrasting against alternative hypotheses. Available null models include:

  • “tipshuffle”: shuffles tip labels multiple times.

  • “rowwise”: shuffles sites (i.e., varying richness) and keeping species occurrence frequency constant.

  • “colwise”: shuffles species occurrence frequency and keeping site richness constant.

reps

Number of replications.

metric

The phylodiversity measure to compute.

...

Further arguments passed to or from other methods.

Value

A data frame of results for each community or grid cell

  • grids: Site identity

  • richness: Number of taxa in community

  • pd_obs: Observed PD in community

  • pd_rand.mean: Mean PD in null communities

  • pd_rand.sd: Standard deviation of PD in null communities

  • pd_obs.rank: Rank of observed PD vs. null communities

  • pd_obs.z: Standardized effect size of PD vs. null communities = (pd_obs - pd_rand.mean) / pd_rand_sd

  • pvalue: P-value (quantile) of observed PD vs. null communities = mpd_obs_rank / iter + 1

  • reps: Number of replicates

  • p_obs_c_lower: Number of times observed value < random value

  • p_obs_c_upper: Number of times observed value > random value

  • p_obs_p_lower: Percentage of times observed value < random value

  • p_obs_p_upper: Percentage of times observed value > random value

  • p_obs_q: Number of the non-NA random values used for comparison

References

Proches, S., Wilson, J.R.U. & Cowling, R.M. (2006) How much evolutionary history in a 10 x 10m plot? Proceedings of Royal Society B 273: 1143-1148.

Examples

library(ape)
library(Matrix)
tree <- read.tree(text ="((t1:1,t2:1)N2:1,(t3:1,t4:1)N3:1)N1;")
com <- sparseMatrix(c(1,3,4,1,4,5,1,2,3,4,5,6,3,4,6),
  c(1,1,1,2,2,2,3,3,3,3,3,3,4,4,4),x=1,
  dimnames = list(paste0("g", 1:6), tree$tip.label))

PD_ses(com, tree, model="rowwise")


phyloregion documentation built on Aug. 15, 2023, 9:07 a.m.