Functions

Index

API documentation

PlantBiophysics.AbstractEnergy_BalanceModelType

energy_balance process abstract model.

All models implemented to simulate the energy_balance process must be a subtype of this type, e.g. struct MyEnergy_BalanceModel <: AbstractEnergy_BalanceModel end.

You can list all models implementing this process using subtypes:

Examples

subtypes(AbstractEnergy_BalanceModel)

Energy balance process. This process computes the energy balance of objects, meaning that it computes the net radiation, the sensible heat flux, and the latent heat flux if necessary. It can be coupled with a photosynthesis model in the case of plants leaves.

At the moment, two models are implemented in the package:

  • Monteith: the model found in Monteith and Unsworth (2013)
  • Missing: if no computation of the energy balance is needed

Note

Some models need input values for some variables. For example Monteith requires a value for Ra_SW_f, d and sky_fraction. If you read the models from a file, you can use init_status! (see examples).

Examples

using PlantMeteo, PlantSimEngine, PlantBiophysics

# ---Simple example---

meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

# Using the model of Monteith and Unsworth (2013) for energy, Farquhar et al. (1980) for
# photosynthesis, and Medlyn et al. (2011) for stomatal conductance:
leaf =
    ModelList(
        energy_balance = Monteith(),
        photosynthesis = Fvcb(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Ra_SW_f = 13.747, sky_fraction = 1.0, aPPFD = 1500.0, d = 0.03)
    )

run!(leaf,meteo)

# ---Using several components---

leaf2 = copy(leaf)
leaf2[:aPPFD] = 800.0

run!([leaf,leaf2],meteo)

# You can use a Dict if you'd like to keep track of the leaf in the returned DataFrame:
run!(Dict(:leaf1 => leaf, :leaf2 => leaf2), meteo)

# ---Using several meteo time-steps---

w = Weather(
    [
        Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65),
        Atmosphere(T = 25.0, Wind = 1.5, P = 101.3, Rh = 0.55)
    ],
    (site = "Test site",)
)

leaf =
    ModelList(
        energy_balance = Monteith(),
        photosynthesis = Fvcb(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Ra_SW_f = [12.0,13.747], sky_fraction = 1.0, aPPFD = 1500.0, d = 0.03)
    )

run!(leaf, w)

# ---Using several meteo time-steps and several components---

leaf2 =
    ModelList(
        energy_balance = Monteith(),
        photosynthesis = Fvcb(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Ra_SW_f = [12.0,13.747], sky_fraction = 1.0, aPPFD = 1500.0, d = 0.01)
    )

run!(Dict(:leaf1 => leaf, :leaf2 => leaf2), w)

# ---Using a model file---

model = read_model(joinpath(dirname(dirname(pathof(PlantBiophysics))),"test","inputs","models","plant_coffee.yml"))

# An example model file is available here:
# "https://raw.githubusercontent.com/VEZY/PlantBiophysics/main/test/inputs/models/plant_coffee.yml"

# Initialising the mandatory variables:
init_status!(model, Ra_SW_f = 13.747, sky_fraction = 1.0, aPPFD = 1500.0, Tₗ = 25.0, d = 0.03)

# NB: To know which variables has to be initialized according to the models used, you can use
# `to_initialize(ComponentModels)`, *e.g.*:
to_initialize(model["Leaf"])

# Running a simulation for all component types in the same scene:
run!(model, meteo)

model["Leaf"].status.Rn
model["Leaf"].status.A
model["Leaf"].status.Cᵢ

# ---Simulation on a full plant using an MTG---

using PlantBiophysics, MultiScaleTreeGraph, PlantGeom, GLMakie, Dates, PlantMeteo

file = joinpath(dirname(dirname(pathof(PlantBiophysics))), "test", "inputs", "scene", "opf", "coffee.opf")
mtg = read_opf(file)

# Import the meteorology:
met_file = joinpath(dirname(dirname(pathof(PlantMeteo))), "test", "data", "meteo.csv")

meteo = read_weather(
    met_file,
    :temperature => :T,
    :relativeHumidity => (x -> x ./ 100) => :Rh,
    :wind => :Wind,
    :atmosphereCO2_ppm => :Cₐ,
    date_format = DateFormat("yyyy/mm/dd")
)

# Make the models:
models = Dict(
    "Leaf" =>
        ModelList(
            energy_balance = Monteith(),
            photosynthesis = Fvcb(),
            stomatal_conductance = Medlyn(0.03, 12.0),
            status = (d = 0.03,)
        )
)

# List the MTG attributes:
names(mtg)
# We have the skyFraction already, but not Ra_SW_f and aPPFD, so we must compute them first.
# Ra_SW_f is the shortwave radiation (or global radiation), so it is the sum of Ra_PAR_f and Ra_NIR_f.
# aPPFD is the PAR in μmol m-2 s-1, so Ra_PAR_f * 4.57.

# We can compute them using the following code (transform! comes from MultiScaleTreeGraph.jl):
transform!(
    mtg,
    [:Ra_PAR_f, :Ra_NIR_f] => ((x, y) -> x + y) => :Ra_SW_f,
    :Ra_PAR_f => (x -> x * 4.57) => :aPPFD,
    ignore_nothing = true
)

# We can now initialize the models in the mtg:
init_mtg_models!(mtg, models, length(meteo))

# Making the simulation:
run!(mtg, meteo)

# Pull the leaf temperature of the first step:
transform!(
    mtg,
    :Tₗ => (x -> x[1]) => :Tₗ_1,
    ignore_nothing = true
)

# Vizualise the output:
viz(mtg, color = :Tₗ_1)

References

Duursma, R. A., et B. E. Medlyn. 2012. « MAESPA: a model to study interactions between water limitation, environmental drivers and vegetation function at tree and stand levels, with an example application to [CO2] × drought interactions ». Geoscientific Model Development 5 (4): 919‑40. https://doi.org/10.5194/gmd-5-919-2012.

Monteith, John L., et Mike H. Unsworth. 2013. « Chapter 13 - Steady-State Heat Balance: (i) Water Surfaces, Soil, and Vegetation ». In Principles of Environmental Physics (Fourth Edition), edited by John L. Monteith et Mike H. Unsworth, 217‑47. Boston: Academic Press.

Schymanski, Stanislaus J., et Dani Or. 2017. « Leaf-Scale Experiments Reveal an Important Omission in the Penman–Monteith Equation ». Hydrology and Earth System Sciences 21 (2): 685‑706. https://doi.org/10.5194/hess-21-685-2017.

Vezy, Rémi, Mathias Christina, Olivier Roupsard, Yann Nouvellon, Remko Duursma, Belinda Medlyn, Maxime Soma, et al. 2018. « Measuring and modelling energy partitioning in canopies of varying complexity using MAESPA model ». Agricultural and Forest Meteorology 253‑254 (printemps): 203‑17. https://doi.org/10.1016/j.agrformet.2018.02.005.

source
PlantBiophysics.AbstractLight_InterceptionModelType

light_interception process abstract model.

All models implemented to simulate the light_interception process must be a subtype of this type, e.g. struct MyLight_InterceptionModel <: AbstractLight_InterceptionModel end.

You can list all models implementing this process using subtypes:

Examples

subtypes(AbstractLight_InterceptionModel)

Light interception process. Available as object.light_interception.

At the moment, two models are implemented in the package:

  • Beer: the Beer-Lambert law of ligth extinction
  • LightIgnore: ignore the computation of light interception (this one is for backward

compatibility with ARCHIMED-ϕ)

Examples

using PlantSimEngine, PlantBiophysics, PlantMeteo
m = ModelList(light_interception=Beer(0.5), status=(LAI=2.0,))

meteo = Atmosphere(T=20.0, Wind=1.0, P=101.3, Rh=0.65, Ri_PAR_f=300.0)

run!(m, meteo)

m[:aPPFD]
source
PlantBiophysics.AbstractPhotosynthesisModelType

photosynthesis process abstract model.

All models implemented to simulate the photosynthesis process must be a subtype of this type, e.g. struct MyPhotosynthesisModel <: AbstractPhotosynthesisModel end.

You can list all models implementing this process using subtypes:

Examples

subtypes(AbstractPhotosynthesisModel)

Photosynthesis process to compute the CO₂ assimilation, and potentially hard-coupled with a stomatal conductance process.

The models used are defined by the types of the photosynthesis and stomatal_conductance fields of the ModelList. For exemple to use the implementation of the Farquhar–von Caemmerer–Berry (FvCB) model, use the type Fvcb (see example below).

Examples

using PlantSimEngine, PlantMeteo, PlantBiophysics

meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

# Using Fvcb model:
leaf =
    ModelList(
        photosynthesis = Fvcb(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Tₗ = 25.0, aPPFD = 1000.0, Cₛ = 400.0, Dₗ = meteo.VPD)
    )

run!(leaf, meteo)

# ---Using several components---

leaf2 = copy(leaf)
leaf2.status.aPPFD = 800.0

run!([leaf,leaf2],meteo)

# ---Using several meteo time-steps---

w = Weather(
        [
            Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65),
            Atmosphere(T = 25.0, Wind = 1.5, P = 101.3, Rh = 0.55)
        ],
        (site = "Test site,)
    )

run!(leaf, w)

# ---Using several meteo time-steps and several components---

run!(Dict(:leaf1 => leaf, :leaf2 => leaf2), w)

# Using a model file:

model = read_model("a-model-file.yml")

# Initialising the mandatory variables:
init_status!(model, Tₗ = 25.0, aPPFD = 1000.0, Cₛ = 400.0, Dₗ = meteo.VPD)

# Running a simulation for all component types in the same scene:
run!(model, meteo)
model["Leaf"].status.A

Note that we use VPD as an approximation of Dₗ here because we don't have the leaf temperature (i.e. Dₗ = VPD when Tₗ = T).

source
PlantBiophysics.AbstractStomatal_ConductanceModelType

stomatal_conductance process abstract model.

All models implemented to simulate the stomatal_conductance process must be a subtype of this type, e.g. struct MyStomatal_ConductanceModel <: AbstractStomatal_ConductanceModel end.

You can list all models implementing this process using subtypes:

Examples

subtypes(AbstractStomatal_ConductanceModel)

Process for the stomatal conductance for CO₂ (mol m⁻² s⁻¹), it takes the form:

leaf.stomatal_conductance.g0 + gs_closure(leaf,meteo) * leaf.status.A

