Title: | Run Allometric Trophic Networks Models |
---|---|
Description: | Implements the differential equations associated to different versions of Allometric Trophic Models (ATN) to estimate the temporal dynamics of species biomasses in food webs. It offers several features to generate synthetic food webs and to parametrise models as well as a wrapper to the ODE solver deSolve. |
Authors: | Benoit Gauzens [cre, aut], Emilio Berti [aut] |
Maintainer: | Benoit Gauzens <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.1.0 |
Built: | 2024-12-24 06:52:48 UTC |
Source: | CRAN |
Filter Extinct Species
.filter_extinct(df, model)
.filter_extinct(df, model)
df |
deSolve matrix as returned from lsoda_wrapper(). |
model |
ATNr model object, from which extinction threshold is extracted. |
Set to zero species biomass that are below the extinction threshold.
df with values below th set to zero.
Make L matrix
create_Lmatrix(BM, nb_b, Ropt = 100, gamma = 2, th = 0.01)
create_Lmatrix(BM, nb_b, Ropt = 100, gamma = 2, th = 0.01)
BM |
float vector, body mass of species. |
nb_b |
integer, number of basal species. |
Ropt |
numeric, consumer/resource optimal body mass ratio. |
gamma |
numeric, code for the width of the Ricker function. |
th |
float, the threshold below which attack rates are considered = 0. |
The L matrix contains the probability for an attack event to be successful based on allometric rules and a Ricker function defined by Ropt and gamma. If at least one species has not resource or consumer (i.e. it is an isolated species), another food web is generated, until a maximum of 100 iterations.
A numeric matrix with the probability for an attack event between two species to be successful.
set.seed(123) mass <- sort(10 ^ runif(30, 2, 6)) L <- create_Lmatrix(mass, nb_b = 10, Ropt = 100) image(t(L))
set.seed(123) mass <- sort(10 ^ runif(30, 2, 6)) L <- create_Lmatrix(mass, nb_b = 10, Ropt = 100) image(t(L))
Make parameter matrix
create_matrix_parameter(BM, b0, bprey, bpred, E, T.K, T0, k)
create_matrix_parameter(BM, b0, bprey, bpred, E, T.K, T0, k)
BM |
float vector, body mass of species. |
b0 |
const |
bprey |
const |
bpred |
const |
E |
const |
T.K |
Celsius to Kelvin conversion |
T0 |
Default temperature in Kelvin |
k |
Boltzmann constant |
Make a parameter matrix that depends on both predators and prey and that is used to define attack rates and handling times based on the general allometric equation:
A matrix filled with estimated values for a model parameter that depends on prey and predator body masses (see details)
Initialize an ATN model, following Delmas et al. 2017, Methods in Ecology and Evolution
create_model_Scaled(nb_s, nb_b, BM, fw)
create_model_Scaled(nb_s, nb_b, BM, fw)
nb_s |
integer, number of total species. |
nb_b |
integer, number of basal species. |
BM |
float vector, body mass of species. |
fw |
binary adjacency matrix of the food web. |
A model is defined by the total number of species (nb_s), the number of basal species (nb_b), the number of nutrients (nb_n), the body masses (BM) of species, and the adjacency matrix (fw) representing species interactions.
An object of class ATN (Rcpp_parameters_prefs).
Delmas, E., Brose, U., Gravel, D., Stouffer, D.B. and Poisot, T. (2017), Simulations of biomass dynamics in community food webs. Methods Ecol Evol, 8: 881-886. https://doi.org/10.1111/2041-210X.12713
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Scaled(n_species, n_basal, masses, fw)
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Scaled(n_species, n_basal, masses, fw)
Initialize an ATN model, following Binzer et al. 2016, Global Change Biology
create_model_Unscaled(nb_s, nb_b, BM, fw)
create_model_Unscaled(nb_s, nb_b, BM, fw)
nb_s |
integer, number of total species. |
nb_b |
integer, number of basal species. |
BM |
float vector, body mass of species. |
fw |
binary adjacency matrix of the food web. |
A model is defined by the total number of species (nb_s), the number of basal species (nb_b), the number of nutrients (nb_n), the body masses (BM) of species, and the adjacency matrix (fw) representing species interactions.
An object of class ATN (Rcpp_parameters_prefs).
Binzer, A., Guill, C., Rall, B.C. and Brose, U. (2016), Interactive effects of warming, eutrophication and size structure: impacts on biodiversity and food-web structure. Glob Change Biol, 22: 220-227. https://doi.org/10.1111/gcb.13086 Gauzens, B., Rall, B.C., Mendonca, V. et al. Biodiversity of intertidal food webs in response to warming across latitudes. Nat. Clim. Chang. 10, 264-269 (2020). https://doi.org/10.1038/s41558-020-0698-z
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 masses <- sort(10^runif(n_species, 1, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Unscaled(n_species, n_basal, masses, fw)
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 masses <- sort(10^runif(n_species, 1, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Unscaled(n_species, n_basal, masses, fw)
Initialize an ATN model, following Schneider et al. 2016, Nature Communication
create_model_Unscaled_nuts(nb_s, nb_b, nb_n = 2, BM, fw)
create_model_Unscaled_nuts(nb_s, nb_b, nb_n = 2, BM, fw)
nb_s |
integer, number of total species. |
nb_b |
integer, number of basal species. |
nb_n |
integer, number of nutrients. |
BM |
float vector, body mass of species. |
fw |
binary adjacency matrix of the food web. |
A model is defined by the total number of species (nb_s), the number of basal species (nb_b), the number of nutrients (nb_n), the body masses (BM) of species, and the adjacency matrix (fw) representing species interactions. Nutrients are not counted as species.
An object of class ATN (Rcpp_parameters_prefs).
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 n_nutrients <- 2 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Unscaled_nuts(n_species, n_basal, n_nutrients, masses, fw)
library(ATNr) set.seed(123) n_species <- 50 n_basal <- 20 n_nutrients <- 2 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species L <- create_Lmatrix(masses, n_basal) fw <- L fw[fw > 0] <- 1 mod <- create_model_Unscaled_nuts(n_species, n_basal, n_nutrients, masses, fw)
Function to generate a food web based on the niche model (Williams and Martinez, 2000) based on the number of species and connectance. Corrections from Allesina et al. (2008) are used.
create_niche_model(S, C)
create_niche_model(S, C)
S |
integer, number of species. |
C |
numeric, connectance i.e. the number of realized links over the all possible links. |
If at least one species has not resource or consumer (i.e. it is an isolated species), another food web is generated, until a maximum of 100 iterations.
A (square) matrix with zeros (no interaction) and ones (species j consume species i).
Williams, R. J., & Martinez, N. D. (2000). Simple rules yield complex food webs. Nature, 404(6774), 180-183.
Allesina, S., Alonso, D., & Pascual, M. (2008). A general model for food web structure. science, 320(5876), 658-661.
set.seed(123) web_niche <- create_niche_model(30, .1) image(t(web_niche))
set.seed(123) web_niche <- create_niche_model(30, .1) image(t(web_niche))
Initialise the default parametrisation for the scaled version of the ATN model as in Delmas et al. (2016).
initialise_default_Scaled(model)
initialise_default_Scaled(model)
model |
an object of class Rcpp_Scaled. |
An object of class Rcpp_Scaled with default parameters as in Delmas et al. (2017).
Delmas, E., Brose, U., Gravel, D., Stouffer, D.B. and Poisot, T. (2017), Simulations of biomass dynamics in community food webs. Methods Ecol Evol, 8: 881-886. https://doi.org/10.1111/2041-210X.12713
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Scaled(20, 10, BM = masses, fw = L) mod <- initialise_default_Scaled(mod)
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Scaled(20, 10, BM = masses, fw = L) mod <- initialise_default_Scaled(mod)
Initialise the default parametrisation for the scaled version of the ATN model as in Binzer et al. (2016), with updates from Gauzens et al. 2020
initialise_default_Unscaled(model, temperature = 20)
initialise_default_Unscaled(model, temperature = 20)
model |
an object of class ATN (Rcpp_Unscaled). |
temperature |
numeric, ambient temperature of the ecosystem in Celsius. |
An object of class ATN (Rcpp_Unscaled) with default parameters as in Delmas et al. (2017).
Binzer, A., Guill, C., Rall, B. C. & Brose, U. Interactive effects of warming, eutrophication and size structure: impacts on biodiversity and food-web structure. Glob. Change Biol. 22, 220-227 (2016). Gauzens, B., Rall, B.C., Mendonca, V. et al. Biodiversity of intertidal food webs in response to warming across latitudes. Nat. Clim. Chang. 10, 264-269 (2020). https://doi.org/10.1038/s41558-020-0698-z
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled(20, 10, masses, L) mod <- initialise_default_Unscaled(mod)
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled(20, 10, masses, L) mod <- initialise_default_Unscaled(mod)
Initialise the default parametrisation for the model for Schneider et al. (2016).
initialise_default_Unscaled_nuts(model, L.mat, temperature = 20)
initialise_default_Unscaled_nuts(model, L.mat, temperature = 20)
model |
an object of class ATN (Rcpp_Unscaled_nuts. |
L.mat |
numeric matrix, probability of a consumer to attack and capture an encountered resource. See |
temperature |
numeric, ambient temperature of the ecosystem in Celsius. |
An object of class ATN (Rcpp_Unscaled_nuts) with default parameters as in Schneider et al. (2016).
Schneider, F. D., Brose, U., Rall, B. C., & Guill, C. (2016). Animal diversity and ecosystem functioning in dynamic food webs. Nature Communications, 7(1), 1-8.
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L)
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L)
Run a deep search first algorithm (DFS)
is_connected(fw)
is_connected(fw)
fw |
binary adjacency matrix of the food web. |
Boolean: TRUE if the food web is connected, FALSE if several disconnected sub-networks are detected.
library(ATNr) set.seed(123) # number of species, nutrients, and body masses n_species <- 20 n_basal <- 5 n_nutrients <- 3 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species # create food web matrix L <- create_Lmatrix(masses, n_basal) L[, 1:n_basal] <- 0 fw <- L fw[fw > 0] <- 1 connected <- is_connected(fw)
library(ATNr) set.seed(123) # number of species, nutrients, and body masses n_species <- 20 n_basal <- 5 n_nutrients <- 3 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species # create food web matrix L <- create_Lmatrix(masses, n_basal) L[, 1:n_basal] <- 0 fw <- L fw[fw > 0] <- 1 connected <- is_connected(fw)
Estimate the Jacobian matrix of a ODE system
jacobian(bioms, ODE, eps = 1e-06)
jacobian(bioms, ODE, eps = 1e-06)
bioms |
float vector, biomass of species. |
ODE |
function that computes the ODEs from one of the model available |
eps |
float, scale precision of the numerical approximation. |
The function provides a numerical estimation of the Jacobian matrix based on the 5 points stencil method. The precision of the method is in
, where
. The choice of eps should ensure that
is always lower to the extinction threshold.
The dimension of the Jacobian matrix are not always matching the number of species in the system. This is because we considered that a perturbation can not correspond to the recolonisation of an extinct species. Therefore, extinct species are removed from the system to calculate the Jacobian matrix.
A matrix corresponding to the Jacobian of the system estimated at the parameter biomasses
library(ATNr) set.seed(123) # first run a model to reach equilibrium masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L) biomasses <- masses ^ -0.75 * 10 ^ 4 #biomasses of species biomasses <- append(runif(3, 20, 30), biomasses) times <- seq(0, 100, 1) sol <- lsoda_wrapper(times, biomasses, mod) # get the final biomasses final.bioms = sol[nrow(sol), -1] # estimate jacobian jacobian(final.bioms, mod$ODE)
library(ATNr) set.seed(123) # first run a model to reach equilibrium masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L) biomasses <- masses ^ -0.75 * 10 ^ 4 #biomasses of species biomasses <- append(runif(3, 20, 30), biomasses) times <- seq(0, 100, 1) sol <- lsoda_wrapper(times, biomasses, mod) # get the final biomasses final.bioms = sol[nrow(sol), -1] # estimate jacobian jacobian(final.bioms, mod$ODE)
This is a wrapper to call lsoda
from
deSolve and solve the ODE.
Package deSolve
needs to be installed to run
this wrapper.
lsoda_wrapper(t, y, model, verbose = FALSE, ...)
lsoda_wrapper(t, y, model, verbose = FALSE, ...)
t |
vector of times. |
y |
vector of biomasses. |
model |
object of class ATN (Rcpp_parameters_prefs). |
verbose |
Boolean, whether a message should be printed when all checks were successful |
... |
additional arguments to pass to 'lsoda' |
A matrix for the ODE solution with species as columns and times as rows.
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L) biomasses <- masses ^ -0.75 * 10 ^ 4 #biomasses of species biomasses <- append(runif(3, 20, 30), biomasses) times <- seq(0, 100, 1) sol <- lsoda_wrapper(times, biomasses, mod) range(sol[, -1]) mod$ext <- 1e-3 sol <- lsoda_wrapper(times, biomasses, mod)
library(ATNr) set.seed(123) masses <- runif(20, 10, 100) #body mass of species L <- create_Lmatrix(masses, 10, Ropt = 10) L[L > 0] <- 1 mod <- create_model_Unscaled_nuts(20, 10, 3, masses, L) mod <- initialise_default_Unscaled_nuts(mod, L) biomasses <- masses ^ -0.75 * 10 ^ 4 #biomasses of species biomasses <- append(runif(3, 20, 30), biomasses) times <- seq(0, 100, 1) sol <- lsoda_wrapper(times, biomasses, mod) range(sol[, -1]) mod$ext <- 1e-3 sol <- lsoda_wrapper(times, biomasses, mod)
Plot solution of the ODE for the food web. Currently only species and not nutrients are plotted.
plot_odeweb(x, nb_s)
plot_odeweb(x, nb_s)
x |
matrix with solutions. First row should be the time vector. |
nb_s |
numeric, number of species as in the model (e.g.,
|
No return value, called for side effects.
## Not run: library(ATNr) library(deSolve) set.seed(123) # number of species, nutrients, and body masses n_species <- 20 n_basal <- 5 n_nutrients <- 3 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species # create food web matrix L <- create_Lmatrix(masses, n_basal) L[, 1:n_basal] <- 0 fw <- L fw[fw > 0] <- 1 model <- create_model_Unscaled_nuts( n_species, n_basal, n_nutrients, masses, fw ) # initialize model as default in Schneider et al. (2016) model <- initialise_default_Unscaled_nuts(model, L) # defining integration time times <- seq(0, 500, 5) biomasses <- runif(n_species + n_nutrients, 2, 3) sol <- lsoda_wrapper(times, biomasses, model, verbose = FALSE) plot_odeweb(sol, model$nb_s) ## End(Not run)
## Not run: library(ATNr) library(deSolve) set.seed(123) # number of species, nutrients, and body masses n_species <- 20 n_basal <- 5 n_nutrients <- 3 masses <- sort(10^runif(n_species, 2, 6)) #body mass of species # create food web matrix L <- create_Lmatrix(masses, n_basal) L[, 1:n_basal] <- 0 fw <- L fw[fw > 0] <- 1 model <- create_model_Unscaled_nuts( n_species, n_basal, n_nutrients, masses, fw ) # initialize model as default in Schneider et al. (2016) model <- initialise_default_Unscaled_nuts(model, L) # defining integration time times <- seq(0, 500, 5) biomasses <- runif(n_species + n_nutrients, 2, 3) sol <- lsoda_wrapper(times, biomasses, model, verbose = FALSE) plot_odeweb(sol, model$nb_s) ## End(Not run)
Function to remove species from a model class
remove_species(species, model, nuts = NULL)
remove_species(species, model, nuts = NULL)
species |
integer vector, the indices of species to remove. |
model |
model object |
nuts |
integer vector, the indices of nutrients to remove. Parameter specific to the Unscaled_nuts model. |
A model object where the data structure has bee updated to remove the species in parameters.
Check if the dimensions of vectors and matrices used in the model are correct. If any dimension is not correct, an error message is returned.
run_checks(model, verbose = TRUE)
run_checks(model, verbose = TRUE)
model |
a model object. |
verbose |
Boolean, whether a message should be printed when all checks were successful |
No return value, only throw an error if parameters are inconsistent.
Type the name of the class to see its methods
nb_s
Total number of species
nb_b
Number of basal species
c
double: interference competition
X
Vector of metabolic rates (length = number of species)
max_feed
Vector of maximum feeding rates (length = number of consumers)
e
Vector of assimilation efficiencies (length = number of species)
r
Vector of producers maximum growth rates (length = number of basal species)
BM
Vector of body masses (length = number of species)
dB
Vector of local derivatives (length = number of species)
B0
Vector of half saturation densities (length = number of consumers)
fw
Adjacency matrix of the food-web (dim = number of species * number of species)
w
Matrix of relative consumption rates (dim = number of species * number of consumers)
F
Matrix of per-capita feeding rates (dim = number of species * number of consumers)
q
hill exponent for the type of functional response
K
Carrying capacity of basal species
ext
Extinction threshold for species
alpha
Plant resource competition
ODE
Calculate the derivatives for the scaled version of the ATN model
Parameter: bioms - Local species biomasses
Parameter: t - Integration time point
Returns a vector of growth rate for each species at time t
To not use. For testing purpose only. Please use Rcpp_Scaled instead.
A dataset containing the default parameters as in the Schneider et al. (2016)
and used to parametrize the default models. See also
create_model_Unscaled_nuts
, create_Lmatrix
,
initialise_default_Unscaled_nuts
.
schneider
schneider
A list with the default parameters:
ambient temperature in Celsius
default temperature, 20 degree Celsius in Kelvin
Boltzmann's constant
20 degree Celsius in Kelvin, used to estimate scaling law of metabolic rates
Hill's exponent of the functional response
consumer/resource optimal body mass ratio
shape of the Ricker function
average predator interference
standard deviation of predator interference
Activation energy for interference
scaling constant of the power-law of handling time with consumer and resource body mass
exponent associated to predator body mass for the allometric scaling of handling time
exponent associated to prey body mass for the allometric scaling of handling time
Activation energy for handling time
normalisation constant for capture coefficient
exponent associated to prey body mass for the allometric scaling of capture coefficient
exponent associated to predator body mass for the allometric scaling of capture coefficient
Activation energy for capture coefficient
Assimilation efficiency associated to the consumption of a plant species
Assimilation efficiency associated to the consumption of an animal species
scaling constant of the power-law of metabolic demand per unit of plant biomass
scaling constant of the power-law of metabolic demand per unit of animal biomass
Activation energy for metabolic rates
TBD
turnover rate of nutrients
Minimum uptake efficiency of plants
Maximum uptake efficiency of plants
Average maximum nutrient concentration
standard deviation of maximum nutrient concentration
relative content of nutrient 1 in plant biomass
Schneider, F. D., Brose, U., Rall, B. C., & Guill, C. (2016). Animal diversity and ecosystem functioning in dynamic food webs. Nature Communications, 7(1), 1-8.
Sort custom input
sort_input(BM, fw)
sort_input(BM, fw)
BM |
numeric vector, body mass of species. |
fw |
adjacency matrix of the food web. |
Body masses and food web matrix should be arranged with the first elements/columns being for basal species. This is a requirement for the Cpp class and must be enforced before initializing the Rcpp_Schneider and Rcpp_Delmas objects.
A list with sorted body masses (body.mass) and food web matrix (food.web).
bm <- runif(10, 10, 50) fw <- matrix(as.numeric(runif(100) > .9), 10, 10) sort_input(bm, fw)
bm <- runif(10, 10, 50) fw <- matrix(as.numeric(runif(100) > .9), 10, 10) sort_input(bm, fw)
Calculate trophic level of species
TroLev(fw)
TroLev(fw)
fw |
numeric matrix, the matrix of the food web. |
A numeric vector of species' trophic level.
Type the name of the class to see its methods
nb_s
Total number of species
nb_b
Number of basal species
c
double: interference competition
X
Vector of metabolic rates (length = number of species)
a
Matrix of attack rates (dim = number of species * number of consumers)
h
Matrix of handling times (dim = number of species * number of consumers)
e
Vector of assimilation efficiencies (length = number of species)
r
Vector of producers maximum growth rates (length = number of basal species)
BM
Vector of body masses (length = number of species)
dB
Vector of local derivatives (length = number of species)
fw
Adjacency matrix of the food-web (dim = number of species * number of species)
F
Matrix of per-capita feeding rates (dim = number of species * number of consumers)
q
hill exponent for the type of functional response
K
Carrying capacity of basal species
alpha
Plant resource competition
ext
Extinction threshold for species
ODE
Calculate the derivatives for the scaled version of the ATN model
Parameter: bioms - Local species biomasses
Parameter: t - Integration time point
Returns a vector of growth rate for each species at time t
To not use. For testing purpose only. Please use Rcpp_Unscaled instead.
Type the name of the class to see its methods
nb_s
Total number of species
nb_b
Number of basal species
nb_n
Number of nutrient pool
c
double: interference competition
b
Matrix of attack rates (dim = number of species * number of consumers)
h
Matrix of handling times (dim = number of species * number of consumers)
X
vector of metabolic rates (length = number of species)
K
matrix of plant nutrient efficiencies (dim = number of nutrients * number of plants)
V
matrix of plant relative nutrient content (dim = number of nutrients * number of plants)
S
Vector of maximum nutrient concentration (length = number of plants)
r
Vector of maximum growth rate of plant species (length = number of plant species)
e
Vector of assimilation efficiencies (length = number of species)
BM
Vector of body masses (length = number of species)
dB
Vector of local derivatives (length = number of species)
fw
Adjacency matrix of the food-web (dim = number of species * number of species)
w
Matrix of relative consumption rates (dim = number of species * number of consumers)
F
Matrix of per-capita feeding rates (dim = number of species * number of consumers)
q
hill exponent for the type of functional response
ext
Extinction threshold for species
ODE
Calculate the derivatives for the scaled version of the ATN model
Parameter: bioms - Local species biomasses
Parameter: t - Integration time point
Returns a vector of growth rate for each species at time t
To not use. For testing purpose only. Please use Rcpp_Unscaled_nuts instead.