Title: | Evolutionary Diversity Metrics for Raster Data |
---|---|
Description: | Phylogenetic Diversity (PD, Faith 1992), Evolutionary Distinctiveness (ED, Isaac et al. 2007), Phylogenetic Endemism (PE, Rosauer et al. 2009; Laffan et al. 2016), and Weighted Endemism (WE, Laffan et al. 2016) for presence-absence raster. Faith, D. P. (1992) <doi:10.1016/0006-3207(92)91201-3> Isaac, N. J. et al. (2007) <doi:10.1371/journal.pone.0000296> Laffan, S. W. et al. (2016) <doi:10.1111/2041-210X.12513> Rosauer, D. et al. (2009) <doi:10.1111/j.1365-294X.2009.04311.x>. |
Authors: | Gabriela Alves-Ferreira [aut, cre, cph]
|
Maintainer: | Gabriela Alves-Ferreira <[email protected]> |
License: | GPL (>= 3) |
Version: | 2.2.0 |
Built: | 2025-02-28 15:25:37 UTC |
Source: | https://github.com/gabferreira/phyloraster |
Check for missing arguments using function call and a provided vector with argument names to check
arg.check( call, arguments = c("LR", "inv.R", "branch.length", "n.descen", "tree") )
arg.check( call, arguments = c("LR", "inv.R", "branch.length", "n.descen", "tree") )
call |
match.call(). To get function call with all of the specified arguments and their full names. |
arguments |
character. Arguments to be checked |
logical
Neander Marcel Heming
geop <- function(x, tree, ...){ f4 <- arg.check(match.call(), c("LR", "inv.R", "branch.length", "n.descen")) f1 <- arg.check(match.call(), c("tree")) c(f1, f4) } geop(1, 1) geop(1) geop(1, LR=1)
geop <- function(x, tree, ...){ f4 <- arg.check(match.call(), c("LR", "inv.R", "branch.length", "n.descen")) f1 <- arg.check(match.call(), c("tree")) c(f1, f4) } geop(1, 1) geop(1) geop(1, LR=1)
Use the results of rast.pe.ses() to identify centers of paleo-, neo-, super-, and mixed- endemism following the CANAPE scheme of Mishler et al., 2014.
canape.rast( pe.obs.p.upper, pe.alt.obs.p.upper, rpe.obs.p.upper, rpe.obs.p.lower, filename = NULL, overwrite = FALSE )
canape.rast( pe.obs.p.upper, pe.alt.obs.p.upper, rpe.obs.p.upper, rpe.obs.p.lower, filename = NULL, overwrite = FALSE )
pe.obs.p.upper |
SpatRaster. Upper p-value comparing the observed phylogenetic endemism and the randomized phylogenetic endemism values |
pe.alt.obs.p.upper |
SpatRaster. Upper p-value comparing the alternate phylogenetic endemism and the randomized alternate phylogenetic endemism |
rpe.obs.p.upper |
SpatRaster. Upper p-value comparing the relative phylogenetic endemism and the randomized relative phylogenetic endemism |
rpe.obs.p.lower |
SpatRaster. Lower p-value comparing the relative phylogenetic endemism and the randomized relative phylogenetic endemism |
filename |
character. Output filename |
overwrite |
logical. If TRUE, filename is overwritten |
SpatRaster
library(SESraster) library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x, tree) ses <- rast.pe.ses(x = data$x, data$tree, aleats = 5, metric = "all") # CANAPE canape <- canape.rast(ses$p.upper.PE, ses$p.upper.PE.alt, ses$p.upper.RPE, ses$p.lower.RPE) unique(canape) plot(canape)
library(SESraster) library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x, tree) ses <- rast.pe.ses(x = data$x, data$tree, aleats = 5, metric = "all") # CANAPE canape <- canape.rast(ses$p.upper.PE, ses$p.upper.PE.alt, ses$p.upper.RPE, ses$p.lower.RPE) unique(canape) plot(canape)
A dataset containing presence-absence of 33 Australian tree frogs. This dataset also provide coordinates x and y for each site.
dataR
dataR
A matrix with 2891 rows and 35 columns.
Rosauer, 2017. Available on: https://github.com/DanRosauer/phylospatial/tree/master/PhyloEndemism_in_R/Tree%20Frog%20Data
Calculates the difference of rasterized diversity metrics (richness, phylogenetic endemism, phylogenetic diversity, weighted endemism, evolutionary distinctiveness) between time periods.
delta.grid(r1, r2, filename = NULL, cores = 1, ...)
delta.grid(r1, r2, filename = NULL, cores = 1, ...)
r1 |
SpatRaster. Rasterized diversity metrics for time 1 (e.g phylogenetic diversity in present). To calculate some diversity metrics for rasters see phyloraster::geo.phylo function. |
r2 |
SpatRaster. Rasterized diversity metrics for time 2 (e.g phylogenetic diversity in future). |
filename |
character. Output filename. |
cores |
positive integer. If cores > 1, a 'parallel' package cluster with that many cores is created and used. |
... |
additional arguments to be passed passed down from a calling function. |
The two input rasters (r1 and r2) must have the same extent and resolution.
SpatRaster
# data x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) # metric SR richness riq.pres <- rast.sr(x) # imagine we lost some species in the future riq.fut <- rast.sr(x[[c(1:15)]]) dg <- delta.grid(riq.pres, riq.fut) terra::plot(dg)
# data x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) # metric SR richness riq.pres <- rast.sr(x) # imagine we lost some species in the future riq.fut <- rast.sr(x[[c(1:15)]]) dg <- delta.grid(riq.pres, riq.fut) terra::plot(dg)
The function transforms a data.frame or a matrix of presence- absence in a raster of distribution.
df2rast(x, CRS = "+proj=longlat +datum=WGS84", ...)
df2rast(x, CRS = "+proj=longlat +datum=WGS84", ...)
x |
data.frame. A data.frame or matrix with species names in columns and sites in rows. The first two columns must provide longitude and latitude, respectively. |
CRS |
character. Description of the Coordinate Reference System (map projection) in PROJ.4. |
... |
additional arguments to be passed passed down from a calling function. |
SpatRaster
dat <- phyloraster::load.data.rosauer() df2rast(dat$presab, crs = "+proj=longlat +datum=WGS84 +ellps=WGS84 +towgs84=0,0,0")
dat <- phyloraster::load.data.rosauer() df2rast(dat$presab, crs = "+proj=longlat +datum=WGS84 +ellps=WGS84 +towgs84=0,0,0")
Calculate species richness, phylogenetic diversity, evolutionary distinctiveness, phylogenetic endemism and weighted endemism using rasters as input.
geo.phylo( x, tree, inv.R, edge.path, branch.length, n.descen, full_tree_metr = TRUE, filename = "", ... )
geo.phylo( x, tree, inv.R, edge.path, branch.length, n.descen, full_tree_metr = TRUE, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
inv.R |
SpatRaster. Inverse of range size. See |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
n.descen |
numeric. A Named numeric vector of number of descendants for
each branch. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
Community metrics calculated:
Phylogenetic diversity (Faith 1992)
Species Richness
Evolutionary distinctiveness by fair-proportion (Isaac et al. 2007)
Phylogenetic endemism (Rosauer et al. 2009)
Weighted endemism (Crisp et al. 2001, Williams et al. 1994)
SpatRaster with one layer for each metric
Neander Marcel Heming
Rosauer, D. A. N., Laffan, S. W., Crisp, M. D., Donnellan, S. C. and Cook, L. G. (2009). Phylogenetic endemism: a new approach for identifying geographical concentrations of evolutionary history. Molecular ecology, 18(19), 4061-4072.
Faith, D. P. (1992). Conservation evaluation and phylogenetic diversity. Biological conservation, 61(1), 1-10.
Williams, P.H., Humphries, C.J., Forey, P.L., Humphries, C.J. and VaneWright, R.I. (1994). Biodiversity, taxonomic relatedness, and endemism in conservation. In: Systematics and Conservation Evaluation (eds Forey PL, Humphries C.J., Vane-Wright RI), p. 438. Oxford University Press, Oxford.
Crisp, M., Laffan, S., Linder, H. and Monro, A. (2001). Endemism in the Australian flora. Journal of Biogeography, 28, 183–198.
Isaac, N. J., Turvey, S. T., Collen, B., Waterman, C. and Baillie, J. E. (2007). Mammals on the EDGE: conservation priorities based on threat and phylogeny. PLoS ONE 2, e296.
Laffan, S. W., Rosauer, D. F., Di Virgilio, G., Miller, J. T., González‐Orozco, C. E., Knerr, N., ... & Mishler, B. D. (2016). Range‐weighted metrics of species and phylogenetic turnover can better resolve biogeographic transition zones. Methods in Ecology and Evolution, 7(5), 580-588.
phylo.pres
, inv.range
,
rast.ed
, rast.pd
,
rast.we
, rast.pe
, rast.sr
,
geo.phylo.ses
,
library(terra) library(phyloraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster"))[[1:10]] tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x, tree) inv.R <- inv.range(data$x) t <- geo.phylo(data$x, inv.R = inv.R, edge.path = data$edge.path, branch.length = data$branch.length, n.descen = data$n.descendants) terra::plot(t)
library(terra) library(phyloraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster"))[[1:10]] tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x, tree) inv.R <- inv.range(data$x) t <- geo.phylo(data$x, inv.R = inv.R, edge.path = data$edge.path, branch.length = data$branch.length, n.descen = data$n.descendants) terra::plot(t)
Calculates the standardized effect size for phylogenetic community metrics. See Details for more information.
geo.phylo.ses( x, tree, inv.R, edge.path, branch.length, n.descen, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, cores = 1, filename = "", ... )
geo.phylo.ses( x, tree, inv.R, edge.path, branch.length, n.descen, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, cores = 1, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
inv.R |
SpatRaster. Inverse of range size. See |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
n.descen |
numeric. A Named numeric vector of number of descendants for
each branch. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
spat_alg |
A function with the algorithm implementing the desired
randomization method. It must work with SpatRaster objects. See examples.
Example of functions that work are: |
spat_alg_args |
List of arguments passed to the randomization method
chosen in 'spat_alg'. See |
aleats |
positive integer. A positive integer indicating how many times the calculation should be repeated. |
cores |
positive integer. If |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells.
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells.
SpatRaster. The function returns the observed value of the metric, the mean of the simulations calculated over n times, the standard deviation of the simulations, the standardized effect size (SES) for the metric, and the p-values.
Neander Marcel Heming
Gotelli, N. J. 2000. Null model analysis of species co-occurrence patterns. – Ecology 81: 2606–2621.
Heming, N. M., Mota, F. M. M. and Alves-Ferreira, G. 2023. SESraster: raster randomization for null hypothesis testing. https://CRAN.R-project.org/package=SESraster.
phylo.pres
,
inv.range
,
geo.phylo
,
rast.ed.ses
,
rast.pd.ses
,
rast.we.ses
,
rast.pe.ses
,
bootspat_str
,
bootspat_naive
,
bootspat_ff
,
SESraster
library(terra) library(phyloraster) require("SESraster") x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) tses <- geo.phylo.ses(x = x, tree = tree, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 2) terra::plot(tses)
library(terra) library(phyloraster) require("SESraster") x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) tses <- geo.phylo.ses(x = x, tree = tree, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 2) terra::plot(tses)
Get range size in square kilometers for all cells that are not NA, the inverse of range size and the inverse of range size multiplied by branch length for multiple species using a raster of presence-absence.
inv.range(x, filename = "", overwrite = FALSE, ...)
inv.range(x, filename = "", overwrite = FALSE, ...)
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
filename |
character. Output filename |
overwrite |
logical. If TRUE, filename is overwritten |
... |
additional arguments to be passed passed down from a calling function. |
SpatRaster and numeric
Neander Marcel Heming and Gabriela Alves-Ferreira
# calculating the inverse of range size x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) inv.range(x[[5]])
# calculating the inverse of range size x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) inv.range(x[[5]])
This function load a phylogenetic tree, a raster and a data.frame with presence-absence of 33 Australian tree frogs from Rosauer (2017). We also provide distribution shapefiles for ten species according to the IUCN.
load.data.rosauer()
load.data.rosauer()
data.frame, SpatRaster, SpatVector and phylo
Rosauer, 2017. Available on: Github
IUCN. 2022. The IUCN Red List of Threatened Species (spatial data). Version 2022-1. IUCN
Reorder a stack of rasters of species distribution to match the order of the tips of the tree, and get branch length and number of descendants for each species to calculate diversity metrics using phyloraster::geo.phylo(). The branch length and the number of descendants can be calculated based on the full tree or the raster based tree subset. The names must be the same in the phylogenetic tree and in the raster for the same species. For example, if you have the name "Leptodactylus_latrans" in the raster and "Leptodactylus latrans" in the tree, the function will not work. The same goes for uppercase and lowercase letters.
phylo.pres(x, tree, full_tree_metr = TRUE, ...)
phylo.pres(x, tree, full_tree_metr = TRUE, ...)
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. |
tree |
phylo. A dated tree. |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
... |
additional arguments to be passed passed down from a calling function. |
Returns a list containing a SpatRaster reordered according to the order that the species appear in the phylogenetic tree, a subtree containing only the species that are in the stack of rasters and finally two named numerical vectors containing the branch length and the number of descendants of each species.
Neander Marcel Heming and Gabriela Alves Ferreira
library(phyloraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) phylo.pres(x[[1:3]], tree, full_tree_metr = TRUE) # using the prunned tree phylo.pres(x[[1:3]], tree, full_tree_metr = FALSE)
library(phyloraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) phylo.pres(x[[1:3]], tree, full_tree_metr = TRUE) # using the prunned tree phylo.pres(x[[1:3]], tree, full_tree_metr = FALSE)
This function calculate range size in square meters (by default) for all cells that are not NA. The size of the cells is constant in degrees but not in square meters, which was considered in the method applied to calculate the area.
range_size(x, cellSz, unit = "m", ...)
range_size(x, cellSz, unit = "m", ...)
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
cellSz |
SpatRaster. A SpatRaster containing cellSize values.
See |
unit |
character. One of "m", "km", or "ha" |
... |
additional arguments to be passed passed down from a calling function. |
vector
Gabriela Alves Ferreira and Neander Marcel Heming
x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) range_size(x[[1:2]], cellSz <- terra::cellSize(x))
x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) range_size(x[[1:2]], cellSz <- terra::cellSize(x))
This function calculates evolutionary distinctiveness according to the fair-proportion index. The values represents the mean ED for species presents in each raster cell.
rast.ed( x, tree, edge.path, branch.length, n.descen, full_tree_metr = TRUE, filename = "", ... )
rast.ed( x, tree, edge.path, branch.length, n.descen, full_tree_metr = TRUE, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
n.descen |
numeric. A Named numeric vector of number of descendants for
each branch. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
SpatRaster
Neander Marcel Heming and Gabriela Alves-Ferreira
Isaac, N. J., Turvey, S. T., Collen, B., Waterman, C. and Baillie, J. E. (2007). Mammals on the EDGE: conservation priorities based on threat and phylogeny. PLoS ONE 2, e296.
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) # phylogenetic tree tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) ed <- rast.ed(data$x, edge.path = data$edge.path, branch.length = data$branch.length, n.descen = data$n.descen) plot(ed)
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) # phylogenetic tree tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) ed <- rast.ed(data$x, edge.path = data$edge.path, branch.length = data$branch.length, n.descen = data$n.descen) plot(ed)
Calculates the standardized effect size for evolutionary distinctiveness. See Details for more information.
rast.ed.ses( x, tree, edge.path, branch.length, n.descen, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
rast.ed.ses( x, tree, edge.path, branch.length, n.descen, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
n.descen |
numeric. A Named numeric vector of number of descendants for
each branch. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
spat_alg |
A function with the algorithm implementing the desired
randomization method. It must work with SpatRaster objects. See examples.
Example of functions that work are: |
spat_alg_args |
List of arguments passed to the randomization method
chosen in 'spat_alg'. See |
aleats |
positive integer. A positive integer indicating how many times the calculation should be repeated. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells.
SpatRaster. The function returns the observed value of the metric, the mean of the simulations calculated over n times, the standard deviation of the simulations, the standardized effect size (SES) for the metric, and the p-values.
Neander M. Heming and Gabriela Alves-Ferreira
Gotelli, N. J. 2000. Null model analysis of species co-occurrence patterns. – Ecology 81: 2606–2621.
Heming, N. M., Mota, F. M. M. and Alves-Ferreira, G. 2023. SESraster: raster randomization for null hypothesis testing. https://CRAN.R-project.org/package=SESraster.
phylo.pres
,
inv.range
,
geo.phylo.ses
,
rast.ed.ses
,
rast.pd.ses
,
rast.we.ses
,
rast.pe.ses
,
bootspat_str
,
bootspat_naive
,
bootspat_ff
,
SESraster
library(phyloraster) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) t <- rast.ed.ses(x[[1:10]], tree, aleats = 3) terra::plot(t)
library(phyloraster) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) t <- rast.ed.ses(x[[1:10]], tree, aleats = 3) terra::plot(t)
Calculate the sum of the branch length for species present in each cell of the raster.
rast.pd( x, tree, edge.path, branch.length, full_tree_metr = TRUE, filename = "", ... )
rast.pd( x, tree, edge.path, branch.length, full_tree_metr = TRUE, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
SpatRaster
Neander Marcel Heming and Gabriela Alves-Ferreira
Faith, D. P. (1992). Conservation evaluation and phylogenetic diversity. Biological conservation, 61(1), 1-10.
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) pd <- rast.pd(data$x, data$tree) plot(pd)
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) pd <- rast.pd(data$x, data$tree) plot(pd)
Calculates the standardized effect size for phylogenetic diversity. See Details for more information.
rast.pd.ses( x, tree, edge.path, branch.length, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
rast.pd.ses( x, tree, edge.path, branch.length, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
edge.path |
matrix. Matrix representing the paths through the tree from
root to each tip. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
spat_alg |
A function with the algorithm implementing the desired
randomization method. It must work with SpatRaster objects. See examples.
Example of functions that work are: |
spat_alg_args |
List of arguments passed to the randomization method
chosen in 'spat_alg'. See |
aleats |
positive integer. A positive integer indicating how many times the calculation should be repeated. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells.
SpatRaster. The function returns the observed value of the metric, the mean of the simulations calculated over n times, the standard deviation of the simulations, the standardized effect size (SES) for the metric, and the p-values.
Gabriela Alves-Ferreira and Neander Heming
Gotelli, N. J. 2000. Null model analysis of species co-occurrence patterns. – Ecology 81: 2606–2621.
Heming, N. M., Mota, F. M. M. and Alves-Ferreira, G. 2023. SESraster: raster randomization for null hypothesis testing. https://CRAN.R-project.org/package=SESraster.
phylo.pres
, inv.range
,
geo.phylo.ses
,
rast.ed.ses
, rast.pd.ses
,
rast.we.ses
, rast.pe.ses
,
bootspat_str
,
bootspat_naive
,
bootspat_ff
,
SESraster
library(terra) library(phyloraster) library(SESraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:10]], tree) t <- rast.pd.ses(data$x, edge.path = data$edge.path, branch.length = data$branch.length, aleats = 3) plot(t)
library(terra) library(phyloraster) library(SESraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:10]], tree) t <- rast.pd.ses(data$x, edge.path = data$edge.path, branch.length = data$branch.length, aleats = 3) plot(t)
Calculate the sum of the inverse of the range size multiplied by the branch length for the species present in raster data.
rast.pe( x, tree, inv.R, branch.length, full_tree_metr = TRUE, filename = "", ... )
rast.pe( x, tree, inv.R, branch.length, full_tree_metr = TRUE, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
inv.R |
SpatRaster. Inverse of range size. See |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
SpatRaster
Gabriela Alves-Ferreira and Neander Marcel Heming
Laffan, S. W., Rosauer, D. F., Di Virgilio, G., Miller, J. T., González‐Orozco, C. E., Knerr, N., ... & Mishler, B. D. (2016). Range‐weighted metrics of species and phylogenetic turnover can better resolve biogeographic transition zones. Methods in Ecology and Evolution, 7(5), 580-588.
Rosauer, D. A. N., Laffan, S. W., Crisp, M. D., Donnellan, S. C. and Cook, L. G. (2009). Phylogenetic endemism: a new approach for identifying geographical concentrations of evolutionary history. Molecular ecology, 18(19), 4061-4072.
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package = "phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package = "phyloraster")) pe <- rast.pe(x = x[[1:3]], tree) plot(pe)
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package = "phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package = "phyloraster")) pe <- rast.pe(x = x[[1:3]], tree) plot(pe)
Calculates the standardized effect size for phylogenetic endemism. See Details for more information.
rast.pe.ses( x, tree, branch.length, branch.length.alt, inv.R, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), metric = c("pe", "pe.alt", "rpe", "all")[4], aleats = 10, cores = 1, filename = "", overwrite = TRUE, ... )
rast.pe.ses( x, tree, branch.length, branch.length.alt, inv.R, full_tree_metr = TRUE, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), metric = c("pe", "pe.alt", "rpe", "all")[4], aleats = 10, cores = 1, filename = "", overwrite = TRUE, ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
tree |
phylo. A dated tree. |
branch.length |
numeric. A Named numeric vector of branch length for
each species. See |
branch.length.alt |
numeric. Branch length calculated by using an alternative phylogeny with non-zero branch lengths converted to a constant value (1) and rescaled so the sum of all branch lengths is 1. |
inv.R |
SpatRaster. Inverse of range size. See |
full_tree_metr |
logical. Whether edge.path, branch length and number of descendants should be calculated with the full (TRUE) or the prunned tree (FALSE). The default is TRUE. |
spat_alg |
A function with the algorithm implementing the desired
randomization method. It must work with SpatRaster objects. See examples.
Example of functions that work are: |
spat_alg_args |
List of arguments passed to the randomization method
chosen in 'spat_alg'. See |
metric |
character. Names of biodiversity metrics to calculate (pe, pe_alt, rpe, all). See details. |
aleats |
positive integer. A positive integer indicating how many times the calculation should be repeated. |
cores |
positive integer. If |
filename |
character. Output filename |
overwrite |
logical. If TRUE, filename is overwritten |
... |
additional arguments passed for terra::app |
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells. Biodiversity metrics available are:
pe: Phylogenetic endemism (Rosauer et al., 2009)
pe.alt: Alternate Phylogenetic endemism (Mishler et al., 2014)
rpe: Relative Phylogenetic endemism (Mishler et al., 2014)
all: Calculate all available metrics Alternate phylogenetic endemism (PE.alt, Mishler et al., 2014) is calculated using an alternate phylogeny with non-zero branch lengths converted to a constant value (here we use 1) and rescaled so the sum of all branch lengths is 1. Relative phylogenetic endemism (RPE, Mishler et al., 2014) is the ratio of phylogenetic endemism (PE, Rosauer et al., 2009) measured on the original tree versus PE measured on a alternate tree (PE.alt).
SpatRaster. The function returns the observed value of the metric, the mean of the simulations calculated over n times, the standard deviation of the simulations, the standardized effect size (SES) for the metric, and the p-values.
Gabriela Alves-Ferreira and Neander Heming
Gotelli, N. J. 2000. Null model analysis of species co-occurrence patterns. Ecology 81: 2606–2621.
Heming, N. M., Mota, F. M. M. and Alves-Ferreira, G. 2023. SESraster: raster randomization for null hypothesis testing. https://CRAN.R-project.org/package=SESraster.
Mishler, B. D., Knerr, N., González-Orozco, C. E., Thornhill, A. H., Laffan, S. W. and Miller, J. T. 2014. Phylogenetic measures of biodiversity and neo- and paleo-endemism in Australian Acacia. – Nat. Commun. 5: 4473.
Rosauer, D. A. N., Laffan, S. W., Crisp, M. D., Donnellan, S. C., & Cook, L. G. (2009). Phylogenetic endemism: a new approach for identifying geographical concentrations of evolutionary history. Molecular ecology, 18(19), 4061-4072.
phylo.pres
,
inv.range
,
geo.phylo.ses
,
rast.ed.ses
,
rast.pd.ses
,
rast.we.ses
,
rast.pe.ses
,
bootspat_str
,
bootspat_naive
,
bootspat_ff
,
SESraster
library(terra) library(phyloraster) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) t <- rast.pe.ses(x = data$x, data$tree, aleats = 99, metric = "all") plot(t)
library(terra) library(phyloraster) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) data <- phylo.pres(x[[1:3]], tree) t <- rast.pe.ses(x = data$x, data$tree, aleats = 99, metric = "all") plot(t)
Calculate the species richness for raster data.
rast.sr(x, filename = "", ...)
rast.sr(x, filename = "", ...)
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. |
filename |
character. Output filename. |
... |
additional arguments to be passed passed down from a calling function. |
SpatRaster
Gabriela Alves Ferreira and Neander Marcel Heming
x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) rse <- phyloraster::rast.sr(x) terra::plot(rse)
x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) rse <- phyloraster::rast.sr(x) terra::plot(rse)
Calculate the weighted endemism for species present in raster data.
rast.we(x, inv.R, filename = "", ...)
rast.we(x, inv.R, filename = "", ...)
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
inv.R |
SpatRaster. Inverse of range size. See |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
SpatRaster
Neander Marcel Heming and Gabriela Alves Ferreira
Laffan, S. W., Rosauer, D. F., Di Virgilio, G., Miller, J. T., González‐Orozco, C. E., Knerr, N., ... & Mishler, B. D. (2016). Range‐weighted metrics of species and phylogenetic turnover can better resolve biogeographic transition zones. Methods in Ecology and Evolution, 7(5), 580-588.
Williams, P.H., Humphries, C.J., Forey, P.L., Humphries, C.J., VaneWright, R.I. (1994). Biodiversity, taxonomic relatedness, and endemism in conservation. In: Systematics and Conservation Evaluation (eds Forey PL, Humphries CJ, Vane-Wright RI), p. 438. Oxford University Press, Oxford.
Crisp, M., Laffan, S., Linder, H., Monro, A. (2001). Endemism in the Australian flora. Journal of Biogeography, 28, 183–198.
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) inv.R <- inv.range(x) we <- rast.we(x, inv.R) plot(we)
library(terra) library(phyloraster) x <- rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) inv.R <- inv.range(x) we <- rast.we(x, inv.R) plot(we)
Calculates the standardized effect size for weighted endemism. See Details for more information.
rast.we.ses( x, inv.R, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
rast.we.ses( x, inv.R, spat_alg = "bootspat_str", spat_alg_args = list(rprob = NULL, rich = NULL, fr_prob = NULL), aleats = 10, filename = "", ... )
x |
SpatRaster. A SpatRaster containing presence-absence data (0 or 1) for a set of species. The layers (species) will be sorted according to the tree order. See the phylo.pres function. |
inv.R |
SpatRaster. Inverse of range size. See |
spat_alg |
A function with the algorithm implementing the desired
randomization method. It must work with SpatRaster objects. See examples.
Example of functions that work are: |
spat_alg_args |
List of arguments passed to the randomization method
chosen in 'spat_alg'. See |
aleats |
positive integer. A positive integer indicating how many times the calculation should be repeated. |
filename |
character. Output filename |
... |
additional arguments passed for terra::app |
The dependency ‘SESraster’ is used to calculate the null models. This package currently implements six algorithms to randomize binary species distribution with several levels of constraints: SIM1, SIM2, SIM3, SIM5, SIM6 and SIM9 (sensu Gotelli 2000). The methods implemented in ‘SESraster’ are based on how species (originally rows) and sites (originally columns) are treated (i.e. fixed, equiprobable, or proportional sums) (Gotelli 2000). By default, the ‘phyloraster’ uses the function bootspat_ str() from the ‘SESraster’ package to conduct the randomizations, but the user is free to choose any of the other methods mentioned above through the spat_alg argument in the *.ses() functions of the ‘phyloraster’ package. The bootspat_str() is equivalent to the SIM5 (proportional-fixed) method of Gotelli (2000), which partially relaxes the spatial structure of species distributions, but keeps the spatial structure of the observed richness pattern across cells.
SpatRaster. The function returns the observed value of the metric, the mean of the simulations calculated over n times, the standard deviation of the simulations, the standardized effect size (SES) for the metric, and the p-values.
Gabriela Alves-Ferreira and Neander Heming
Gotelli, N. J. 2000. Null model analysis of species co-occurrence patterns. – Ecology 81: 2606–2621.
Heming, N. M., Mota, F. M. M. and Alves-Ferreira, G. 2023. SESraster: raster randomization for null hypothesis testing. https://CRAN.R-project.org/package=SESraster.
phylo.pres
,
inv.range
,
geo.phylo.ses
,
rast.ed.ses
,
rast.pd.ses
,
rast.we.ses
,
rast.pe.ses
,
bootspat_str
,
bootspat_naive
,
bootspat_ff
,
SESraster
library(terra) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) t <- rast.we.ses(x[[1:10]], aleats = 3) plot(t)
library(terra) library(SESraster) x <- terra::rast(system.file("extdata", "rast.presab.tif", package="phyloraster")) t <- rast.we.ses(x[[1:10]], aleats = 3) plot(t)
The function will rasterize the shapefile using the parameters of y, a SpatRaster. When the argument y is provided, the resolution parameter is ignored. When the argument ymask is TRUE, y is used as a mask for x.
shp2rast( x, y = NULL, sps.col, ymask = FALSE, background = NA, touches = TRUE, resolution, values = 1, filename = NULL, ... )
shp2rast( x, y = NULL, sps.col, ymask = FALSE, background = NA, touches = TRUE, resolution, values = 1, filename = NULL, ... )
x |
SpatVector or a two-column matrix (point coordinates) |
y |
SpatRaster |
sps.col |
character. It should be a variable name in x. |
ymask |
logical. If TRUE, y will be used as a mask for x. |
background |
numeric. Value to put in the cells that are not covered by any of the features of |
touches |
logical. If |
resolution |
numeric vector of length 1 or 2 to set the spatial resolution (see |
values |
typically a numeric vector of length |
filename |
character. Output filename |
... |
additional arguments passed to |
SpatRaster
library(terra) shp <- terra::vect(system.file("extdata", "shps_iucn_spps_rosauer.shp", package="phyloraster")) # create a polygon to use as mask with an extent e <- terra::ext(113, 123, -43.64, -33.90) p <- terra::as.polygons(e, crs="") coun.crop <- terra::crop(p, terra::ext(shp)) coun.rast <- terra::rasterize(coun.crop, terra::rast(terra::ext(shp), resolution = 0.5)) plot(coun.rast, col = "green") # rasterizing with the mask of the polygon shp.t <- shp2rast(shp, y = coun.rast, sps.col = "BINOMIAL", ymask = TRUE, background = 0) plot(shp.t, col = c("grey", "green")) # rasterizing without using mask shp.t2 <- shp2rast(shp, sps.col = "BINOMIAL", ymask = FALSE, background = NA, resolution = 0.1) plot(shp.t2[[9]], col = c("grey", "green"))
library(terra) shp <- terra::vect(system.file("extdata", "shps_iucn_spps_rosauer.shp", package="phyloraster")) # create a polygon to use as mask with an extent e <- terra::ext(113, 123, -43.64, -33.90) p <- terra::as.polygons(e, crs="") coun.crop <- terra::crop(p, terra::ext(shp)) coun.rast <- terra::rasterize(coun.crop, terra::rast(terra::ext(shp), resolution = 0.5)) plot(coun.rast, col = "green") # rasterizing with the mask of the polygon shp.t <- shp2rast(shp, y = coun.rast, sps.col = "BINOMIAL", ymask = TRUE, background = 0) plot(shp.t, col = c("grey", "green")) # rasterizing without using mask shp.t2 <- shp2rast(shp, sps.col = "BINOMIAL", ymask = FALSE, background = NA, resolution = 0.1) plot(shp.t2[[9]], col = c("grey", "green"))
This function calculates evolutionary distinctiveness according to the fair-proportion index for each species.
species.ed(tree)
species.ed(tree)
tree |
phylo. A dated tree. |
data.frame
Neander Marcel Heming and Gabriela Alves-Ferreira
Isaac, N. J., Turvey, S. T., Collen, B., Waterman, C. and Baillie, J. E. (2007). Mammals on the EDGE: conservation priorities based on threat and phylogeny. PLoS ONE 2, e296.
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) plot(tree) ed <- species.ed(tree) ed
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) plot(tree) ed <- species.ed(tree) ed
Computation of species tip length using a phylogeny.
species.tip.length(tree = NULL, edge.info = NULL, ...)
species.tip.length(tree = NULL, edge.info = NULL, ...)
tree |
phylo. A dated tree. |
edge.info |
Object returned by |
... |
additional arguments to be passed passed down from a calling function. |
Calculates tip lengths for all species in a phylogeny
returns a numeric vector giving the length of species branch.
Neander M. Heming
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) species.tip.length(tree) library(ape) set.seed(1) tree <- rtree(n=40) plot(tree) species.tip.length(tree) edge.info <- tip.root.path(tree) species.tip.length(edge.info = edge.info)
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) species.tip.length(tree) library(ape) set.seed(1) tree <- rtree(n=40) plot(tree) species.tip.length(tree) edge.info <- tip.root.path(tree) species.tip.length(edge.info = edge.info)
Computation of tree edge lengths and node paths from root to each tip to calculate PD for a entire phylogeny (= sum of all edge or branch lengths)
tip.root.path(tree)
tip.root.path(tree)
tree |
phylo. A dated tree. |
Based on the algorithm FastXtreePhylo of Peter D. Wilson
returns a list with two components: matrix H1 representing the paths through the tree from root to each tip, and edge.length a numeric vector giving the length of each branch in the tree. Some matrix algebra and a summation of the resulting vector gives the whole-tree PD value.
Peter Wilson
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) fxtp <- tip.root.path(tree) H1 <- fxtp$H1 edge.length <- fxtp$edge.length # PD for the whole community pres <- rep(1, nrow(H1)) sum((crossprod(H1, pres)>0) * edge.length) # PD for a random subset of the community pres <- sample(c(1, 0), nrow(H1), TRUE) sum((crossprod(H1, pres)>0) * edge.length)
library(phyloraster) tree <- ape::read.tree(system.file("extdata", "tree.nex", package="phyloraster")) fxtp <- tip.root.path(tree) H1 <- fxtp$H1 edge.length <- fxtp$edge.length # PD for the whole community pres <- rep(1, nrow(H1)) sum((crossprod(H1, pres)>0) * edge.length) # PD for a random subset of the community pres <- sample(c(1, 0), nrow(H1), TRUE) sum((crossprod(H1, pres)>0) * edge.length)