where gsclosure(leaf,meteo) computes the stomatal closure, and must be implemented for the type of `leaf.stomatalconductance. The stomatal conductance is not allowed to go belowleaf.stomatalconductance.gsmin`.

Arguments

  • Gs::Gsm: a stomatal conductance model, usually the leaf model (i.e. leaf.stomatal_conductance)
  • models::ModelList: A leaf struct holding the parameters for the model. See

ModelList, and Medlyn or ConstantGs for the conductance models.

  • status::Status: A status, usually the leaf status (i.e. leaf.status)
  • gs_mod: the output from a gs_closure implementation (the conductance models

generally only implement this function)

  • meteo<:PlantMeteo.AbstractAtmosphere: meteo data, see Atmosphere

Examples

using PlantMeteo, PlantSimEngine, PlantBiophysics
meteo = Atmosphere(T = 22.0, Wind = 0.8333, P = 101.325, Rh = 0.4490995)

# Using a constant value for Gs:

leaf =
    ModelList(
        stomatal_conductance = Medlyn(0.03,12.0), # Instance of a Medlyn type
        status = (A = 20.0, Cₛ = 380.0, Dₗ = meteo.VPD)
    )

# Computing the stomatal conductance using the Medlyn et al. (2011) model:
run!(leaf,meteo)
source
PlantBiophysics.BeerType
Beer(k)

Beer-Lambert law for light interception.

Required inputs: LAI in m² m⁻². Required meteorology data: Ri_PAR_f, the incident flux of atmospheric radiation in the PAR, in W m[soil]⁻² (== J m[soil]⁻² s⁻¹).

Output: aPPFD, the absorbed Photosynthetic Photon Flux Density in μmol[PAR] m[leaf]⁻² s⁻¹.

source
PlantBiophysics.BeerShortwaveType
BeerShortwave(k, f)
BeerShortwave(k)

The Beer-Lambert law for light interception for the shortwave radiation.

Arguments

  • k_PAR: extinction coefficient for the PAR
  • k_NIR: extinction coefficient for the NIR, taken equal to the one for the PAR by default

Required inputs

  • LAI: the leaf area index (m² m⁻²)
  • Ri_PAR_f: (from meteorology) the incident flux of atmospheric radiation in the PAR, in W m[soil]⁻² (== J m[soil]⁻² s⁻¹).
  • Ri_NIR_f: (from meteorology) the incident flux of atmospheric radiation in the NIR, in W m[soil]⁻² (== J m[soil]⁻² s⁻¹).

Outputs

  • aPPFD: the absorbed Photosynthetic Photon Flux Density in μmol[PAR] m[leaf]⁻² s⁻¹.
  • Ra_PAR_f: the absorbed PAR in W m[leaf]⁻².
  • Ra_NIR_f: the absorbed NIR in W m[leaf]⁻².
  • Ra_SW_f: the absorbed shortwave radiation in W m[leaf]⁻².

Examples

using PlantSimEngine, PlantBiophysics, PlantMeteo

m = ModelList(light_interception=BeerShortwave(0.5), status=(LAI=2.0,))
meteo = Atmosphere(T=20.0, Wind=1.0, P=101.3, Rh=0.65, Ri_PAR_f=300.0, Ri_NIR_f=280.0)

run!(m, meteo)

m[:aPPFD]
m[:Ra_SW_f]
m[:Ra_PAR_f]
m[:Ra_NIR_f]
source
PlantBiophysics.ConstantAGsType

Constant (forced) assimilation, given in $μmol\ m^{-2}\ s^{-1}$, coupled with a stomatal conductance model that helps computing Cᵢ.

Examples

ConstantAGs(30.0)
source
PlantBiophysics.ConstantGsType

Constant stomatal conductance for CO₂ struct.

Arguments

  • g0: intercept (only used when calling from a photosynthesis model, e.g. Fvcb).
  • Gₛ: stomatal conductance.

Then used as follows: Gs = ConstantGs(0.0,0.1)

source
PlantBiophysics.FvcbType

Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981) coupled with a conductance model.

Parameters

  • Tᵣ: the reference temperature (°C) at which other parameters were measured
  • VcMaxRef: maximum rate of Rubisco activity ($μmol\ m^{-2}\ s^{-1}$)
  • JMaxRef: potential rate of electron transport ($μmol\ m^{-2}\ s^{-1}$)
  • RdRef: mitochondrial respiration in the light at reference temperature ($μmol\ m^{-2}\ s^{-1}$)
  • TPURef: triose phosphate utilization-limited photosynthesis rate ($μmol\ m^{-2}\ s^{-1}$)
  • Eₐᵣ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for Rd.
  • O₂: intercellular dioxygen concentration ($ppm$)
  • Eₐⱼ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for JMax.
  • Hdⱼ: rate of decrease of the function above the optimum (also called EDVJ) for JMax.
  • Δₛⱼ: entropy factor for JMax.
  • Eₐᵥ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for VcMax.
  • Hdᵥ: rate of decrease of the function above the optimum (also called EDVC) for VcMax.
  • Δₛᵥ: entropy factor for VcMax.
  • α: quantum yield of electron transport ($mol_e\ mol^{-1}_{quanta}$). See also eq. 4 of

Medlyn et al. (2002), equation 9.16 from von Caemmerer et al. (2009) ((1-f)/2) and its implementation in get_J

  • θ: determines the curvature of the light response curve for J~aPPFD. See also eq. 4 of

Medlyn et al. (2002) and its implementation in get_J

Note on parameters

The default values of the temperature correction parameters are taken from plantecophys. If there is no negative effect of high temperatures on the reaction (Jmax or VcMax), then Δₛ can be set to 0.0.

θ is taken at 0.7 according to (Von Caemmerer, 2000) but it can be modified to 0.9 as in (Su et al., 2009). The larger it is, the lower the smoothing.

α is taken at 0.425 as proposed in von Caemmerer et al. (2009) eq. 9.16, where α = (1-f)/2.

Medlyn et al. (2002) found relatively low influence ("a slight effect") of α and θ. They also say that Kc, Ko and Γ* "are thought to be intrinsic properties of the Rubisco enzyme and are generally assumed constant among species".

See also

References

Caemmerer, S. von, et G. D. Farquhar. 1981. « Some Relationships between the Biochemistry of Photosynthesis and the Gas Exchange of Leaves ». Planta 153 (4): 376‑87. https://doi.org/10.1007/BF00384257.

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Medlyn, B. E., E. Dreyer, D. Ellsworth, M. Forstreuter, P. C. Harley, M. U. F. Kirschbaum, X. Le Roux, et al. 2002. « Temperature response of parameters of a biochemically based model of photosynthesis. II. A review of experimental data ». Plant, Cell & Environment 25 (9): 1167‑79. https://doi.org/10.1046/j.1365-3040.2002.00891.x.

Su, Y., Zhu, G., Miao, Z., Feng, Q. and Chang, Z. 2009. « Estimation of parameters of a biochemically based model of photosynthesis using a genetic algorithm ». Plant, Cell & Environment, 32: 1710-1723. https://doi.org/10.1111/j.1365-3040.2009.02036.x.

Von Caemmerer, Susanna. 2000. Biochemical models of leaf photosynthesis. Csiro publishing.

Duursma, R. A. 2015. « Plantecophys - An R Package for Analysing and Modelling Leaf Gas Exchange Data ». PLoS ONE 10(11): e0143346. https://doi:10.1371/journal.pone.0143346.

Examples

Get the fieldnames:
fieldnames(Fvcb)
# Using default values for the model:
A = Fvcb()

A.Eₐᵥ
source
PlantBiophysics.FvcbIterType

Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981).

Iterative implementation, i.e. the assimilation is computed iteratively over Cᵢ.

Parameters

  • Tᵣ: the reference temperature (°C) at which other parameters were measured
  • VcMaxRef: maximum rate of Rubisco activity ($μmol\ m^{-2}\ s^{-1}$)
  • JMaxRef: potential rate of electron transport ($μmol\ m^{-2}\ s^{-1}$)
  • RdRef: mitochondrial respiration in the light at reference temperature ($μmol\ m^{-2}\ s^{-1}$)
  • TPURef: triose phosphate utilization-limited photosynthesis rate ($μmol\ m^{-2}\ s^{-1}$)
  • Eₐᵣ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for Rd.
  • O₂: intercellular dioxygen concentration ($ppm$)
  • Eₐⱼ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for JMax.
  • Hdⱼ: rate of decrease of the function above the optimum (also called EDVJ) for JMax.
  • Δₛⱼ: entropy factor for JMax.
  • Eₐᵥ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for VcMax.
  • Hdᵥ: rate of decrease of the function above the optimum (also called EDVC) for VcMax.
  • Δₛᵥ: entropy factor for VcMax.
  • α: quantum yield of electron transport ($mol_e\ mol^{-1}_{quanta}$). See also eq. 4 of

Medlyn et al. (2002), equation 9.16 from von Caemmerer et al. (2009) ((1-f)/2) and its implementation in get_J

  • θ: determines the curvature of the light response curve for J~aPPFD. See also eq. 4 of

Medlyn et al. (2002) and its implementation in get_J

  • iter_A_max::Int: maximum number of iterations allowed for the iteration on the assimilation.
  • ΔT_A::T = 1: threshold bellow which the assimilation is considered constant. Given in

percent of change, i.e. 1% means that two successive assimilations with less than 1% difference in value are considered the same value.

Note on parameters

The default values of the temperature correction parameters are taken from plantecophys. If there is no negative effect of high temperatures on the reaction (Jmax or VcMax), then Δₛ can be set to 0.0.

θ is taken at 0.7 according to (Von Caemmerer, 2000) but it can be modified to 0.9 as in (Su et al., 2009). The larger it is, the lower the smoothing.

α is taken at 0.425 as proposed in von Caemmerer et al. (2009) eq. 9.16, where α = (1-f)/2.

Medlyn et al. (2002) found relatively low influence ("a slight effect") of α and θ. They also say that Kc, Ko and Γ* "are thought to be intrinsic properties of the Rubisco enzyme and are generally assumed constant among species".

source
PlantBiophysics.FvcbRawType

Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981). Direct implementation of the model.

