R/phylobeta_ses.R
phylobeta_ses.Rd
This function computes the standard effect size of phylogenetic beta diversity by correcting for changes in species beta diversity. 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.
phylobeta_ses(
x,
phy,
index.family = "simpson",
model = c("tipshuffle", "rowwise", "colwise"),
reps = 1000,
...
)
a (sparse) community matrix, i.e., an object of class matrix or Matrix.
a phylogenetic tree (object of class phylo).
the family of dissimilarity indices including “simpson”, “sorensen” and “jaccard”.
The null model for separating patterns from processes and for contrasting against alternative hypotheses. Available null models include:
“tipshuffle”: shuffles phylogenetic 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.
Number of replications.
Further arguments passed to or from other methods.
A data frame of results for each community or grid cell
phylobeta_obs: Observed phylobeta in community
phylobeta_rand_mean: Mean phylobeta in null communities
phylobeta_rand_sd: Standard deviation of phylobeta in null communities
phylobeta_obs_z: Standardized effect size of phylobeta vs. null communities \(= (phylobeta_obs - phylobeta_rand_mean) / phylobeta_rand_sd\)
reps: Number of replicates
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.
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))
phylobeta_ses(com, tree, model="rowwise")
#> $phylobeta_obs
#> g1 g2 g3 g4 g5
#> g2 0.0000000
#> g3 0.2000000 0.0000000
#> g4 0.0000000 0.0000000 0.0000000
#> g5 0.0000000 0.0000000 0.2500000 0.0000000
#> g6 0.3333333 0.0000000 0.0000000 0.0000000 0.3333333
#>
#> $phylobeta_rand_mean
#> g1 g2 g3 g4 g5
#> g2 0.0000000
#> g3 0.2000000 0.0000000
#> g4 0.0000000 0.0000000 0.0000000
#> g5 0.0000000 0.0000000 0.2500000 0.0000000
#> g6 0.3333333 0.0000000 0.0000000 0.0000000 0.3333333
#>
#> $phylobeta_rand_sd
#> g1 g2 g3 g4 g5
#> g2 1
#> g3 1 1
#> g4 1 1 1
#> g5 1 1 1 1
#> g6 1 1 1 1 1
#>
#> $phylobeta_obs_z
#> g1 g2 g3 g4 g5
#> g2 0
#> g3 0 0
#> g4 0 0 0
#> g5 0 0 0 0
#> g6 0 0 0 0 0
#>
#> $reps
#> [1] 1000
#>