Title: | Spatial Simulation and Scale-Dependent Analysis of Biodiversity Changes |
---|---|
Description: | Simulation, analysis and sampling of spatial biodiversity data (May, Gerstner, McGlinn, Xiao & Chase 2017) <doi:10.1111/2041-210x.12986>. In the simulation tools user define the numbers of species and individuals, the species abundance distribution and species aggregation. Functions for analysis include species rarefaction and accumulation curves, species-area relationships and the distance decay of similarity. |
Authors: | Felix May [aut, cre, cph] , Alban Sagouis [aut] |
Maintainer: | Felix May <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.3.1 |
Built: | 2024-11-20 05:06:38 UTC |
Source: | https://github.com/mobiodiv/mobsim |
Get local abundance distribution in rectangle bounded by x0, y0, x0 + xsize, y0 + ysize
abund_rect(x0, y0, xsize, ysize, comm)
abund_rect(x0, y0, xsize, ysize, comm)
x0 |
x-coordinate of lower left corner |
y0 |
y-coordinate of lower left corner |
xsize |
Size of the subplot in x-direction |
ysize |
Size of the subplot in y-direction |
comm |
|
Integer vector with local species abundances
Creates a spatial community object with defined extent and with coordinates and species identities of all individuals in the community.
community(x, y, spec_id, xrange = c(0, 1), yrange = c(0, 1))
community(x, y, spec_id, xrange = c(0, 1), yrange = c(0, 1))
x , y
|
Coordinates of individuals (numeric) |
spec_id |
Species names or IDs; can be integers, characters or factors |
xrange |
Extent of the community in x-direction (numeric vector of length 2) |
yrange |
Extent of the community in y-direction (numeric vector of length 2) |
Community object which includes three items:
census: data.frame with three columns: x, y, and species names for each individual
x_min_max: extent of the community in x-direction
y_min_max: extent of the community in y-direction
x <- runif(100) y <- runif(100) species_names <- rep(paste("species",1:10, sep = ""), each = 10) com1 <- community(x,y, species_names) plot(com1) summary(com1)
x <- runif(100) y <- runif(100) species_names <- rep(paste("species",1:10, sep = ""), each = 10) com1 <- community(x,y, species_names) plot(com1) summary(com1)
Get species abundance distribution from community object
community_to_sad(comm)
community_to_sad(comm)
comm |
Community object |
Object of class sad
, which contains a named integer vector
with species abundances
sim1 <- sim_poisson_community(s_pool = 200, n_sim = 20000, sad_type = "lnorm", sad_coef = list("cv_abund" = 2)) sad1 <- community_to_sad(sim1) plot(sad1, method = "rank") plot(sad1, method = "octave")
sim1 <- sim_poisson_community(s_pool = 200, n_sim = 20000, sad_type = "lnorm", sad_coef = list("cv_abund" = 2)) sad1 <- community_to_sad(sim1) plot(sad1, method = "rank") plot(sad1, method = "octave")
Estimate pairwise similarities of communities in subplots as function of distance
dist_decay( comm, prop_area = 0.005, n_samples = 20, method = "bray", binary = FALSE )
dist_decay( comm, prop_area = 0.005, n_samples = 20, method = "bray", binary = FALSE )
comm |
|
prop_area |
Subplot size as proportion of the total area |
n_samples |
Number of randomly located subplots |
method |
Choice of (dis)similarity index. See |
binary |
Perform presence/absence standardization before analysis?
See |
Object of class dist_decay
: a dataframe with distances between
subplot pairs and the respective similarity indices.
sim_com1 <- sim_thomas_community(100, 10000, sigma = 0.1, mother_points = 2) dd1 <- dist_decay(sim_com1, prop_area = 0.005, n_samples = 20) plot(dd1)
sim_com1 <- sim_thomas_community(100, 10000, sigma = 0.1, mother_points = 2) dd1 <- dist_decay(sim_com1, prop_area = 0.005, n_samples = 20) plot(dd1)
Estimate pairwise similarities of communities in quadrats as function of distance. The function allows the user to compute distance decay between the quadrats of his/her choice.
dist_decay_quadrats(samples, method = "bray", binary = FALSE)
dist_decay_quadrats(samples, method = "bray", binary = FALSE)
samples |
A list given by |
method |
Choice of (dis)similarity index. See |
binary |
Perform presence/absence standardization before analysis?
See |
Object of class dist_decay
: a dataframe with distances between
subplot pairs and the respective similarity indices.
sim_com1 <- sim_thomas_community(100, 10000, sigma = 0.1, mother_points = 2) oldpar<- par(mfrow=c(1,2)) samples <- sample_quadrats(sim_com1, avoid_overlap = TRUE, quadrat_area=.005, n_quadrats = 50, plot = TRUE) dd_quadrats <- dist_decay_quadrats(samples) plot(dd_quadrats) par(oldpar)
sim_com1 <- sim_thomas_community(100, 10000, sigma = 0.1, mother_points = 2) oldpar<- par(mfrow=c(1,2)) samples <- sample_quadrats(sim_com1, avoid_overlap = TRUE, quadrat_area=.005, n_quadrats = 50, plot = TRUE) dd_quadrats <- dist_decay_quadrats(samples) plot(dd_quadrats) par(oldpar)
Get mean and standard deviation of diversity indices in several equally sized subplots of a community
div_rand_rect(prop_area = 0.25, comm, n_rect = 100, exclude_zeros = FALSE)
div_rand_rect(prop_area = 0.25, comm, n_rect = 100, exclude_zeros = FALSE)
prop_area |
Size of subplots as proportion of the total area |
comm |
|
n_rect |
Number of randomly located subplots |
exclude_zeros |
Should subplots without individuals be excluded? (logical) |
Vector with mean and standard deviation of the following diversity indices:
Number of species
Number of endemics
Shannon index
Effective number of species (ENS) based on Shannon index
Simpson index
Effective number of species (ENS) based on Simpson index
See the documentation of div_rect
for detailed information on the
definition of the diversity indices.
sim1 <- sim_poisson_community(100,1000) div_rand_rect(prop_area = 0.1, comm = sim1)
sim1 <- sim_poisson_community(100,1000) div_rand_rect(prop_area = 0.1, comm = sim1)
Get diversity indices including species richness, no. of endemics, Shannon and Simpson diversity for one rectangle subplot in the community.
div_rect(x0, y0, xsize, ysize, comm)
div_rect(x0, y0, xsize, ysize, comm)
x0 |
x-coordinate of lower left corner |
y0 |
y-coordinate of lower left corner |
xsize |
Size of the subplot in x-direction |
ysize |
Size of the subplot in y-direction |
comm |
|
The effective number of species is defined as the number of equally abundant species that produce the same value of a certain diversity index as an observed community (Jost 2006). According to Chao et al. 2014 and Chiu et al. 20 ENS_shannon can be interpreted as the number of common species and ENS_simpson as the number of dominant species in a community.
Named vector with six diversity indices
n_species: Number of species
n_endemics: Number of endemics
shannon: Shannon index index defined as ,
where
is the relative abundance of species i:
ens_shannon: Effective number of species (ENS) based on the Shannon index exp(H)
simpson: Simpson index index (= probability of interspecific encounter PIE)
defined as
ens_simpson: Effective number of species (ENS) based on the Simpson index
Jost 2006. Entropy and diversity. Oikos, 113, 363-375.
Chao et al. 2014. Rarefaction and extrapolation with Hill numbers: a framework for sampling and estimation in species diversity studies. Ecological Monographs, 84, 45-67.
Hsieh et al. 2016. iNEXT: an R package for rarefaction and extrapolation of species diversity (Hill numbers). Methods Ecol Evol, 7, 1451-1456.
sim1 <- sim_poisson_community(100,1000) div_rect(0, 0, 0.3, 0.3, sim1)
sim1 <- sim_poisson_community(100,1000) div_rect(0, 0, 0.3, 0.3, sim1)
Estimate diversity indices in subplots of different sizes. This includes the well-known species-area and endemics-area relationships.
divar( comm, prop_area = seq(0.1, 1, by = 0.1), n_samples = 100, exclude_zeros = TRUE )
divar( comm, prop_area = seq(0.1, 1, by = 0.1), n_samples = 100, exclude_zeros = TRUE )
comm |
|
prop_area |
Subplot sizes as proportion of the total area (numeric) |
n_samples |
Number of randomly located subplots per subplot size (single integer) |
exclude_zeros |
Should subplots without individuals be excluded? (logical) |
Dataframe with the proportional area of the subplots and mean and standard deviation of the following diversity indices:
Number of species
Number of endemics
Shannon index
Effective number of species (ENS) based on Shannon index
Simpson index
Effective number of species (ENS) based on Simpson index
See the documentation of div_rect
for detailed information on the
definition of the diversity indices.
sim1 <- sim_thomas_community(100, 1000) divar1 <- divar(sim1, prop_area = seq(0.01, 1.0, length = 20)) plot(divar1)
sim1 <- sim_thomas_community(100, 1000) divar1 <- divar(sim1, prop_area = seq(0.01, 1.0, length = 20)) plot(divar1)
Plot positions and species identities of all individuals in a community object.
## S3 method for class 'community' plot(x, ..., col = NULL, pch = NULL)
## S3 method for class 'community' plot(x, ..., col = NULL, pch = NULL)
x |
Community object |
... |
Other parameters to |
col |
Colour vector to mark species identities |
pch |
Plotting character to mark species identities. pch 16 is advised for large datasets |
This function is called for its side effects and has no return value.
sim1 <- sim_thomas_community(30, 500) plot(sim1)
sim1 <- sim_thomas_community(30, 500) plot(sim1)
Plot distance decay of similarity
## S3 method for class 'dist_decay' plot(x, ...)
## S3 method for class 'dist_decay' plot(x, ...)
x |
Dataframe generated by |
... |
Additional graphical parameters used in |
The function plots the similarity indices between all pairs of
subplots as function of distance. To indicate the relationship a
stats::loess
smoother is added to the plot.
This function is called for its side effects and has no return value.
sim_com1 <- sim_thomas_community(100, 10000) dd1 <- dist_decay(sim_com1) plot(dd1)
sim_com1 <- sim_thomas_community(100, 10000) dd1 <- dist_decay(sim_com1) plot(dd1)
—————————————————————————– Plot diversity-area relationships
## S3 method for class 'divar' plot(x, ...)
## S3 method for class 'divar' plot(x, ...)
x |
Dataframe generated by the function |
... |
Additional graphical parameters used in |
This function is called for its side effects and has no return value.
Plot species abundance distributions
## S3 method for class 'sad' plot(x, ..., method = c("octave", "rank"))
## S3 method for class 'sad' plot(x, ..., method = c("octave", "rank"))
x |
Vector with species abundances (integer vector) |
... |
Additional graphical parameters used in |
method |
Plotting method, partial match to |
With method = "octave"
a histogram showing the number
species in several abundance classes is generated. The abundance class
are a simplified version of the "octaves" suggested by Preston (1948), which
are based on log2-binning. The first abundance class includes species
with 1 individual, the second with 2, the third with 3-4, the fourth with 5-8, etc.
With method = "rank"
rank-abundance curve is generated with
species abundance rank on the x-axis (descending) and species abundance on
the y-axis (Hubbell 2001).
This function is called for its side effects and has no return value.
Preston 1948. The Commonness, and rarity, of species. Ecology 29(3):254-283.
Hubbell 2001. The unified neutral theory of biodiversity and biogeography. Princeton University Press.
abund1 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 1)) plot(abund1, method = "octave") plot(abund1, method = "rank")
abund1 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 1)) plot(abund1, method = "octave") plot(abund1, method = "rank")
Plot species sampling curves
## S3 method for class 'spec_sample_curve' plot(x, ...)
## S3 method for class 'spec_sample_curve' plot(x, ...)
x |
Species sampling curve generated by |
... |
Additional graphical parameters used in |
This function is called for its side effects and has no return value.
sim_com1 <- sim_thomas_community(s_pool = 100, n_sim = 1000) sac1 <- spec_sample_curve(sim_com1, method = c("rare","acc")) plot(sac1)
sim_com1 <- sim_thomas_community(s_pool = 100, n_sim = 1000) sac1 <- spec_sample_curve(sim_com1, method = c("rare","acc")) plot(sac1)
Expected species richness as a function of sample size
rare_curve(abund_vec)
rare_curve(abund_vec)
abund_vec |
Species abundance distribution of the community (integer vector) |
This function essentially evaluates spec_sample
for
sample sizes from 1 to sum(abund_vec)
. It is similar to the function
vegan:rarecurve
in the R package vegan
.
Numeric Vector with expected species richness in samples of 1, 2, 3 ... n individuals
Gotelli & Colwell 2001. Quantifying biodiversity: procedures and pitfalls in the measurement and comparison of species richness. Ecology Letters 4, 379–391.
sad1 <- sim_sad(100, 2000, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1)) rc1 <- rare_curve(sad1) plot(rc1, type = "l", xlab = "Sample size", ylab = "Expected species richness")
sad1 <- sim_sad(100, 2000, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1)) rc1 <- rare_curve(sad1) plot(rc1, type = "l", xlab = "Sample size", ylab = "Expected species richness")
Usually used internally inside sim_thomas_coords
This function randomly draws points (individuals) around one or several mother points using Rcpp.
The function is an efficient re-implementation of the rThomas function from the spatstat package.
n_points |
The total number of points (individuals). |
n_mother_points |
Number of mother points (= cluster centres). |
xmother |
Vector of |
ymother |
Vector of |
sigma |
Mean displacement (along each coordinate axes) of a point from its mother point (= cluster centre). |
xmin |
Left limit, |
xmax |
Right limit, |
ymin |
Bottom limit, |
ymax |
Top limit, |
A dataframe with x and y coordinates.
Felix May, Alban Sagouis
This function allows to sample quadratic subplots from a spatially-explicit
community. The output format are a sites x species abundance table and a
sites x xy-coordinates table. The sites x species abundance is
a classical data format used in community ecology. The table generated
can be for instance be further analysed with the package vegan
.
sample_quadrats( comm, n_quadrats = 20, quadrat_area = 0.01, plot = TRUE, method = "random", avoid_overlap = TRUE, x0 = 0, y0 = 0, delta_x = 0.1, delta_y = 0.1, seed = NULL )
sample_quadrats( comm, n_quadrats = 20, quadrat_area = 0.01, plot = TRUE, method = "random", avoid_overlap = TRUE, x0 = 0, y0 = 0, delta_x = 0.1, delta_y = 0.1, seed = NULL )
comm |
Community object from which the samples are generated |
n_quadrats |
(integer) Number of sampling quadrats |
quadrat_area |
(numeric) Area of the sampling quadrats |
plot |
(logical) Should the sampling design be plotted? default to TRUE. |
method |
(character) Available methods are |
avoid_overlap |
(logical) For the random sampling try to generate a design without overlap of quadrats . Default is TRUE. |
x0 , y0
|
(numeric value) Lower left corner of the first quadrat in transect and grid sampling |
delta_x |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in x-direction (the distance between the left sides is measured) |
delta_y |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in y-direction (the distance between the lower sides is measured) |
seed |
(integer) Any integer passed to |
A list with two items, spec_dat
and xy_dat
.
spec_dat
is a data.frame with sampling quadrats in rows and species abundances
in columns, and xy_dat
is a data.frame with sampling quadrats in rows
and the xy-coordinates of the quadrats (lower left corner) in columns.
library(vegan) sim_com1 <- sim_poisson_community(100, 10000) comm_mat1 <- sample_quadrats(sim_com1, n_quadrats = 100, quadrat_area = 0.002, method = "grid") specnumber(comm_mat1$spec_dat) diversity(comm_mat1$spec_dat, index = "shannon")
library(vegan) sim_com1 <- sim_poisson_community(100, 10000) comm_mat1 <- sample_quadrats(sim_com1, n_quadrats = 100, quadrat_area = 0.002, method = "grid") specnumber(comm_mat1$spec_dat) diversity(comm_mat1$spec_dat, index = "shannon")
Creates square quadrats aligned on a regular grid
sampling_grids( n_quadrats, xmin, xmax, ymin, ymax, x0, y0, delta_x, delta_y, quadrat_size )
sampling_grids( n_quadrats, xmin, xmax, ymin, ymax, x0, y0, delta_x, delta_y, quadrat_size )
n_quadrats |
(integer) Number of sampling quadrats |
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
x0 , y0
|
(numeric value) Lower left corner of the first quadrat in transect and grid sampling |
delta_x |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in x-direction (the distance between the left sides is measured) |
delta_y |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in y-direction (the distance between the lower sides is measured) |
quadrat_size |
(numeric) width of the quadrats. |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrats.
Creates one square quadrat randomly located in the landscape
sampling_one_quadrat(xmin, xmax, ymin, ymax, seed = NULL)
sampling_one_quadrat(xmin, xmax, ymin, ymax, seed = NULL)
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
seed |
(integer) Any integer passed to |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrat.
This function works without having the spatstat.random
package install.
sampling_random_bruteforce( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
sampling_random_bruteforce( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
n_quadrats |
Number of sampling quadrats |
min_dist |
(numeric) minimal distance between two points to avoid overlap. Equal to the length of a quadrat diagonal |
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
seed |
(integer) Any integer passed to |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrats.
Creates coordinates (lower left corner of a quadrat) randomly distributed that may overlap each other
sampling_random_overlap( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
sampling_random_overlap( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
n_quadrats |
Number of sampling quadrats |
min_dist |
(numeric) minimal distance between two points to avoid overlap. Equal to the length of a quadrat diagonal |
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
seed |
(integer) Any integer passed to |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrats.
Efficient algorithm from package spatstat.random
is used.
Produces similar results as sampling_random_bruteforce
.
sampling_random_spatstat( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
sampling_random_spatstat( n_quadrats, min_dist, xmin, xmax, ymin, ymax, seed = NULL )
n_quadrats |
Number of sampling quadrats |
min_dist |
(numeric) minimal distance between two points to avoid overlap. Equal to the length of a quadrat diagonal |
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
seed |
(integer) Any integer passed to |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrats.
Creates square quadrats aligned along a transect
sampling_transects( n_quadrats, xmin, xmax, ymin, ymax, x0, y0, delta_x, delta_y, quadrat_size )
sampling_transects( n_quadrats, xmin, xmax, ymin, ymax, x0, y0, delta_x, delta_y, quadrat_size )
n_quadrats |
(integer) Number of sampling quadrats |
xmin |
(numeric) minimum possible value on the x axis a quadrat can cover. |
xmax |
(numeric) maximum possible value on the x axis a quadrat can cover. |
ymin |
(numeric) minimum possible value on the y axis a quadrat can cover. |
ymax |
(numeric) maximum possible value on the y axis a quadrat can cover. |
x0 , y0
|
(numeric value) Lower left corner of the first quadrat in transect and grid sampling |
delta_x |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in x-direction (the distance between the left sides is measured) |
delta_y |
(numeric value) Distance between consecutive quadrats in transect and grid sampling in y-direction (the distance between the lower sides is measured) |
quadrat_size |
(numeric) width of the quadrats. |
a data.frame with 2 columns x and y giving the coordinates of the lower left corner of the square quadrats.
This function simulates a community with a certain abundance distribution and
and random spatial coordinates. This function consecutively calls
sim_sad
and sim_poisson_coords
sim_poisson_community( s_pool, n_sim, sad_type = "lnorm", sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
sim_poisson_community( s_pool, n_sim, sad_type = "lnorm", sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
s_pool |
Number of species in the pool (integer) |
n_sim |
Number of individuals in the simulated community (integer) |
sad_type |
Root name of the species abundance distribution model of the
species pool (character) - e.g., "lnorm" for the lognormal distribution
( See the table in Details below, or |
sad_coef |
List with named arguments to be passed to the distribution
function defined by the argument In Please note that the parameters mu and sigma are not equal to the mean and standard deviation of the log-normal distribution. |
fix_s_sim |
Should the simulation constrain the number of species in the simulated local community? (logical) |
xrange |
Extent of the community in x-direction (numeric vector of length 2) |
yrange |
Extent of the community in y-direction (numeric vector of length 2) |
seed |
Integer. Any integer passed to |
A community object as defined by community
.
Felix May
com1 <- sim_poisson_community(s_pool = 20, n_sim = 500, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1)) plot(com1)
com1 <- sim_poisson_community(s_pool = 20, n_sim = 500, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1)) plot(com1)
Add random spatial positions to a species abundance distribution.
sim_poisson_coords(abund_vec, xrange = c(0, 1), yrange = c(0, 1), seed = NULL)
sim_poisson_coords(abund_vec, xrange = c(0, 1), yrange = c(0, 1), seed = NULL)
abund_vec |
Species abundance vector (integer) |
xrange |
Extent of the community in x-direction (numeric vector of length 2) |
yrange |
Extent of the community in y-direction (numeric vector of length 2) |
seed |
Integer. Any integer passed to |
A community object as defined by community
.
Felix May
abund <- sim_sad(s_pool = 100, n_sim = 1000) sim_com1 <- sim_poisson_coords(abund) plot(sim_com1) summary(sim_com1)
abund <- sim_sad(s_pool = 100, n_sim = 1000) sim_com1 <- sim_poisson_coords(abund) plot(sim_com1) summary(sim_com1)
Simulate species abundance distribution (SAD) of a local community with user-defined number of species and relative abundance distribution in the pool, and user-defined number of individuals in the simulated local community.
sim_sad( s_pool = NULL, n_sim = NULL, sad_type = c("lnorm", "bs", "gamma", "geom", "ls", "mzsm", "nbinom", "pareto", "poilog", "power", "powbend", "weibull"), sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, drop_zeros = TRUE, seed = NULL )
sim_sad( s_pool = NULL, n_sim = NULL, sad_type = c("lnorm", "bs", "gamma", "geom", "ls", "mzsm", "nbinom", "pareto", "poilog", "power", "powbend", "weibull"), sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, drop_zeros = TRUE, seed = NULL )
s_pool |
Number of species in the pool (integer) |
n_sim |
Number of individuals in the simulated community (integer) |
sad_type |
Root name of the species abundance distribution model of the
species pool (character) - e.g., "lnorm" for the lognormal distribution
( See the table in Details below, or |
sad_coef |
List with named arguments to be passed to the distribution
function defined by the argument In Please note that the parameters mu and sigma are not equal to the mean and standard deviation of the log-normal distribution. |
fix_s_sim |
Should the simulation constrain the number of species in the simulated local community? (logical) |
drop_zeros |
Should the function remove species with abundance zero from the output? (logical) |
seed |
Integer. Any integer passed to |
The function sim_sad
was built using code of the function
sads::rsad
from the R package sads
. However, in
contrast to sads::rsad
, the function sim_sad
allows to
define the number of individuals in the simulated local community. This is
implemented by converting the abundance distribution simulated based on
sads::rsad
into a relative abundance distribution. This
relative abundance distribution is considered as the species pool for the
local community. In a second step the required no. of individuals (n_sim)
is sampled (with replacement) from this relative abundance distribution.
Please note that this might effect the interpretation of the parameters of
the underlying statistical distribution, e.g. the mean abundance will always
be n_sim/s_pool
irrespective of the settings of sad_coef
.
When fix_s_sim = FALSE
the species number in the local
community might deviate from s_pool
due to stochastic sampling. When
fix_s_sim = TRUE
the local number of species will equal
s_pool
, but this constraint can result in systematic biases from the
theoretical distribution parameters. Generally, with fix_s_sim = TRUE
additional very rare species will be added to the community, while the abundance
of the most common ones is reduced to keep the defined number of individuals.
Here is an overview of all available models (sad_type
) and their
respective coefficients (sad_coef
). Further information is provided
by the documentation of the specific functions that can be accesses by the
links. Please note that the coefficient cv_abund
for the log-normal
and Poisson log-normal model are only available within mobsim
.
SAD function | Distribution name | coef #1 | coef #2 | coef #3 |
sads::rbs |
Mac-Arthur's brokenstick | N | S | |
stats:rgamma |
Gamma distribution | shape | rate | scale |
rgeom |
Geometric distribution | prob | ||
rlnorm |
Log-normal distributions | meanlog | sdlog | cv_abund |
rls |
Fisher's log-series distribution | N | alpha | |
sads::rmzsm |
Metacommunity zero-sum multinomial | J | theta | |
stats::rnbinom |
Negative binomial distribution | size | prob | mu |
sads::rpareto |
Pareto distribution | shape | scale | |
sads::rpoilog |
Poisson-lognormal distribution | mu | sigma | cv_abund |
sads::rpower |
Power discrete distributions | s | ||
sads::rpowbend |
Puyeo's Power-bend discrete distribution | s | omega | |
stats::rweibull |
Weibull distribution | shape | scale | |
Object of class sad
, which contains a named integer vector
with species abundances
Felix May
#Simulate log-normal species abundance distribution sad_lnorm1 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("meanlog" = 5, "sdlog" = 0.5)) plot(sad_lnorm1, method = "octave") plot(sad_lnorm1, method = "rank") # Alternative parameterization of the log-normal distribution sad_lnorm2 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 0.5)) plot(sad_lnorm2, method = "octave") # Fix species richness in the simulation by adding rare species sad_lnorm3a <- sim_sad(s_pool = 500, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 5), fix_s_sim = TRUE) sad_lnorm3b <- sim_sad(s_pool = 500, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 5)) plot(sad_lnorm3a, method = "rank") points(1:length(sad_lnorm3b), sad_lnorm3b, type = "b", col = 2) legend("topright", c("fix_s_sim = TRUE","fix_s_sim = FALSE"), col = 1:2, pch = 1) # Different important SAD models # Fisher's log-series sad_logseries <- sim_sad(s_pool = NULL, n_sim = NULL, sad_type = "ls", sad_coef = list("N" = 1e5, "alpha" = 20)) # Poisson log-normal sad_poilog <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "poilog", sad_coef = list("mu" = 5, "sig" = 0.5)) # Mac-Arthur's broken stick sad_broken_stick <- sim_sad(s_pool = NULL, n_sim = NULL, sad_type = "bs", sad_coef = list("N" = 1e5, "S" = 100)) # Plot all SADs together as rank-abundance curves plot(sad_logseries, method = "rank") lines(1:length(sad_lnorm2), sad_lnorm2, type = "b", col = 2) lines(1:length(sad_poilog), sad_poilog, type = "b", col = 3) lines(1:length(sad_broken_stick), sad_broken_stick, type = "b", col = 4) legend("topright", c("Log-series","Log-normal","Poisson log-normal","Broken stick"), col = 1:4, pch = 1)
#Simulate log-normal species abundance distribution sad_lnorm1 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("meanlog" = 5, "sdlog" = 0.5)) plot(sad_lnorm1, method = "octave") plot(sad_lnorm1, method = "rank") # Alternative parameterization of the log-normal distribution sad_lnorm2 <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 0.5)) plot(sad_lnorm2, method = "octave") # Fix species richness in the simulation by adding rare species sad_lnorm3a <- sim_sad(s_pool = 500, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 5), fix_s_sim = TRUE) sad_lnorm3b <- sim_sad(s_pool = 500, n_sim = 10000, sad_type = "lnorm", sad_coef = list("cv_abund" = 5)) plot(sad_lnorm3a, method = "rank") points(1:length(sad_lnorm3b), sad_lnorm3b, type = "b", col = 2) legend("topright", c("fix_s_sim = TRUE","fix_s_sim = FALSE"), col = 1:2, pch = 1) # Different important SAD models # Fisher's log-series sad_logseries <- sim_sad(s_pool = NULL, n_sim = NULL, sad_type = "ls", sad_coef = list("N" = 1e5, "alpha" = 20)) # Poisson log-normal sad_poilog <- sim_sad(s_pool = 100, n_sim = 10000, sad_type = "poilog", sad_coef = list("mu" = 5, "sig" = 0.5)) # Mac-Arthur's broken stick sad_broken_stick <- sim_sad(s_pool = NULL, n_sim = NULL, sad_type = "bs", sad_coef = list("N" = 1e5, "S" = 100)) # Plot all SADs together as rank-abundance curves plot(sad_logseries, method = "rank") lines(1:length(sad_lnorm2), sad_lnorm2, type = "b", col = 2) lines(1:length(sad_poilog), sad_poilog, type = "b", col = 3) lines(1:length(sad_broken_stick), sad_broken_stick, type = "b", col = 4) legend("topright", c("Log-series","Log-normal","Poisson log-normal","Broken stick"), col = 1:4, pch = 1)
This function simulates a community with a certain abundance distribution and with intraspecific aggregation, i.e. individuals of the same species are distributed in clusters.
sim_thomas_community( s_pool, n_sim, sad_type = "lnorm", sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, sigma = 0.02, cluster_points = NA, mother_points = NA, xmother = NA, ymother = NA, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
sim_thomas_community( s_pool, n_sim, sad_type = "lnorm", sad_coef = list(cv_abund = 1), fix_s_sim = FALSE, sigma = 0.02, cluster_points = NA, mother_points = NA, xmother = NA, ymother = NA, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
s_pool |
Number of species in the pool (integer) |
n_sim |
Number of individuals in the simulated community (integer) |
sad_type |
Root name of the species abundance distribution model of the
species pool (character) - e.g., "lnorm" for the lognormal distribution
( See the table in Details below, or |
sad_coef |
List with named arguments to be passed to the distribution
function defined by the argument In Please note that the parameters mu and sigma are not equal to the mean and standard deviation of the log-normal distribution. |
fix_s_sim |
Should the simulation constrain the number of species in the simulated local community? (logical) |
sigma |
Mean displacement (along each coordinate axes) of a point from
its mother point (= cluster centre). |
cluster_points |
Mean number of points per cluster. If this is
a single value, species have the same average number of points per cluster.
If this is a vector of the same length as |
mother_points |
Number of mother points (= cluster centres).
If this is a single value, all species have the same number of clusters.
For example |
xmother |
List of length equal to the number of species. Each list element is a vector of x coordinates for every mother points. If one element is NA, the the corresponding species is not clustered. |
ymother |
List of length equal to the number of species. Each list element is a vector of y coordinates for every mother points. If one element is NA, the the corresponding species is not clustered. |
xrange |
Extent of the community in x-direction. If this a numeric vector
of length 2, all species share the same range. To specify different x ranges for
all species, |
yrange |
Extent of the community in y-direction. If this a numeric vector
of length 2, all species share the same range. To specify different y ranges for
all species, |
seed |
Integer. Any integer passed to |
This function consecutively calls sim_sad
and
sim_thomas_coords
See the documentations of sim_sad
and
sim_thomas_coords
for details.
A community object as defined by community
Felix May
com1 <- sim_thomas_community(s_pool = 20, n_sim = 500, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1), sigma = 0.01) plot(com1)
com1 <- sim_thomas_community(s_pool = 20, n_sim = 500, sad_type = "lnorm", sad_coef = list("meanlog" = 2, "sdlog" = 1), sigma = 0.01) plot(com1)
Add clumped (aggregated) positions to a species abundance distribution. Clumping is simulated using a Thomas cluster process, also known as Poisson cluster process (Morlon et al. 2008, Wiegand & Moloney 2014)
sim_thomas_coords( abund_vec, sigma = 0.02, mother_points = NA, xmother = NA, ymother = NA, cluster_points = NA, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
sim_thomas_coords( abund_vec, sigma = 0.02, mother_points = NA, xmother = NA, ymother = NA, cluster_points = NA, xrange = c(0, 1), yrange = c(0, 1), seed = NULL )
abund_vec |
Species abundance vector (integer) |
sigma |
Mean displacement (along each coordinate axes) of a point from
its mother point (= cluster centre). |
mother_points |
Number of mother points (= cluster centres).
If this is a single value, all species have the same number of clusters.
For example |
xmother |
List of length equal to the number of species. Each list element is a vector of x coordinates for every mother points. If one element is NA, the the corresponding species is not clustered. |
ymother |
List of length equal to the number of species. Each list element is a vector of y coordinates for every mother points. If one element is NA, the the corresponding species is not clustered. |
cluster_points |
Mean number of points per cluster. If this is
a single value, species have the same average number of points per cluster.
If this is a vector of the same length as |
xrange |
Extent of the community in x-direction. If this a numeric vector
of length 2, all species share the same range. To specify different x ranges for
all species, |
yrange |
Extent of the community in y-direction. If this a numeric vector
of length 2, all species share the same range. To specify different y ranges for
all species, |
seed |
Integer. Any integer passed to |
To generate a Thomas cluster process of a single species this
function uses a C++ re-implementation of the function
rThomas
in the package
spatstat.random.
There is an inherent link between the parameters abund_vec
,
mother_points
, and cluster_points
. For every species the
abundance has to be equal to the number of clusters
(mother_points
) times the number of points per cluster
(cluster_points
).
Accordingly, if one of the parameters is provided, the other one is directly
calculated from the abundance. Values for mother_points
override values
for cluster_points
. If none of the parameters is specified, it is assumed
that for every species there is a similar number of clusters and of points
per cluster.
In this case rare species have few clusters with few points per cluster, while abundant species have many clusters with many points per cluster.
A community object as defined by community
.
Felix May, Alban Sagouis
Morlon et al. 2008. A general framework for the distance-decay of similarity in ecological communities. Ecology Letters 11, 904-917.
Wiegand and Moloney 2014. Handbook of Spatial Point-Pattern Analysis in Ecology. CRC Press
abund <- c(10,20,50,100) sim1 <- sim_thomas_coords(abund, sigma = 0.02) plot(sim1) # Simulate species "ranges" sim2 <- sim_thomas_coords(abund, sigma = 0.02, mother_points = 1) plot(sim2) # Equal numbers of points per cluster sim3 <- sim_thomas_coords(abund, sigma = 0.02, cluster_points = 5) plot(sim3) # With large sigma the distribution will be essentially random (see Details) sim4 <- sim_thomas_coords(abund, sigma = 10) plot(sim4) # Some random and some clustered species with different numbers of mother points. mother_points <- sample(0:3, length(abund), replace = TRUE) sim5 <- sim_thomas_coords(abund, mother_points = mother_points, sigma=0.01) plot(sim5) # Specifying mother point coordinates or no-clustering (\code{NA}). mother_points <- sample(1:3, length(abund), replace = TRUE) xmother <- lapply(1:length(abund), function(i) runif(mother_points[i], 0, 1)) ymother <- lapply(1:length(abund), function(i) runif(mother_points[i], 0, 1)) xmother[[1]] <- NA ymother[[1]] <- NA sim6 <- sim_thomas_coords(abund, xmother=xmother, ymother=ymother, sigma=0.01) plot(sim6) # Species having different ranges. xrange <- data.frame(t(sapply(1:length(abund), function(i) sort(runif(2, 0, 1))))) yrange <- data.frame(t(sapply(1:length(abund), function(i) sort(runif(2, 0, 1))))) sim7 <- sim_thomas_coords(abund, mother_points=1, sigma=1, xrange=xrange, yrange=yrange) plot(sim7)
abund <- c(10,20,50,100) sim1 <- sim_thomas_coords(abund, sigma = 0.02) plot(sim1) # Simulate species "ranges" sim2 <- sim_thomas_coords(abund, sigma = 0.02, mother_points = 1) plot(sim2) # Equal numbers of points per cluster sim3 <- sim_thomas_coords(abund, sigma = 0.02, cluster_points = 5) plot(sim3) # With large sigma the distribution will be essentially random (see Details) sim4 <- sim_thomas_coords(abund, sigma = 10) plot(sim4) # Some random and some clustered species with different numbers of mother points. mother_points <- sample(0:3, length(abund), replace = TRUE) sim5 <- sim_thomas_coords(abund, mother_points = mother_points, sigma=0.01) plot(sim5) # Specifying mother point coordinates or no-clustering (\code{NA}). mother_points <- sample(1:3, length(abund), replace = TRUE) xmother <- lapply(1:length(abund), function(i) runif(mother_points[i], 0, 1)) ymother <- lapply(1:length(abund), function(i) runif(mother_points[i], 0, 1)) xmother[[1]] <- NA ymother[[1]] <- NA sim6 <- sim_thomas_coords(abund, xmother=xmother, ymother=ymother, sigma=0.01) plot(sim6) # Species having different ranges. xrange <- data.frame(t(sapply(1:length(abund), function(i) sort(runif(2, 0, 1))))) yrange <- data.frame(t(sapply(1:length(abund), function(i) sort(runif(2, 0, 1))))) sim7 <- sim_thomas_coords(abund, mother_points=1, sigma=1, xrange=xrange, yrange=yrange) plot(sim7)
Expected species richness in a random sample of fixed size.
spec_sample(abund_vec, n)
spec_sample(abund_vec, n)
abund_vec |
Species abundance distribution of the community (integer vector) |
n |
Sample size in terms of number of individuals (integer) |
The expected number of species is calculated after Hurlbert 1971, Equation 3.
spec_sample
is similar to the function rarefy
in the R package vegan
.
Expected number of species in a sample of n individuals
Hurlbert, S.H. 1971. The nonconcept of species diversity: a critique and + alternative parameters. Ecology 52, 577-586.
sad1 <- sim_sad(100, 1000) spec_sample(abund_vec = sad1, n = 20)
sad1 <- sim_sad(100, 1000) spec_sample(abund_vec = sad1, n = 20)
Expected species richness as function of sample size (no. of individuals), when individuals are sampled randomly (rarefaction) or when nearest-neighbours are samples (accumulation).
spec_sample_curve(comm, method = c("accumulation", "rarefaction"))
spec_sample_curve(comm, method = c("accumulation", "rarefaction"))
comm |
Community object |
method |
Partial match to |
Non-spatial sampling corresponds to the species rarefaction curve, which only
depends on the species abundance distribution and can thus be also calculated
from abundance data (see rare_curve
).
In contrast the species-accumulation curve starts from a focal individual and only samples the nearest neighbours of the focal individual. The final species accumulation curves is calculated as the mean over the accumulation curves starting from all individuals.
In contrast to the rarefaction curve the accumulation curve is not only influenced by the species abundance distribution, but also by the spatial distribution of individuals.
A dataframe with 2-3 columns. The first column indicates the sample size (numbers of individuals), and the second and third column indicate the expected species richness for spatial sampling (column: "spec_accum") and/or random sampling (column "spec_rarefied")
sim_com1 <- sim_thomas_community(s_pool = 100, n_sim = 1000) sac1 <- spec_sample_curve(sim_com1, method = c("rare","acc")) head(sac1) plot(sac1)
sim_com1 <- sim_thomas_community(s_pool = 100, n_sim = 1000) sac1 <- spec_sample_curve(sim_com1, method = c("rare","acc")) head(sac1) plot(sac1)
Print summary of spatial community object
## S3 method for class 'community' summary(object, digits = 2, ...)
## S3 method for class 'community' summary(object, digits = 2, ...)
object |
Community object of class |
digits |
Integer. Number of digits to print |
... |
Additional arguments passed to |
This function is called for its side effects and has no return value.
Print summary of species abundance distribution object
## S3 method for class 'sad' summary(object, ...)
## S3 method for class 'sad' summary(object, ...)
object |
Community object of class |
... |
Additional arguments passed to |
This function is called for its side effects and has no return value.