Parameters

  • Tᵣ: the reference temperature (°C) at which other parameters were measured
  • VcMaxRef: maximum rate of Rubisco activity ($μmol\ m^{-2}\ s^{-1}$)
  • JMaxRef: potential rate of electron transport ($μmol\ m^{-2}\ s^{-1}$)
  • RdRef: mitochondrial respiration in the light at reference temperature ($μmol\ m^{-2}\ s^{-1}$)
  • TPURef: triose phosphate utilization-limited photosynthesis rate ($μmol\ m^{-2}\ s^{-1}$)
  • Eₐᵣ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for Rd.
  • O₂: intercellular dioxygen concentration ($ppm$)
  • Eₐⱼ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for JMax.
  • Hdⱼ: rate of decrease of the function above the optimum (also called EDVJ) for JMax.
  • Δₛⱼ: entropy factor for JMax.
  • Eₐᵥ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for VcMax.
  • Hdᵥ: rate of decrease of the function above the optimum (also called EDVC) for VcMax.
  • Δₛᵥ: entropy factor for VcMax.
  • α: quantum yield of electron transport ($mol_e\ mol^{-1}_{quanta}$). See also eq. 4 of

Medlyn et al. (2002), equation 9.16 from von Caemmerer et al. (2009) ((1-f)/2) and its implementation in get_J

  • θ: determines the curvature of the light response curve for J~aPPFD. See also eq. 4 of

Medlyn et al. (2002) and its implementation in get_J

See also

References

Caemmerer, S. von, et G. D. Farquhar. 1981. « Some Relationships between the Biochemistry of Photosynthesis and the Gas Exchange of Leaves ». Planta 153 (4): 376‑87. https://doi.org/10.1007/BF00384257.

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Examples

Get the fieldnames:
fieldnames(FvcbRaw)
# Using default values for the model:
A = FvcbRaw()

A.Eₐᵥ
source
PlantBiophysics.LightIgnoreType

Ignore model for light interception, see here. Make the mesh invisible, and not computed. Can save a lot of time for the computations when there are components types that are not visible anyway (e.g. inside others).

source
PlantBiophysics.MedlynType

Medlyn et al. (2011) stomatal conductance model for CO₂.

Arguments

  • g0: intercept.
  • g1: slope.
  • gs_min = 0.001: residual conductance. We consider the residual conductance being different from g0 because in practice g0 can be negative when fitting real-world data.

Examples

using PlantMeteo, PlantSimEngine, PlantBiophysics
meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

leaf =
    ModelList(
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (A = A, Cₛ = 380.0, Dₗ = meteo.VPD)
    )
run!(leaf,meteo)

Note that we use VPD as an approximation of Dₗ here because we don't have the leaf temperature (i.e. Dₗ = VPD when Tₗ = T).

References

Medlyn, Belinda E., Remko A. Duursma, Derek Eamus, David S. Ellsworth, I. Colin Prentice, Craig V. M. Barton, Kristine Y. Crous, Paolo De Angelis, Michael Freeman, et Lisa Wingate.

  1. « Reconciling the optimal and empirical approaches to modelling stomatal conductance ».

Global Change Biology 17 (6): 2134‑44. https://doi.org/10.1111/j.1365-2486.2010.02375.x.

source
PlantBiophysics.MonteithType

Struct to hold parameter and values for the energy model close to the one in Monteith and Unsworth (2013)

Arguments

  • aₛₕ = 2: number of faces of the object that exchange sensible heat fluxes
  • aₛᵥ = 1: number of faces of the object that exchange latent heat fluxes (hypostomatous => 1)
  • ε = 0.955: emissivity of the object
  • maxiter = 10: maximal number of iterations allowed to close the energy balance
  • ΔT = 0.01 (°C): maximum difference in object temperature between two iterations to consider convergence

Examples

energy_model = Monteith() # a leaf in an illuminated chamber
source
PlantBiophysics.TranslucentType

Translucent model for light interception, see here. This model is not yet implemented in PlantBiophysics, so it takes the values from the MTG nodes, which means all nodes with the model should provide the necessary attributes:

  • Ra_SW_f the absorbed flux of atmospheric radiation in the short wave bandwidth (PAR+NIR), in W m[object]⁻² (== J m[object]⁻² s⁻¹).
  • Ra_PAR_f the absorbed flux of atmospheric radiation in the PAR bandwidth, in W m[object]⁻² (== J m[object]⁻² s⁻¹).
  • sky_fraction the sky fraction seen by the the node, in [0, 1].
source
PlantBiophysics.Fvcb_net_assimiliationMethod
Fvcb_net_assimiliation(Cᵢ,Vⱼ,Γˢ,VcMax,Km,Rd)

Net assimilation following the Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981)

source
PlantBiophysics.arrheniusFunction
arrhenius(A,Eₐ,Tₖ,Tᵣₖ,R = PlantMeteo.Constants().R)

The Arrhenius function for dependence of the rate constant of a chemical reaction.

Arguments

  • A: pre-exponential factor, a constant for each chemical reaction
  • Eₐ: activation energy for the reaction ($J\ mol^{-1}$)
  • Tₖ: temperature (Kelvin)
  • Tᵣₖ: reference temperature (Kelvin) at which A was measured
  • R: universal gas constant ($J\ mol^{-1}\ K^{-1}$)

Examples

using PlantBiophysics, PlantMeteo
# Importing physical constants
constants = PlantMeteo.Constants()
# Using default values for the model:
A = Fvcb()

# Computing Jmax:
arrhenius(A.JMaxRef,A.Eₐⱼ,28.0-constants.K₀,A.Tᵣ-constants.K₀,constants.R)
# ! Warning: temperatures must be given in Kelvin

# Computing Vcmax:
arrhenius(A.VcMaxRef,A.Eₐᵥ,28.0-constants.K₀,A.Tᵣ-constants.K₀,constants.R)
source
PlantBiophysics.arrheniusFunction
arrhenius(A,Eₐ,Tₖ,Tᵣₖ,Hd,Δₛ,R = Constants().R)

The Arrhenius function for dependence of the rate constant of a chemical reaction, modified following equation (17) from Medlyn et al. (2002) to consider the negative effect of very high temperatures.

Arguments

  • A: the pre-exponential factor, a constant for each chemical reaction
  • Eₐ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise

of the function (Ha in the equation of Medlyn et al. (2002))

  • Tₖ: current temperature (Kelvin)
  • Tᵣₖ: reference temperature (Kelvin) at which A was measured
  • Hd: rate of decrease of the function above the optimum (called EDVJ in

MAESPA and plantecophys)

  • Δₛ: entropy factor
  • R: is the universal gas constant ($J\ mol^{-1}\ K^{-1}$)

References

Medlyn, B. E., E. Dreyer, D. Ellsworth, M. Forstreuter, P. C. Harley, M. U. F. Kirschbaum, X. Le Roux, et al. 2002. « Temperature response of parameters of a biochemically based model of photosynthesis. II. A review of experimental data ». Plant, Cell & Environment 25 (9): 1167‑79. https://doi.org/10.1046/j.1365-3040.2002.00891.x.

Examples

using PlantBiophysics, PlantMeteo
# Importing physical constants
constants = Constants()
# Using default values for the model:
A = Fvcb()

# Computing Jmax:
PlantBiophysics.arrhenius(A.JMaxRef,A.Eₐⱼ,28.0-constants.K₀,A.Tᵣ-constants.K₀,A.Hdⱼ,A.Δₛⱼ)
# ! Warning: temperatures must be given in Kelvin

# Computing Vcmax:
PlantBiophysics.arrhenius(A.VcMaxRef,A.Eₐᵥ,28.0-constants.K₀,A.Tᵣ-constants.K₀,A.Hdᵥ,A.Δₛᵥ)
source
PlantBiophysics.black_bodyMethod
black_body(T, K₀, σ)
black_body(T)

Thermal infrared, i.e. longwave radiation emitted from a black body at temperature T.

Note

K₀ and σ are taken from PlantMeteo.Constants if not provided.

source
PlantBiophysics.gbh_to_gbwFunction
gbh_to_gbw(gbh, Gbₕ_to_Gbₕ₂ₒ = PlantMeteo.Constants().Gbₕ_to_Gbₕ₂ₒ)
gbw_to_gbh(gbh, Gbₕ_to_Gbₕ₂ₒ = PlantMeteo.Constants().Gbₕ_to_Gbₕ₂ₒ)

Boundary layer conductance for water vapor from boundary layer conductance for heat.

Arguments

  • gbh (m s-1): boundary layer conductance for heat under mixed convection.
  • Gbₕ_to_Gbₕ₂ₒ: conversion factor.

Note

Gbₕ is the sum of free and forced convection. See gbₕ_free and gbₕ_forced.

source
PlantBiophysics.gbₕ_forcedMethod
gbₕ_forced(Wind,d)

Boundary layer conductance for heat under forced convection (m s-1). See eq. E1 from Leuning et al. (1995) for more details.

Arguments

  • Wind (m s-1): wind speed
  • d (m): characteristic dimension, e.g. leaf width (see eq. 10.9 from Monteith and Unsworth, 2013).

Notes

d is the minimal dimension of the surface of an object in contact with the air.

References

Leuning, R., F. M. Kelliher, DGG de Pury, et E.-D. SCHULZE. 1995. « Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

source
PlantBiophysics.gbₕ_freeFunction
gbₕ_free(Tₐ,Tₗ,d,Dₕ₀)
gbₕ_free(Tₐ,Tₗ,d)

Leaf boundary layer conductance for heat under free convection (m s-1).

Arguments

  • Tₐ (°C): air temperature
  • Tₗ (°C): leaf temperature
  • d (m): characteristic dimension, e.g. leaf width (see eq. 10.9 from Monteith and Unsworth, 2013).
  • Dₕ₀ = 21.5e-6: molecular diffusivity for heat at base temperature. Use value from

PlantMeteo.Constants if not provided.

Note

R and Dₕ₀ can be found using PlantMeteo.Constants. To transform in $mol\ m^{-2}\ s^{-1}$, use ms_to_mol.

References

Leuning, R., F. M. Kelliher, DGG de Pury, et E.-D. SCHULZE. 1995. « Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

Monteith, John, et Mike Unsworth. 2013. Principles of environmental physics: plants, animals, and the atmosphere. Academic Press. Paragraph 10.1.3, eq. 10.9.

source
PlantBiophysics.get_CᵢᵥMethod

Analytic resolution of Cᵢ when the RuBisCo activity is limiting ($μmol\ mol^{-1}$)

Arguments

  • VcMAX: maximum rate of RuBisCo activity($μmol\ m^{-2}\ s^{-1}$)
  • Γˢ: CO2 compensation point $Γ^⋆$ ($μmol\ mol^{-1}$)
  • Cₛ: Air CO₂ concentration at the leaf surface ($μmol\ mol^{-1}$)
  • Rd: day respiration ($μmol\ m^{-2}\ s^{-1}$)
  • g0: residual stomatal conductance ($μmol\ m^{-2}\ s^{-1}$)
  • st_closure: stomatal conductance term computed from a given implementation of a Gs model,

e.g. Medlyn.

  • Km: effective Michaelis–Menten coefficient for CO2 ($μ mol\ mol^{-1}$)
source
PlantBiophysics.get_CᵢⱼMethod

Analytic resolution of Cᵢ when the rate of electron transport is limiting ($μmol\ mol^{-1}$)

Arguments

  • Vⱼ: RuBP regeneration (J/4.0, $μmol\ m^{-2}\ s^{-1}$)
  • Γˢ: CO2 compensation point $Γ^⋆$ ($μmol\ mol^{-1}$)
  • Cₛ: Air CO₂ concentration at the leaf surface ($μmol\ mol^{-1}$)
  • Rd: day respiration ($μmol\ m^{-2}\ s^{-1}$)
  • g0: residual stomatal conductance ($μmol\ m^{-2}\ s^{-1}$)
  • st_closure: stomatal conductance term computed from a given implementation of a Gs model,

e.g. Medlyn.

References

Duursma, R. A., et B. E. Medlyn. 2012. « MAESPA: a model to study interactions between water limitation, environmental drivers and vegetation function at tree and stand levels, with an example application to [CO2] × drought interactions ». Geoscientific Model Development 5 (4): 919‑40. https://doi.org/10.5194/gmd-5-919-2012.

Wang and Leuning, 1998

source
PlantBiophysics.get_DₕFunction
get_Dₕ(T,Dₕ₀)
get_Dₕ(T)

Dₕ -molecular diffusivity for heat at base temperature- from Dₕ₀ (corrected by temperature). See Monteith and Unsworth (2013, eq. 3.10).

Arguments

  • Tₐ (°C): temperature
  • Dₕ₀: molecular diffusivity for heat at base temperature. Use value from PlantMeteo.Constants

if not provided.

References

Monteith, John, et Mike Unsworth. 2013. Principles of environmental physics: plants, animals, and the atmosphere. Academic Press. Paragraph 10.1.3.

source
PlantBiophysics.get_JMethod

Rate of electron transport J ($μmol\ m^{-2}\ s^{-1}$), computed using the smaller root of the quadratic equation (eq. 4 from Medlyn et al., 2002):

θ * J² - (α * aPPFD + JMax) * J + α * aPPFD * JMax

NB: we use the smaller root because considering the range of values for θ and α (quite stable), and aPPFD and JMax, the function always tends to JMax with high aPPFD with the smaller root (behavior we are searching), and the opposite with the larger root.

Returns

A tuple with (A, Gₛ, Cᵢ):

  • A: carbon assimilation (μmol m-2 s-1)
  • Gₛ: stomatal conductance (mol m-2 s-1)
  • Cᵢ: intercellular CO₂ concentration (ppm)

Arguments

  • aPPFD: absorbed photon irradiance ($μmol_{quanta}\ m^{-2}\ s^{-1}$)
  • α: quantum yield of electron transport ($mol_e\ mol^{-1}_{quanta}$)
  • JMax: maximum rate of electron transport ($μmol\ m^{-2}\ s^{-1}$)
  • θ: determines the shape of the non-rectangular hyperbola (-)

References

Medlyn, B. E., E. Dreyer, D. Ellsworth, M. Forstreuter, P. C. Harley, M. U. F. Kirschbaum, X. Le Roux, et al. 2002. « Temperature response of parameters of a biochemically based model of photosynthesis. II. A review of experimental data ». Plant, Cell & Environment 25 (9): 1167‑79. https://doi.org/10.1046/j.1365-3040.2002.00891.x.

Von Caemmerer, Susanna. 2000. Biochemical models of leaf photosynthesis. Csiro publishing.

Examples

# Using default values for the model:
julia> A = Fvcb()
Fvcb{Float64}(25.0, 200.0, 250.0, 0.6, 9999.0, 46390.0, 210.0, 29680.0, 200000.0, 631.88, 58550.0, 200000.0, 629.26, 0.425, 0.7)

julia> PlantBiophysics.get_J(1500, A.JMaxRef, A.α, A.θ)
216.5715752671342
source
PlantBiophysics.get_kmFunction

Compute the effective Michaelis–Menten coefficient for CO2 $Km$ ($μ mol\ mol^{-1}$) according to Medlyn et al. (2002), equations (5) and (6).

References

Medlyn, B. E., E. Dreyer, D. Ellsworth, M. Forstreuter, P. C. Harley, M. U. F. Kirschbaum, X. Le Roux, et al. 2002. « Temperature response of parameters of a biochemically based model of photosynthesis. II. A review of experimental data ». Plant, Cell & Environment 25 (9): 1167‑79. https://doi.org/10.1046/j.1365-3040.2002.00891.x.

Examples

# computing the temperature dependence of γˢ:
get_km(28,25,210.0)
source
PlantBiophysics.grey_bodyMethod

Thermal infrared, i.e. longwave radiation emitted from an object at temperature T.

  • T: temperature of the object in Celsius degree
  • ε object emissivity (not to confuse with ε the

ratio of molecular weights from PlantMeteo.Constants). A typical value for a leaf is 0.955.

Note

K₀ and σ are taken from PlantMeteo.Constants if not provided.

Examples

# Thermal infrared radiation of water at 25 °C:
grey_body(25.0, 0.96)
source
PlantBiophysics.gs_closureFunction
gs_closure(::Medlyn, models, status, meteo, constants=nothing, extra=nothing)

Stomatal closure for CO₂ according to Medlyn et al. (2011). Carefull, this is just a part of the computation of the stomatal conductance.

The result of this function is then used as:

gs_mod = gs_closure(leaf,meteo)

# And then stomatal conductance (μmol m-2 s-1) calling [`stomatal_conductance`](@ref):
Gₛ = leaf.stomatal_conductance.g0 + gs_mod * leaf.status.A

Arguments

  • ::Medlyn: an instance of the Medlyn model type
  • models::ModelList: A ModelList struct holding the parameters for the models.
  • status: A status struct holding the variables for the models.
  • meteo: meteorology structure, see Atmosphere. Is not used in this model.
  • constants: A constants struct holding the constants for the models. Is not used in this model.
  • extra: A struct holding the extra variables for the models. Is not used in this model.

Details

Use variables() on Medlyn to get the variables that must be instantiated in the ModelList struct.

Notes

  • Cₛ is used instead of Cₐ because Gₛ is between the surface and the intercellular space. The conductance

between the atmosphere and the surface is accounted for using the boundary layer conductance (Gbc in Monteith). Medlyn et al. (2011) uses Cₐ in their paper because they relate their models to the measurements made at leaf level, with a well-mixed chamber whereCₛ ≈ Cₐ.

  • Dₗ is forced to be >= 1e-9 because it is used in a squared root. It is prefectly acceptable to

get a negative Dₗ when leaves are re-hydrating from air. Cloud forests are the perfect example. See e.g.: Guzmán‐Delgado, P, Laca, E, Zwieniecki, MA. Unravelling foliar water uptake pathways: The contribution of stomata and the cuticle. Plant Cell Environ. 2021; 1– 13. https://doi.org/10.1111/pce.14041

Examples

meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

leaf =
    ModelList(
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Cₛ = 380.0, Dₗ = meteo.VPD)
    )

gs_mod = gs_closure(leaf, meteo)

A = 20 # example assimilation (μmol m-2 s-1)
Gs = leaf.stomatal_conductance.g0 + gs_mod * A

# Or more directly using `run!()`:

leaf =
    ModelList(
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (A = A, Cₛ = 380.0, Dₗ = meteo.VPD)
    )
run!(leaf,meteo)

Note that we use VPD as an approximation of Dₗ here because we don't have the leaf temperature (i.e. Dₗ = VPD when Tₗ = T).

References

Medlyn, Belinda E., Remko A. Duursma, Derek Eamus, David S. Ellsworth, I. Colin Prentice, Craig V. M. Barton, Kristine Y. Crous, Paolo De Angelis, Michael Freeman, et Lisa Wingate.

  1. « Reconciling the optimal and empirical approaches to modelling stomatal conductance ».

Global Change Biology 17 (6): 2134‑44. https://doi.org/10.1111/j.1365-2486.2010.02375.x.

source
PlantBiophysics.gs_closureFunction

Constant stomatal closure. Usually called from a photosynthesis model.

Note

meteo is just declared here for compatibility with other formats of calls.

source
PlantBiophysics.gsc_to_gswFunction
gsc_to_gsw(Gₛ, Gsc_to_Gsw = PlantMeteo.Constants().Gsc_to_Gsw)

Conversion of a stomatal conductance for CO₂ into stomatal conductance for H₂O.

source
PlantBiophysics.gsw_to_gscFunction
gsw_to_gsc(Gₛ, Gsc_to_Gsw = PlantMeteo.Constants().Gsc_to_Gsw)

Conversion of a stomatal conductance for H₂O into stomatal conductance for CO₂.

source
PlantBiophysics.instantiateMethod
instantiate(x)

Instantiate a model given its parameter names, considering that parameter names can be different compared to the model fields (used to insure compatibility with Archimed).

source
PlantBiophysics.is_modelMethod
is_model(model)

Check if a model object has the"Group" and "Type" keys as the first level of a Dict type object. But the function is generic as long as the input struct has a keys() method.

Examples

models = read_model("path_to_a_model_file.yaml")
is_model(models)
source
PlantBiophysics.latent_heatMethod
latent_heat(Rn, VPD, γˢ, Rbₕ, Δ, ρ, aₛₕ, Cₚ)
latent_heat(Rn, VPD, γˢ, Rbₕ, Δ, ρ, aₛₕ)

λE -the latent heat flux (W m-2)- using the Monteith and Unsworth (2013) definition corrected by Schymanski et al. (2017), eq.22.

  • Rn (W m-2): net radiation. Carefull: not the isothermal net radiation
  • VPD (kPa): air vapor pressure deficit
  • γˢ (kPa K−1): apparent value of psychrometer constant (see PlantMeteo.γ_star)
  • Rbₕ (s m-1): resistance for heat transfer by convection, i.e. resistance to sensible heat
  • Δ (KPa K-1): rate of change of saturation vapor pressure with temperature (see PlantMeteo.e_sat_slope)
  • ρ (kg m-3): air density of moist air.
  • aₛₕ (1,2): number of sides that exchange energy for heat (2 for leaves)
  • Cₚ (J K-1 kg-1): specific heat of air for constant pressure

References

Monteith, J. and Unsworth, M., 2013. Principles of environmental physics: plants, animals, and the atmosphere. Academic Press. See eq. 13.33.

Schymanski et al. (2017), Leaf-scale experiments reveal an important omission in the Penman–Monteith equation, Hydrology and Earth System Sciences. DOI: https://doi.org/10.5194/hess-21-685-2017. See equ. 22.

Examples

Tₐ = 20.0 ; P = 100.0 ;
ρ = air_density(Tₐ, P) # in kg m-3
Δ = e_sat_slope(Tₐ)

latent_heat(300.0, 2.0, 0.1461683, 50.0, Δ, ρ, 2.0)
source
PlantBiophysics.mol_to_msMethod
ms_to_mol(G,T,P,R,K₀)
ms_to_mol(G,T,P)

Conversion of a conductance G from $mol\ m^{-2}\ s^{-1}$ to $m\ s^{-1}$.

Arguments

  • G ($m\ s^{-1}$): conductance
  • T (°C): air temperature
  • P (kPa): air pressure
  • R ($J\ mol^{-1}\ K^{-1}$): universal gas constant.
  • K₀ (°C): absolute zero

See also

ms_to_mol for the inverse process.

source
PlantBiophysics.ms_to_molMethod
ms_to_mol(G,T,P,R,K₀)
ms_to_mol(G,T,P)

Conversion of a conductance G from $m\ s^{-1}$ to $mol\ m^{-2}\ s^{-1}$.

Arguments

  • G ($m\ s^{-1}$): conductance
  • T (°C): air temperature
  • P (kPa): air pressure
  • R ($J\ mol^{-1}\ K^{-1}$): universal gas constant.
  • K₀ (°C): absolute zero

See also

mol_to_ms for the inverse process.

source
PlantBiophysics.net_longwave_radiationMethod
net_longwave_radiation(T₁,T₂,ε₁,ε₂,F₁,K₀,σ)
net_longwave_radiation(T₁,T₂,ε₁,ε₂,F₁)

Net longwave radiation fluxes (i.e. thermal radiation, W m-2) between an object and another. The object of interest is at temperature T₁ and has an emissivity ε₁, and the object with which it exchanges energy is at temperature T₂ and has an emissivity ε₂.

If the result is positive, then the object of interest gain energy.

Arguments

  • T₁ (Celsius degree): temperature of the target object (object 1)
  • T₂ (Celsius degree): temperature of the object with which there is potential exchange (object 2)
  • ε₁: object 1 emissivity
  • ε₂: object 2 emissivity
  • F₁: view factor (0-1), i.e. visible fraction of object 2 from object 1 (see note)
  • K₀: absolute zero (°C)
  • σ ($W\ m^{-2}\ K^{-4}$) Stefan-Boltzmann constant

Note

F₁, the view factor (also called shape factor) is a coefficient applied to the semi-hemisphere field of view of object 1 that "sees" object 2. E.g. a leaf can be viewed as a plane. If one side of the leaf sees only object 2 in its field of view (e.g. the sky), then F₁ = 1. Then the net longwave radiation flux for this part of the leaf is multiplied by its actual surface to get the exchange. Note that we apply reciprocity between the two objects for the view factor (they have the same value), i.e.: A₁F₁₂ = A₂F₂₁.

Then, if we take a leaf as object 1, and the sky as object 2, the visible fraction of sky viewed by the leaf would be:

  • 0.5 if the leaf is on top of the canopy, i.e. the upper side of the leaf sees the sky,

the side bellow sees other leaves and the soil.

  • between 0 and 0.5 if it is within the canopy and partly shaded by other objects.

Note that A₁ for a leaf is twice its common used leaf area, because A₁ is the total leaf area of the object that exchange energy.

# Net thermal radiation fluxes between a leaf and the sky considering the leaf at the top of
# the canopy:
Tₗ = 25.0 ; Tₐ = 20.0
ε₁ = 0.955 ; ε₂ = 1.0
Ra_LW_f = net_longwave_radiation(Tₗ,Tₐ,ε₁,ε₂,1.0)
Ra_LW_f

# Ra_LW_f is the net longwave radiation flux between the leaf and the atmosphere per surface area.
# To get the actual net longwave radiation flux we need to multiply by the surface of the
# leaf, e.g. for a leaf of 2cm²:
leaf_area = 2e-4 # in m²
Ra_LW_f * leaf_area

# The leaf lose ~0.0055 W towards the atmosphere.

References

Cengel, Y, et Transfer Mass Heat. 2003. A practical approach. New York, NY, USA: McGraw-Hill.

source
PlantBiophysics.read_licor6400Method
read_licor6400(file; abs=0.85)

Import Licor6400 data (such as Medlyn 2001 data) with the units and names corresponding to the ones used in PlantBiophysics.jl.

Arguments

  • file: a string or a vector of strings containing the path to the file(s) to read.
  • abs=0.85: the absorptance of the leaf. Default is 0.85 (von Caemmerer et al., 2009).
source
PlantBiophysics.read_modelMethod
read_model(file)

Read a model file. The model file holds the choice and the parameterization of the models.

Arguments

  • file::String: path to a model file

Examples

file = joinpath(dirname(dirname(pathof(PlantBiophysics))), "test", "inputs", "models", "plant_coffee.yml")
models = read_model(file)
source
PlantBiophysics.read_walzMethod
read_walz(file; abs=0.85)

Import a Walz GFS-3000 output file, perform variables conversion and rename according to package conventions.

Arguments

  • file: a string or a vector of strings containing the path to the file(s) to read.
  • abs=0.85: the absorptance of the leaf. Default is 0.85 (von Caemmerer et al., 2009).

Returns

A DataFrame containing the data read and transformed from the file(s). The units are the same than in the Walz output file, except for:

  • Dₗ: kPa
  • Rh: fraction (0-1)
  • VPD: kPa
  • Gₛ: mol[CO₂] m⁻² s⁻¹

Examples

Reading one file:

file = joinpath(dirname(dirname(pathof(PlantBiophysics))),"test","inputs","data","P1F20129.csv")
read_walz(file)

We can also read multiple files at once:

files = readdir(joinpath(dirname(dirname(pathof(PlantBiophysics))),"test","inputs","data"), join = true)
df = read_walz(files)

In this case, the source of the data is added as the source columns into the DataFrame:

df.source
source
PlantBiophysics.sensible_heatMethod
sensible_heat(Rn, VPD, γˢ, Rbₕ, Δ, ρ, aₛₕ, Cₚ)
sensible_heat(Rn, VPD, γˢ, Rbₕ, Δ, ρ, aₛₕ)

H -the sensible heat flux (W m-2)- using the Monteith and Unsworth (2013) definition corrected by Schymanski et al. (2017), eq.22.

  • Rn (W m-2): net radiation. Carefull: not the isothermal net radiation
  • VPD (kPa): air vapor pressure deficit
  • γˢ (kPa K−1): apparent value of psychrometer constant (see PlantMeteo.γ_star)
  • Rbₕ (s m-1): resistance for heat transfer by convection, i.e. resistance to sensible heat
  • Δ (KPa K-1): rate of change of saturation vapor pressure with temperature (see PlantMeteo.e_sat_slope)
  • ρ (kg m-3): air density of moist air.
  • aₛₕ (1,2): number of sides that exchange energy for heat (2 for leaves)
  • Cₚ (J K-1 kg-1): specific heat of air for constant pressure

References

Monteith, J. and Unsworth, M., 2013. Principles of environmental physics: plants, animals, and the atmosphere. Academic Press. See eq. 13.33.

Schymanski et al. (2017), Leaf-scale experiments reveal an important omission in the Penman–Monteith equation, Hydrology and Earth System Sciences. DOI: https://doi.org/10.5194/hess-21-685-2017. See equ. 22.

Examples

Tₐ = 20.0 ; P = 100.0 ;
ρ = air_density(Tₐ, P) # in kg m-3
Δ = PlantMeteo.e_sat_slope(Tₐ)

sensible_heat(300.0, 2.0, 0.1461683, 50.0, Δ, ρ, 2.0)
source
PlantBiophysics.Γ_starFunction
Γ_star(Tₖ,Tᵣₖ,R = PlantMeteo.Constants().R)

CO₂ compensation point $Γ^⋆$ ($μ mol\ mol^{-1}$) according to equation (12) from Medlyn et al. (2002).

$Γ^⋆$ is the [CO₂] at which oxygenation proceeds at

twice the rate of carboxylation causing photosynthetic uptake of CO2 to be exactly compensated by photorespiratory CO₂ release (Sharkey et al., 2007).

Notes

Could be replaced by equation (38) from Farquhar et al. (1980), but Medlyn et al. (2002) states that $Γ^⋆$ as a relatively low effect on the model outputs.

Arguments

  • Tₖ (Kelvin): current temperature
  • Tᵣₖ (Kelvin): reference temperature at which A was measured
  • R ($J\ mol^{-1}\ K^{-1}$): is the universal gas constant

Examples

using PlantBiophysics, PlantMeteo
# Importing the physical constants:
constants = Constants()
# computing the temperature dependence of γˢ:
Γ_star(28-constants.K₀,25-constants.K₀,constants.R)

References

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Medlyn, B. E., E. Dreyer, D. Ellsworth, M. Forstreuter, P. C. Harley, M. U. F. Kirschbaum, X. Le Roux, et al. 2002. « Temperature response of parameters of a biochemically based model of photosynthesis. II. A review of experimental data ». Plant, Cell & Environment 25 (9): 1167‑79. https://doi.org/10.1046/j.1365-3040.2002.00891.x.

Sharkey, Thomas D., Carl J. Bernacchi, Graham D. Farquhar, et Eric L. Singsaas. 2007. « Fitting Photosynthetic Carbon Dioxide Response Curves for C3 Leaves ». Plant, Cell & Environment 30 (9): 1035‑40. https://doi.org/10.1111/j.1365-3040.2007.01710.x.

source

Un-exported functions

PlantBiophysics.f_ms_to_molMethod

Conversion factor between conductance in $m\ s^{-1}$ to $mol\ m^{-2}\ s^{-1}$.

Arguments

  • T (°C): air temperature
  • P (kPa): air pressure
  • R ($J\ mol^{-1}\ K^{-1}$): universal gas constant.
  • K₀ (°C): absolute zero
source
PlantBiophysics.locf!Method
locf!(var)

Last observation carried forward (LOCF) iterates forwards var and fills missing data with the last existing observation.

This function is heavily inspired (i.e. copied) by the function locf in the package Impute (MIT licence). See here for more details.

source
PlantBiophysics.negative_rootMethod

Negative root of a quadratic equation, but returns 0 if Δ is negative. Careful, this is not right mathematically, but biologically OK because used in the computation of Cᵢ (gives A = 0 in this case).

source
PlantBiophysics.positive_rootMethod

Positive root of a quadratic equation, but returns 0 if Δ is negative. Careful, this is not right mathematically, but biologically OK because used in the computation of Cᵢ (gives A = 0 in this case).

source
PlantBiophysics.γ_starMethod

γstar(γ, ash, a_s, rbv, Rsᵥ, Rbₕ)

γ∗, the apparent value of psychrometer constant (kPa K−1).

Arguments

  • γ (kPa K−1): psychrometer constant
  • aₛₕ (1,2): number of faces exchanging heat fluxes (see Schymanski et al., 2017)
  • aₛᵥ (1,2): number of faces exchanging water fluxes (see Schymanski et al., 2017)
  • Rbᵥ (s m-1): boundary layer resistance to water vapor
  • Rsᵥ (s m-1): stomatal resistance to water vapor
  • Rbₕ (s m-1): boundary layer resistance to heat

Note

Using the corrigendum from Schymanski et al. (2017) in here so the definition of latent_heat remains generic.

Not to be confused with Γ_star the CO₂ compensation point.

References

Monteith, John L., et Mike H. Unsworth. 2013. « Chapter 13 - Steady-State Heat Balance: (i) Water Surfaces, Soil, and Vegetation ». In Principles of Environmental Physics (Fourth Edition), edited by John L. Monteith et Mike H. Unsworth, 217‑47. Boston: Academic Press.

Schymanski, Stanislaus J., et Dani Or. 2017. Leaf-Scale Experiments Reveal an Important Omission in the Penman–Monteith Equation ». Hydrology and Earth System Sciences 21 (2): 685‑706. https://doi.org/10.5194/hess-21-685-2017.

source
PlantBiophysics.λE_to_EFunction
λE_to_E(λE, λ, Mₕ₂ₒ=PlantMeteo.Constants().Mₕ₂ₒ)
E_to_λE(E, λ, Mₕ₂ₒ=PlantMeteo.Constants().Mₕ₂ₒ)

Conversion from latent heat (W m-2) to evaporation (mol[H₂O] m-2 s-1) or the opposite (E_to_λE).

Arguments

  • λE: latent heat flux (W m-2)
  • E: water evaporation (mol[H₂O] m-2 s-1)
  • λ (J kg-1): latent heat of vaporization
  • Mₕ₂ₒ = 18.0e-3 (kg mol-1): Molar mass for water.

Note

λ can be computed using:

λ = latent_heat_vaporization(T, constants.λ₀)

It is also directly available from the Atmosphere structure, and by extention in Weather.

To convert E from mol[H₂O] m-2 s-1 to mm s-1 you can simply do:

E_mms = E_mol / constants.Mₕ₂ₒ

mm[H₂O] s-1 is equivalent to kg[H₂O] m-2 s-1, wich is equivalent to l[H₂O] m-2 s-1.

source
PlantSimEngine.fitMethod
fit(::Type{Beer}, df; J_to_umol=PlantMeteo.Constants().J_to_umol)

Optimize k, the coefficient of the Beer-Lambert law of light extinction.

Arguments

  • df: a DataFrame with columns RiPARf (Incoming light flux in the PAR, W m⁻²),

aPPFD (μmol m⁻² s⁻¹) and LAI (m² m⁻²), where each row is an observation. The column names should match exactly.

Examples

using PlantSimEngine, PlantBiophysics, DataFrames, PlantMeteo

# Defining dummy data:
df = DataFrame(
    Ri_PAR_f = [200.0, 250.0, 300.0], 
    aPPFD = [548.4, 685.5, 822.6], 
    LAI = [1.0, 1.0, 1.0],
    T = [20.0, 20.0, 20.0],
    Rh = [0.5, 0.5, 0.5],
    Wind = [10.0, 10.0, 10.0],
)

# Fit the parameters values:
k = fit(Beer, df)

# Re-simulating aPPFD using the newly fitted parameters:
w = Weather(df)
leaf = ModelList(
        Beer(k.k),
        status = (LAI = df.LAI,)
    )
run!(leaf, w)

leaf
source
PlantSimEngine.fitMethod
fit(::Type{Medlyn}, df)

Optimize the parameters of the Medlyn model. Note that here Gₛ is stomatal conductance for CO2, not H2O.

Careful, g0 and g1 are fitted on CO₂ stomatal conductance, and so are given for CO₂, not H₂O. If you want to compare the values with the ones usually given in the literature, you need to multiply g0 and g1 by 1.57.

Arguments

  • df: a DataFrame with columns A, Dₗ, Cₐ and Gₛ, where each row is an observation. The column

names should match exactly.

Note that Dₗ is the leaf to air vapour pressure deficit, computed as follows:

Dₗ = PlantMeteo.e_sat(Tₗ) - PlantMeteo.e_sat(T) * Rh

with Tₗ the leaf temperature, T the air temperature and Rh the air relative humidity. If Tₗ is not available, you can use Tₗ = T as a crude approximation, which gives Dₗ = VPD.

Examples

using PlantBiophysics, PlantMeteo, Plots, DataFrames

file = joinpath(dirname(dirname(pathof(PlantBiophysics))),"test","inputs","data","P1F20129.csv")
df = read_walz(file)
# Removing the CO2 and ligth Curve, we fit the parameters on the Rh curve:
filter!(x -> x.curve != "ligth Curve" && x.curve != "CO2 Curve", df)

# Fit the parameters values:
g0, g1 = fit(Medlyn, df)

# Re-simulating Gₛ using the newly fitted parameters:
w = Weather(select(df, :T, :P, :Rh, :Cₐ, :VPD, :T => (x -> 10) => :Wind))
leaf = ModelList(
        stomatal_conductance = Medlyn(g0, g1),
        status = (A = df.A, Cₛ = df.Cₐ, Dₗ = df.Dₗ)
    )
run!(leaf, w)

# Visualising the results:
gsAvpd = PlantBiophysics.GsADₗ(g0, g1, df.Gₛ, df.Dₗ, df.A, df.Cₐ, leaf[:Gₛ])
plot(gsAvpd,leg=:bottomright)
# As in [`Medlyn`](@ref) reference paper, linear regression is also plotted.
source
PlantSimEngine.fitMethod
fit(
    ::Type{Fvcb}, df; 
    Tᵣ = nothing, 
    VcMaxRef = 0.0, JMaxRef = 0.0, RdRef = 0.0, TPURef = 0.0, 
    VcMaxRef_bound=[0.0, Inf], JMaxRef_bound=[0.0, Inf], RdRef_bound=[0.0, Inf], TPURef_bound=[0.0, Inf],
    verbose = true,
    Eₐᵣ=46390.0, O₂=210.0, Eₐⱼ=29680.0, Hdⱼ=200000.0, Δₛⱼ=631.88, Eₐᵥ=58550.0, Hdᵥ=200000.0, Δₛᵥ=629.26, α=0.425, θ=0.7
)

Optimize the parameters of the Fvcb model. Also works for FvcbIter.

Arguments

  • df: a DataFrame with columns A, aPPFD, Tₗ and Cᵢ, where each row is an observation. The column

names should match exactly

  • Tᵣ: reference temperature for the optimized parameter values. If not provided, use the average Tₗ.
  • VcMaxRef, JMaxRef, RdRef, TPURef: initialisation values for the parameter optimisation
  • VcMaxRefbound, JMaxRefbound, RdRefbound, TPURefbound: boundary values for the parameter optimisation
  • verbose: if true, print the optimisation results
  • Eₐᵣ, O₂, Eₐⱼ, Hdⱼ, Δₛⱼ, Eₐᵥ, Hdᵥ, Δₛᵥ, α, θ: parameters for the FvCB model

FvCB model parameters:

Parameters

  • Tᵣ: the reference temperature (°C) at which other parameters were measured
  • VcMaxRef: maximum rate of Rubisco activity ($μmol\ m^{-2}\ s^{-1}$)
  • JMaxRef: potential rate of electron transport ($μmol\ m^{-2}\ s^{-1}$)
  • RdRef: mitochondrial respiration in the light at reference temperature ($μmol\ m^{-2}\ s^{-1}$)
  • TPURef: triose phosphate utilization-limited photosynthesis rate ($μmol\ m^{-2}\ s^{-1}$)
  • Eₐᵣ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for Rd.
  • O₂: intercellular dioxygen concentration ($ppm$)
  • Eₐⱼ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for JMax.
  • Hdⱼ: rate of decrease of the function above the optimum (also called EDVJ) for JMax.
  • Δₛⱼ: entropy factor for JMax.
  • Eₐᵥ: activation energy ($J\ mol^{-1}$), or the exponential rate of rise for VcMax.
  • Hdᵥ: rate of decrease of the function above the optimum (also called EDVC) for VcMax.
  • Δₛᵥ: entropy factor for VcMax.
  • α: quantum yield of electron transport ($mol_e\ mol^{-1}_{quanta}$). See also eq. 4 of

Medlyn et al. (2002), equation 9.16 from von Caemmerer et al. (2009) ((1-f)/2) and its implementation in get_J

  • θ: determines the curvature of the light response curve for J~aPPFD. See also eq. 4 of

Medlyn et al. (2002) and its implementation in get_J

Note that boundary values are set to [0.0, Inf] by default. You should adapt them to your use case. Note that no boundary can be set using [-Inf, Inf].

Examples

using PlantBiophysics, PlantSimEngine, PlantMeteo, Plots, DataFrames

file = joinpath(dirname(dirname(pathof(PlantBiophysics))),"test","inputs","data","P1F20129.csv")
df = read_walz(file)
# Removing the Rh and light curves for the fitting because temperature varies
filter!(x -> x.curve != "Rh Curve" && x.curve != "ligth Curve", df)

# Fit the parameter values:
VcMaxRef, JMaxRef, RdRef, TPURef = fit(Fvcb, df; Tᵣ = 25.0)
# Note that Tᵣ was set to 25 °C in our response curve. You should adapt its value to what you
# had during the response curves

# Checking the results:
filter!(x -> x.curve == "CO2 Curve", df)

# Sort the DataFrame by :Cᵢ to get ordered data point
sort!(df, :Cᵢ)

# Re-simulating A using the newly fitted parameters:
leaf =
    ModelList(
        photosynthesis = FvcbRaw(VcMaxRef = VcMaxRef, JMaxRef = JMaxRef, RdRef = RdRef, TPURef = TPURef),
        status = (Tₗ = df.Tₗ, aPPFD = df.aPPFD, Cᵢ = df.Cᵢ)
    )
run!(leaf)
df_sim = DataFrame(leaf)

# Visualising the results:
ACi_struct = PlantBiophysics.ACi(VcMaxRef, JMaxRef, RdRef, df.A, df_sim.A, df[:,:Cᵢ], df_sim.Cᵢ)
plot(ACi_struct,leg=:bottomright)

# Note that we can also simulate the results using the full photosynthesis model too (Fvcb):
# Adding the windspeed to simulate the boundary-layer conductance (we put a high value):
df[!, :Wind] .= 10.0

leaf = ModelList(
        photosynthesis = Fvcb(VcMaxRef = VcMaxRef, JMaxRef = JMaxRef, RdRef = RdRef, Tᵣ = 25.0, TPURef = TPURef),
        # stomatal_conductance = ConstantGs(0.0, df[i,:Gₛ]),
        stomatal_conductance = Medlyn(0.03, 12.),
        status = (Tₗ = df.Tₗ, aPPFD = df.aPPFD, Cₛ = df.Cₐ, Dₗ = 0.1)
    )

w = Weather(select(df, :T, :P, :Rh, :Cₐ, :T => (x -> 10) => :Wind))
run!(leaf, w)
df_sim2 = DataFrame(leaf)

# And finally we plot the results:
ACi_struct_full = PlantBiophysics.ACi(VcMaxRef, JMaxRef, RdRef, df.A, df_sim2.A, df[:,:Cᵢ], df_sim2.Cᵢ)
plot(ACi_struct_full,leg=:bottomright)
# Note that the results differ a bit because there are more variables that are re-simulated (e.g. Cᵢ)
source
PlantSimEngine.run!Function
run!(::LightIgnore, models::ModelList, status, meteo,constants = Constants())

Method for when light interception should be explicitely ignored (do nothing).

Arguments

  • ::LightIgnore: an Ignore model.
  • models: a ModelList struct with a missing energy model.
  • status: the status of the model, usually the one from the models (i.e. models.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details
source
PlantSimEngine.run!Function
run!(::Monteith, models, status, meteo, constants=Constants())

Leaf energy balance according to Monteith and Unsworth (2013), and corrigendum from Schymanski et al. (2017). The computation is close to the one from the MAESPA model (Duursma et al., 2012, Vezy et al., 2018) here. The leaf temperature is computed iteratively to close the energy balance using the mass flux (~ Rn - λE).

Arguments

  • ::Monteith: a Monteith model, usually from a model list (i.e. m.energy_balance)
  • models: A ModelList struct holding the parameters for the model with

initialisations for: - Ra_SW_f (W m-2): net shortwave radiation (PAR + NIR). Often computed from a light interception model - sky_fraction (0-2): view factor between the object and the sky for both faces (see details). - d (m): characteristic dimension, e.g. leaf width (see eq. 10.9 from Monteith and Unsworth, 2013).

  • status: the status of the model, usually the model list status (i.e. leaf.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Details

The sky_fraction in the variables is equal to 2 if all the leaf is viewing is sky (e.g. in a controlled chamber), 1 if the leaf is e.g. up on the canopy where the upper side of the leaf sees the sky, and the side bellow sees soil + other leaves that are all considered at the same temperature than the leaf, or less than 1 if it is partly shaded.

Notes

If you want the algorithm to print a message whenever it does not reach convergence, use the debugging mode by executing this in the REPL: ENV["JULIA_DEBUG"] = PlantBiophysics.

More information here.

Examples

meteo = Atmosphere(T = 22.0, Wind = 0.8333, P = 101.325, Rh = 0.4490995)

# Using a constant value for Gs:
leaf = ModelList(
    energy_balance = Monteith(),
    photosynthesis = Fvcb(),
    stomatal_conductance = ConstantGs(0.0, 0.0011),
    status = (Ra_SW_f = 13.747, sky_fraction = 1.0, d = 0.03)
)

run!(leaf,meteo)
leaf.status.Rn
julia> 12.902547446281233

# Using the model from Medlyn et al. (2011) for Gs:
leaf = ModelList(
    energy_balance = Monteith(),
    photosynthesis = Fvcb(),
    stomatal_conductance = Medlyn(0.03, 12.0),
    status = (Ra_SW_f = 13.747, sky_fraction = 1.0, aPPFD = 1500.0, d = 0.03)
)

run!(leaf,meteo)
leaf[:Rn]
leaf[:Ra_LW_f]
leaf[:A]

DataFrame(leaf)

References

Duursma, R. A., et B. E. Medlyn. 2012. « MAESPA: a model to study interactions between water limitation, environmental drivers and vegetation function at tree and stand levels, with an example application to [CO2] × drought interactions ». Geoscientific Model Development 5 (4): 919‑40. https://doi.org/10.5194/gmd-5-919-2012.

Monteith, John L., et Mike H. Unsworth. 2013. « Chapter 13 - Steady-State Heat Balance: (i) Water Surfaces, Soil, and Vegetation ». In Principles of Environmental Physics (Fourth Edition), edited by John L. Monteith et Mike H. Unsworth, 217‑47. Boston: Academic Press.

Schymanski, Stanislaus J., et Dani Or. 2017. « Leaf-Scale Experiments Reveal an Important Omission in the Penman–Monteith Equation ». Hydrology and Earth System Sciences 21 (2): 685‑706. https://doi.org/10.5194/hess-21-685-2017.

Vezy, Rémi, Mathias Christina, Olivier Roupsard, Yann Nouvellon, Remko Duursma, Belinda Medlyn, Maxime Soma, et al. 2018. « Measuring and modelling energy partitioning in canopies of varying complexity using MAESPA model ». Agricultural and Forest Meteorology 253‑254 (printemps): 203‑17. https://doi.org/10.1016/j.agrformet.2018.02.005.

source
PlantSimEngine.run!Function
run!(::Fvcb, models, status, meteo, constants=Constants())

Coupled photosynthesis and conductance model using the Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981) that models the assimilation as the most limiting factor between three processes:

  • RuBisCo-limited photosynthesis, when the kinetics of the RuBisCo enzyme for fixing

CO₂ is at its maximum (RuBisCo = Ribulose-1,5-bisphosphate carboxylase-oxygenase). It happens mostly when the CO₂ concentration in the stomata is too low. The main parameter is VcMaxRef, the maximum rate of RuBisCo activity at reference temperature. See get_Cᵢᵥ for the computation.

  • RuBP-limited photosynthesis, when the rate of RuBP (ribulose-1,5-bisphosphate) regeneration

associated with electron transport rates on the thylakoid membrane (RuBP) is limiting. It happens mostly when light is limiting, or when CO₂ concentration is rather high. It is parameterized using JMaxRef, the potential rate of electron transport. See get_Cᵢⱼ for the computation.

  • TPU-limited photosynthesis, when the rate at which inorganic phosphate is released for

regenerating ATP from ADP during the utilization of triose phosphate (TPU) is limiting. It happens at very high assimilation rate, when neither light or CO₂ are limiting factors. The parameter is TPURef.

The computation in this function is made following Farquhar & Wong (1984), Leuning et al. (1995), and the MAESPA model (Duursma et al., 2012).

The resolution is analytical as first presented in Baldocchi (1994), and needs Cₛ as input.

Triose phosphate utilization (TPU) limitation is taken into account as proposed in Lombardozzi (2018) (i.e. Aₚ = 3 * TPURef, making the assumption that glycolate recycling is set to 0). TPURef is set at 9999.0 by default, meaning there is no limitation of photosynthesis by TPU. Note that TPURef can be (badly) approximated using the simple equation TPURef = 0.167 * VcMaxRef as presented in Lombardozzi (2018).

If you prefer to use Gbc, you can use the iterative implementation of the Fvcb model FvcbIter

If you want a version that is de-coupled from the stomatal conductance use FvcbRaw, but you'll need Cᵢ as input of the model.

Returns

Modify the first argument in place for A, Gₛ and Cᵢ:

  • A: carbon assimilation (μmol[CO₂] m-2 s-1)
  • Gₛ: stomatal conductance for CO₂ (mol[CO₂] m-2 s-1)
  • Cᵢ: intercellular CO₂ concentration (ppm)

Arguments

  • ::Fvcb: the Farquhar–von Caemmerer–Berry (FvCB) model
  • models: a ModelList struct holding the parameters for the model with

initialisations for: - Tₗ (°C): leaf temperature - aPPFD (μmol m-2 s-1): absorbed Photosynthetic Photon Flux Density - Cₛ (ppm): Air CO₂ concentration at the leaf surface - Dₗ (kPa): vapour pressure difference between the surface and the saturated air vapour pressure in case you're using the stomatal conductance model of Medlyn.

  • status: A status, usually the leaf status (i.e. leaf.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Note

Tₗ, aPPFD, Cₛ (and Dₗ if you use Medlyn) must be initialized by providing them as keyword arguments (see examples). If in doubt, it is simpler to compute the energy balance of the leaf with the photosynthesis to get those variables. See AbstractEnergy_BalanceModel for more details.

Examples

using PlantBiophysics, PlantMeteo
meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

leaf =
    ModelList(
        photosynthesis = Fvcb(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Tₗ = 25.0, aPPFD = 1000.0, Cₛ = 400.0, Dₗ = meteo.VPD)
    )
# NB: we need to initalise Tₗ, aPPFD and Cₛ.
# NB2: we provide the name of the process before the model but it is not mandatory.

run!(leaf,meteo,PlantMeteo.Constants())
leaf.status.A
leaf.status.Cᵢ

Note that we use VPD as an approximation of Dₗ here because we don't have the leaf temperature (i.e. Dₗ = VPD when Tₗ = T).

References

Baldocchi, Dennis. 1994. « An analytical solution for coupled leaf photosynthesis and stomatal conductance models ». Tree Physiology 14 (7-8‑9): 1069‑79. https://doi.org/10.1093/treephys/14.7-8-9.1069.

Duursma, R. A., et B. E. Medlyn. 2012. « MAESPA: a model to study interactions between water limitation, environmental drivers and vegetation function at tree and stand levels, with an example application to [CO2] × drought interactions ». Geoscientific Model Development 5 (4): 919‑40. https://doi.org/10.5194/gmd-5-919-2012.

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Leuning, R., F. M. Kelliher, DGG de Pury, et E.D. Schulze. 1995. « Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

Lombardozzi, L. D. et al. 2018.« Triose phosphate limitation in photosynthesis models reduces leaf photosynthesis and global terrestrial carbon storage ». Environmental Research Letters 13.7: 1748-9326. https://doi.org/10.1088/1748-9326/aacf68.

source
PlantSimEngine.run!Function
run!(::ConstantAGs, models, status, meteo, constants=Constants())

Constant photosynthesis coupled with a stomatal conductance model.

Returns

Modify the leaf status in place for A, Gₛ and Cᵢ:

  • A: carbon assimilation, set to leaf.photosynthesis.A (μmol[CO₂] m-2 s-1)
  • Gₛ: stomatal conductance for CO₂ (mol[CO₂] m-2 s-1)
  • Cᵢ: intercellular CO₂ concentration (ppm)

Arguments

  • ::ConstantAGs: a constant assimilation model coupled to a stomatal conductance model
  • models: a ModelList struct holding the parameters for the model with

initialisations for: - Cₛ (mol m-2 s-1): surface CO₂ concentration. - any other value needed by the stomatal conductance model.

  • status: A status, usually the leaf status (i.e. leaf.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Examples

using PlantBiophysics, PlantMeteo
meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)
leaf = ModelList(
    photosynthesis = ConstantAGs(),
    stomatal_conductance = Medlyn(0.03, 12.0),
    status = (Cₛ = 400.0, Dₗ = 2.0)
)

run!(leaf,meteo,PlantMeteo.Constants())

status(leaf, :A)
status(leaf, :Cᵢ)
source
PlantSimEngine.run!Function
run!(::FvcbIter, models, status, meteo, constants=Constants())

Photosynthesis using the Farquhar–von Caemmerer–Berry (FvCB) model for C3 photosynthesis (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981). Computation is made following Farquhar & Wong (1984), Leuning et al. (1995), and the Archimed model.

Iterative implementation, i.e. the assimilation is computed iteratively over Cᵢ. For the analytical resolution, see Fvcb.

Returns

Modify the first argument in place for A, Gₛ and Cᵢ:

  • A: carbon assimilation (μmol[CO₂] m-2 s-1)
  • Gₛ: stomatal conductance for CO₂ (mol[CO₂] m-2 s-1)
  • Cᵢ: intercellular CO₂ concentration (ppm)

Arguments

  • ::FvcbIter: Farquhar–von Caemmerer–Berry (FvCB) model with iterative resolution.
  • models: a ModelList struct holding the parameters for the model with

initialisations for: - Tₗ (°C): leaf temperature - aPPFD (μmol m-2 s-1): absorbed Photosynthetic Photon Flux Density - Gbc (mol m-2 s-1): boundary conductance for CO₂ - Dₗ (kPa): is the difference between the vapour pressure at the leaf surface and the saturated air vapour pressure in case you're using the stomatal conductance model of Medlyn.

  • status: A status, usually the leaf status (i.e. leaf.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Note

Tₗ, aPPFD, Gbc (and Dₗ if you use Medlyn) must be initialized by providing them as keyword arguments (see examples). If in doubt, it is simpler to compute the energy balance of the leaf with the photosynthesis to get those variables. See AbstractEnergy_BalanceModel for more details.

Examples

using PlantBiophysics, PlantMeteo
meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)

leaf =
    ModelList(
        photosynthesis = FvcbIter(),
        stomatal_conductance = Medlyn(0.03, 12.0),
        status = (Tₗ = 25.0, aPPFD = 1000.0, Gbc = 0.67, Dₗ = meteo.VPD)
    )
# NB: we need  to initalise Tₗ, aPPFD and Gbc.

run!(leaf,meteo,PlantMeteo.Constants())
leaf.status.A
leaf.status.Cᵢ

Note that we use VPD as an approximation of Dₗ here because we don't have the leaf temperature (i.e. Dₗ = VPD when Tₗ = T).

References

Baldocchi, Dennis. 1994. « An analytical solution for coupled leaf photosynthesis and stomatal conductance models ». Tree Physiology 14 (7-8‑9): 1069‑79. https://doi.org/10.1093/treephys/14.7-8-9.1069.

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Leuning, R., F. M. Kelliher, DGG de Pury, et E.D. Schulze. 1995. Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

source
PlantSimEngine.run!Function
run!(object, meteo, constants = Constants())

Computes the light interception of an object using the Beer-Lambert law.

Arguments

  • ::Beer: a Beer model, from the model list (i.e. m.light_interception)
  • models: A ModelList struct holding the parameters for the model with

initialisations for LAI (m² m⁻²): the leaf area index.

  • status: the status of the model, usually the model list status (i.e. m.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Examples

using PlantSimEngine, PlantBiophysics, PlantMeteo
m = ModelList(light_interception=Beer(0.5), status=(LAI=2.0,))

meteo = Atmosphere(T=20.0, Wind=1.0, P=101.3, Rh=0.65, Ri_PAR_f=300.0)

run!(m, meteo)

m[:aPPFD]
source
PlantSimEngine.run!Function
run!(::ConstantA; models, status, meteo, constants=Constants())

Constant photosynthesis (forcing the value).

Returns

Modify the leaf status in place for A with a constant value:

  • A: carbon assimilation, set to leaf.photosynthesis.A (μmol[CO₂] m-2 s-1)

Arguments

  • ::ConstantA: a constant assimilation model
  • models: a ModelList struct holding the parameters for the model.
  • status: A status, usually the leaf status (i.e. leaf.status)
  • meteo: meteorology structure, see Atmosphere
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Examples

meteo = Atmosphere(T = 20.0, Wind = 1.0, P = 101.3, Rh = 0.65)
leaf = ModelList(photosynthesis = ConstantA(26.0))

run!(leaf,meteo,Constants())

leaf.status.A
source
PlantSimEngine.run!Function
run!(::FvcbRaw, models, status, meteo=nothing, constants=Constants())

Direct implementation of the photosynthesis model for C3 photosynthesis from Farquhar–von Caemmerer–Berry (Farquhar et al., 1980; von Caemmerer and Farquhar, 1981).

Returns

Modify the first argument in place for A, the carbon assimilation (μmol[CO₂] m-2 s-1).

Arguments

  • ::FvcbRaw: the Farquhar–von Caemmerer–Berry (FvCB) model (not coupled)
  • models: a ModelList struct holding the parameters for the model with

initialisations for: - Tₗ (°C): leaf temperature - aPPFD (μmol m-2 s-1): absorbed Photosynthetic Photon Flux Density - Cₛ (ppm): Air CO₂ concentration at the leaf surface - Dₗ (kPa): vapour pressure difference between the surface and the saturated air vapour pressure in case you're using the stomatal conductance model of Medlyn.

  • status: A status, usually the leaf status (i.e. leaf.status)
  • constants = PlantMeteo.Constants(): physical constants. See PlantMeteo.Constants for more details

Note

Tₗ, aPPFD, Cₛ (and Dₗ if you use Medlyn) must be initialized by providing them as keyword arguments (see examples). If in doubt, it is simpler to compute the energy balance of the leaf with the photosynthesis to get those variables. See AbstractEnergy_BalanceModel for more details.

Examples

using PlantSimEngine
leaf = ModelList(photosynthesis = FvcbRaw(), status = (Tₗ = 25.0, aPPFD = 1000.0, Cᵢ = 400.0))
# NB: we need Tₗ, aPPFD and Cᵢ as inputs (see [`inputs`](@ref))

run!(leaf)
leaf.status.A
leaf.status.Cᵢ

# using several time-steps:
leaf =
    ModelList(
        photosynthesis = FvcbRaw(),
        status = (Tₗ = [20., 25.0], aPPFD = 1000.0, Cᵢ = [380.,400.0])
    )
# NB: we need Tₗ, aPPFD and Cᵢ as inputs (see [`inputs`](@ref))

run!(leaf)
DataFrame(leaf) # fetch the leaf status as a DataFrame

References

Baldocchi, Dennis. 1994. « An analytical solution for coupled leaf photosynthesis and stomatal conductance models ». Tree Physiology 14 (7-8‑9): 1069‑79. https://doi.org/10.1093/treephys/14.7-8-9.1069.

Duursma, R. A., et B. E. Medlyn. 2012. « MAESPA: a model to study interactions between water limitation, environmental drivers and vegetation function at tree and stand levels, with an example application to [CO2] × drought interactions ». Geoscientific Model Development 5 (4): 919‑40. https://doi.org/10.5194/gmd-5-919-2012.

Farquhar, G. D., S. von von Caemmerer, et J. A. Berry. 1980. « A biochemical model of photosynthetic CO2 assimilation in leaves of C3 species ». Planta 149 (1): 78‑90.

Leuning, R., F. M. Kelliher, DGG de Pury, et E.D. Schulze. 1995. « Leaf nitrogen, photosynthesis, conductance and transpiration: scaling from leaves to canopies ». Plant, Cell & Environment 18 (10): 1183‑1200.

Lombardozzi, L. D. et al. 2018.« Triose phosphate limitation in photosynthesis models reduces leaf photosynthesis and global terrestrial carbon storage ». Environmental Research Letters 13.7: 1748-9326. https://doi.org/10.1088/1748-9326/aacf68.

source
PlantSimEngine.run!Method

Constant stomatal conductance for CO₂ (mol m-2 s-1).

Note

meteo or gs_mod are just declared here for compatibility with the call from photosynthesis (need a constant way of calling the functions).

source