Photosynthesis

Models and Parameters

ClimaLand.Canopy.PModelType
PModel{FT,
            OPFT <: PModelParameters{FT},
            OPCT <: PModelConstants{FT}
            } <: AbstractPhotosynthesisModel{FT}

An implementation of the optimality photosynthesis model "P-model v1.0" of Stocker et al. (2020).

Stocker, B. D., Wang, H., Smith, N. G., Harrison, S. P., Keenan, T. F., Sandoval, D., Davis, T.,
    and Prentice, I. C.: P-model v1.0: an optimality-based light use efficiency model for simulating
    ecosystem gross primary production, Geosci. Model Dev., 13, 1545–1581,
    https://doi.org/10.5194/gmd-13-1545-2020, 2020.

The P-model computes photosynthesis rates at the canopy level, and ci, Γstar, Ko, Kc are in
units of Pa.
source
ClimaLand.Canopy.PModelMethod
function PModel{FT}(;
    cstar = FT(0.41),
    β = FT(146),
    ϕc = FT(0.087),
    ϕ0 = FT(NaN),
    ϕa0 = FT(0.352),
    ϕa1 = FT(0.022),
    ϕa2 = FT(-0.00034),
    α = FT(0.933),
    sc = LP.get_default_parameter(FT, :low_water_pressure_sensitivity),
    pc = LP.get_default_parameter(FT, :moisture_stress_ref_water_pressure),
) where {FT <: AbstractFloat}

Constructs a P-model (an optimality model for photosynthesis) using default parameters.

The following default parameters are used:

  • cstar = 0.41 (unitless) - 4 * dA/dJmax, assumed to be a constant marginal cost (Wang 2017, Stocker 2020)
  • β = 146 (unitless) - Unit cost ratio of Vcmax to transpiration (Stocker 2020)
  • ϕc = 0.087 (unitless) - constant linear scaling factor for the intrinsic quantum yield (hat{c}_L in Stocker 2020)
  • ϕ0 = NaN (unitless) - constant intrinsic quantum yield. If set to NaN, intrinsic quantum yield is computed from the temperature-dependent form (Stocker 2020)
  • ϕa0 = 0.352 (unitless) - constant term in quadratic intrinsic quantum yield (Stocker 2020)
  • ϕa1 = 0.022 (K^-1) - first order term in quadratic intrinsic quantum yield (Stocker 2020)
  • ϕa2 = -0.00034 (K^-2) - second order term in quadratic intrinsic quantum yield (Stocker 2020)
  • α = 0.933 (unitless) - 1 - 1/T where T is the timescale of Vcmax, Jmax acclimation. Here T = 15 days. (Mengoli 2022)
  • sc = 5e-6 (Pa^{-1}) - sensitivity to low water pressure in the moisture stress factor [Tuzet et al. (2003)]
  • pc = -2e6 (Pa) - reference water pressure for the moisture stress factor [Tuzet et al. (2003)]
source
ClimaLand.Canopy.PModelParametersType
PModelParameters{FT<:AbstractFloat}

The required parameters for P-model (Stocker et al. 2020). Parameters are typically tunable with considerable uncertainty.

  • cstar: Constant describing cost of maintaining electron transport (unitless) Typical value = 0.41

  • β: Ratio of unit costs of transpiration and carboxylation (unitless) Typical value = 146

  • ϕc: Scaling parameter for temp-dependent intrinsic quantum yield (unitless) Typical value = 0.087

  • ϕ0: Temp-independent intrinsic quantum yield. If provided, overrides ϕc. (unitless) Typical value = 0.05

  • ϕa0: Constant term in temp-dependent intrinsic quantum yield (unitless).

  • ϕa1: First order term in temp-dependent intrinsic quantum yield (K^-1).

  • ϕa2: Second order term in temp-dependent intrinsic quantum yield (K^-2).

  • α: Timescale parameter used in EMA for acclimation of optimal photosynthetic capacities (unitless). Setting this to 0 represents no incorporation of past values. Since we update the EMA equation once per day, α = 1 - 1 day/τ where τ is the acclimation timescale in days.

  • sc: Sensitivity to low water pressure, in the moisture stress factor, (Pa^{-1}) [Tuzet et al. (2003)]

  • pc: Reference water pressure for the moisture stress factor (Pa) [Tuzet et al. (2003)]

source
ClimaLand.Canopy.PModelConstantsType
PModelConstants{FT<:AbstractFloat}

The required constants for P-model (Stocker et al. 2020). These are physical or biochemical constants that are not expected to change with time or space.

  • R: Gas constant (J mol^-1 K^-1)

  • Kc25: Michaelis-Menten parameter for carboxylation at 25°C (Pa)

  • Ko25: Michaelis-Menten parameter for oxygenation at 25°C (Pa)

  • To: Reference temperature equal to 25˚C (K)

  • ΔHkc: Energy of activation for Kc (J mol^-1)

  • ΔHko: Energy of activation for Ko (J mol^-1)

  • Drel: Relative diffusivity of CO2 in the stomatal pores, equal to 1.6.

  • ΔHΓstar: Effective energy of activation for Γstar (J mol^-1)

  • Γstar25: Γstar at 25 °C (Pa)

  • Ha_Vcmax: Effective energy of activation for Vcmax (J mol^-1)

  • Hd_Vcmax: Effective energy of deactivation for Vcmax (J mol^-1)

  • aS_Vcmax: Intercept term for dS in Vcmax deactivation factor (J K^-1 mol^-1)

  • bS_Vcmax: Slope term for dS in Vcmax deactivation factor (J K^-2 mol^-1)

  • Ha_Jmax: Effective energy of activation for Jmax (J mol^-1)

  • Hd_Jmax: Effective energy of deactivation for Jmax (J mol^-1)

  • aS_Jmax: Intercept term for dS in Jmax deactivation factor (J K^-1 mol^-1)

  • bS_Jmax: Slope term for dS in Jmax deactivation factor (J K^-2 mol^-1)

  • Mc: Molar mass of carbon (kg mol^-1)

  • oi: Intercellular O2 mixing ratio (unitless)

  • aRd: First order coefficient for temp-dependent Rd (K^-1)

  • bRd: Second order coefficient for temp-dependent Rd (K^-2)

  • fC3: Constant factor appearing the dark respiration term for C3 plants (unitless)

  • planck_h: Planck constant (J s)

  • lightspeed: Speed of light (m s^-1)

  • N_a: Avogadro constant (mol^-1)

  • ρ_water: Density of water (kg m^-3)

source
ClimaLand.Canopy.FarquharModelType

FarquharModel{FT, FP <: FarquharParameters{FT}}

A photosynthesis model taking leaf-level Vcmax25 and multiple other constants and predicting dark respiration at the leaf level, net photosynthesis at the leaf level, and GPP at the canopy level.

The Farquhar model computes photosynthetic rates using leaf-level properties, and so Vcmax25, Rd, An, etc are per unit area leaf. This version also works with ci, Γstar, Ko, Ki, in units of mol/mol, and not Pa.

source
ClimaLand.Canopy.FarquharModelMethod
FarquharModel{FT}(
    domain;
    photosynthesis_parameters = clm_photosynthesis_parameters(
        domain.space.surface,
    ),
    sc::FT = LP.get_default_parameter(FT, :low_water_pressure_sensitivity),
    pc::FT = LP.get_default_parameter(FT, :moisture_stress_ref_water_pressure),
) where {
    FT <: AbstractFloat,
    MECH <: Union{FT, ClimaCore.Fields.Field},
    VC <: Union{FT, ClimaCore.Fields.Field},
}

Creates a FarquharModel using default parameters of type FT.

The photosynthesis_parameters argument is a NamedTuple that contains

  • is_c3: a Float or Field indicating if plants are C3 (1) or C4 (0) (unitless)
  • Vcmax25: a Float or Field representing the maximum carboxylation rate at 25C (mol m^-2 s^-1)

By default, these parameters are set by the clm_photosynthesis_parameters function, which reads in CLM data onto the surface space as ClimaUtilities SpaceVaryingInputs.

The following additional default parameters are used:

  • sc = 5e-6 (Pa^{-1}) - sensitivity to low water pressure in the moisture stress factor [Tuzet et al. (2003)]
  • pc = -2e6 (Pa) - reference water pressure for the moisture stress factor [Tuzet et al. (2003)]
source
ClimaLand.Canopy.FarquharParametersType
FarquharParameters{
    FT<:AbstractFloat,
    MECH <: Union{FT, ClimaCore.Fields.Field},
    VC <: Union{FT, ClimaCore.Fields.Field},
}

The required parameters for the Farquhar photosynthesis model.

  • Vcmax25: Vcmax at 25 °C (mol CO2/m^2/s); leaf level

  • Γstar25: Γstar at 25 °C (mol/mol)

  • Kc25: Michaelis-Menten parameter for CO2 at 25 °C (mol/mol)

  • Ko25: Michaelis-Menten parameter for O2 at 25 °C (mol/mol)

  • ΔHkc: Energy of activation for CO2 (J/mol)

  • ΔHko: Energy of activation for oxygen (J/mol)

  • ΔHVcmax: Energy of activation for Vcmax (J/mol)

  • ΔHΓstar: Energy of activation for Γstar (J/mol)

  • ΔHJmax: Energy of activation for Jmax (J/mol)

  • ΔHRd: Energy of activation for Rd (J/mol)

  • To: Reference temperature equal to 25 degrees Celsius (K)

  • oi: Intercelluar O2 concentration (mol/mol); taken to be constant

  • ϕ: Quantum yield of photosystem II (Bernacchi, 2003; unitless)

  • θj: Curvature parameter, a fitting constant to compute J, unitless

  • fC3: Constant factor appearing the dark respiration term for C3, equal to 0.015.

  • fC4: Constant factor appearing the dark respiration term for C4, equal to 0.025.

  • sc: Sensitivity to low water pressure, in the moisture stress factor, (Pa^{-1}) [Tuzet et al. (2003)]

  • pc: Reference water pressure for the moisture stress factor (Pa) [Tuzet et al. (2003)]

  • Q10: Q10 temperature parameter for Vcmax and Rd for C4 photosynthesis; unitless

  • s1: Parameter appearing in temperature dependence of C4 Vcmax; K^{-1}

  • s2: Parameter appearing in temperature dependence of C4 Vcmax; K

  • s3: Parameter appearing in temperature dependence of C4 Vcmax; K^{-1}

  • s4: Parameter appearing in temperature dependence of C4 Vcmax; K

  • s5: Parameter appearing in temperature dependence of C4 Rd; K^{-1}

  • s6: Parameter appearing in temperature dependence of C4 Rd; K

  • E: Quantum yield for C4 photosynthesis; mol/mol

  • is_c3: Photosynthesis mechanism: 1.0 indicates C3, 0.0 indicates C4

source
ClimaLand.Canopy.FarquharParametersMethod
function FarquharParameters(
    ::Type{FT},
    is_c3::Union{FT, ClimaCore.Fields.Field};
    Vcmax25 = FT(5e-5),
    kwargs...  # For individual parameter overrides
)

function FarquharParameters(
    toml_dict::CP.AbstractTOMLDict,
    is_c3::Union{AbstractFloat, ClimaCore.Fields.Field};
    Vcmax25 = FT(5e-5),
    kwargs...  # For individual parameter overrides
)

Constructors for the FarquharParameters struct. Two variants:

  1. Pass in the float-type and retrieve parameter values from the default TOML dict.
  2. Pass in a TOML dictionary to retrieve parameter values.Possible calls:
ClimaLand.Canopy.FarquharParameters(Float64, 1.0)
# Kwarg overrides
ClimaLand.Canopy.FarquharParameters(Float64, 1.0; Vcmax25 = 99999999, pc = 444444444)
# TOML Dictionary:
import ClimaParams as CP
toml_dict = CP.create_toml_dict(Float32);
ClimaLand.Canopy.FarquharParameters(toml_dict, 1.0f0; Vcmax25 = 99999999, pc = 444444444)
source

FarquharModel Methods

ClimaLand.Canopy.arrhenius_functionFunction
arrhenius_function(T::FT, To::FT, R::FT, ΔH::FT)

Computes the Arrhenius function at temperature T given the reference temperature To=298.15K, the universal gas constant R, and the energy activation ΔH.

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.intercellular_co2_farquharFunction
intercellular_co2_farquhar(ca::FT, Γstar::FT, medlyn_factor::FT) where{FT}

Computes the intercellular CO2 concentration (mol/mol) given the atmospheric concentration (ca, mol/mol), the CO2 compensation (Γstar, mol/mol), and the Medlyn factor (unitless).

source
ClimaLand.Canopy.co2_compensation_farquharFunction
co2_compensation_farquhar(Γstar25::FT,
                 ΔHΓstar::FT,
                 T::FT,
                 To::FT,
                 R::FT) where {FT}

Computes the CO2 compensation point (Γstar), in units of mol/mol, as a function of its value at 25 °C (Γstar25), a constant energy of activation (ΔHΓstar), a standard temperature (To), the unversal gas constant (R), and the temperature (T).

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.rubisco_assimilationFunction
rubisco_assimilation(is_c3::AbstractFloat, args...)

Calls the correct rubisco assimilation function based on the is_c3; the returned quantity may be canopy-level (p-model) or leaf-level (standard Farquhar).

A is_c3 value of 1.0 corresponds to C3 photosynthesis and calls c3_rubisco_assimilation, while 0.0 corresponds to C4 photsynthesis and calls c4_rubisco_assimilation.

source
ClimaLand.Canopy.light_assimilationFunction
light_assimilation(is_c3::AbstractFloat, args...)

Calls the correct light assimilation function based on the is_c3; the returned quantity may be canopy-level (p-model) or leaf-level (standard Farquhar).

A is_c3 value of 1.0 corresponds to C3 photosynthesis and calls c3_light_assimilation, while 0.0 corresponds to C4 photsynthesis and calls c4_light_assimilation.

source
ClimaLand.Canopy.max_electron_transport_farquharFunction
max_electron_transport_farquhar(
    Vcmax25::FT,
    ΔHJmax::FT,
    T::FT,
    To::FT,
    R::FT,
) where {FT}

Computes the maximum potential rate of electron transport (Jmax), in units of mol/m^2/s, as a function of Vcmax at 25 °C (Vcmax25), a constant (ΔHJmax), a standard temperature (To), the unversal gas constant (R), and the temperature (T).

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.electron_transport_farquharFunction
electron_transport_farquhar(APAR::FT,
                   Jmax::FT,
                   θj::FT,
                   ϕ::FT) where {FT}

Computes the rate of electron transport (J), in units of mol/m^2/s, as a function of the maximum potential rate of electron transport (Jmax), absorbed photosynthetically active radiation (APAR), an empirical "curvature parameter" (θj; Bonan Eqn 11.21) and the quantum yield of photosystem II (ϕ).

In this function, all fluxes are per unit leaf area, including APAR.

See Ch 11, G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.net_photosynthesisFunction
net_photosynthesis(A::FT,
                   Rd::FT) where {FT}

Computes the total net carbon assimilation (An), in units of mol CO2/m^2/s, as a function of the gross assimilation (A) and dark respiration (Rd).

Note that if moisture stress is factored into the model, it must be applied correctly to Rd and A prior to this function.

Both A and Rd must either be leaf level (standard Farquhar), or canopy level (p-model), but not mixed. The returned quantity is similarily either leaf or canopy level.

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.moisture_stressFunction
moisture_stress(pl::FT,
                sc::FT,
                pc::FT) where {FT}

Computes the moisture stress factor (β), which is unitless, as a function of a constant (sc, 1/Pa), a reference pressure (pc, Pa), and the leaf water pressure (pl, Pa) .

See Eqn 12.57 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.dark_respiration_farquharFunction
dark_respiration_farquhar(is_c3::AbstractFloat, args...)

Calls the correct dark respiration function based on is_c3.

A is_c3 value of 1.0 corresponds to C3 photosynthesis and calls c3_dark_respiration, while 0.0 corresponds to C4 photsynthesis and calls c4_dark_respiration.

TODO: Generalize to p-model.

source
ClimaLand.Canopy.MM_KcFunction
MM_Kc(Kc25::FT,
      ΔHkc::FT,
      T::FT,
      To::FT,
      R::FT) where {FT}

Computes the Michaelis-Menten coefficient for CO2 (Kc), in units of mol/mol or Pa (depending on the units of Kc25), as a function of its value at 25 °C (Kc25), a constant (ΔHkc), a standard temperature (To), the unversal gas constant (R), and the temperature (T).

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.MM_KoFunction
MM_Ko(Ko25::FT,
      ΔHko::FT,
      T::FT,
      To::FT,
      R::FT) where {FT}

Computes the Michaelis-Menten coefficient for O2 (Ko), in units of mol/mol or Pa (depending on the units of Ko25), as a function of its value at 25 °C (Ko25), a constant (ΔHko), a standard temperature (To), the universal gas constant (R), and the temperature (T).

See Table 11.5 of G. Bonan's textbook, Climate Change and Terrestrial Ecosystem Modeling (2019).

source
ClimaLand.Canopy.compute_Vcmax_farquharFunction
compute_Vcmax_farquhar(is_c3::AbstractFloat, args...)

Calls the correct Vcmax function based on is_c3.

A is_c3 value of 1.0 corresponds to C3 photosynthesis and calls c3_compute_Vcmax, while 0.0 corresponds to C4 photsynthesis and calls c4_compute_Vcmax.

source

PModel Methods

ClimaLand.Canopy.compute_full_pmodel_outputsFunction
compute_full_pmodel_outputs(
    parameters::PModelParameters{FT},
    constants::PModelConstants{FT},
    T_canopy::FT,
    P_air::FT,
    VPD::FT,
    ca::FT,
    βm::FT,
    APAR::FT,
) where {FT}

Performs the P-model computations as defined in Stocker et al. (2020) and returns a dictionary of full outputs. See https://github.com/geco-bern/rpmodel for a code reference. This should replicate the behavior of the rpmodel package.

Args:

  • parameters: PModelParameters object containing the model parameters.
  • constants: PModelConstants object containing the model constants.
  • T_canopy: Canopy temperature (K).
  • P_air: Ambient air pressure (Pa).
  • VPD: Vapor pressure deficit (Pa).
  • ca: Ambient CO2 concentration (mol/mol).
  • βm: Soil moisture stress factor (unitless).
  • APAR: Absorbed photosynthetically active radiation (mol photons m^-2 s^-1).

Returns: named tuple with the following keys and descriptions: Output name Description (units) "gpp" Gross primary productivity (kg m^-2 s^-1) "gammastar" CO2 compensation point (Pa) "kmm" Effective MM coefficient for Rubisco-limited photosynthesis (Pa) "ns_star" Viscosity of water normalized to 25 deg C (unitless) "chi" Optimal ratio of intercellular to ambient CO2 (unitless) "xi" Sensitivity of χ to VPD (Pa^1/2) "mj" CO2 limitation factor for light-limited photosynthesis (unitless) "mc" CO2 limitation factor for Rubisco-limited photosynthesis (unitless) "ci" Intercellular CO2 concentration (Pa) "iwue" Intrinsic water use efficiency (Pa) "gs" Stomatal conductance (mol m^-2 s^-1 Pa^-1) "vcmax" Maximum rate of carboxlation (mol m^-2 s^-1) "vcmax25" Vcmax normalized to 25°C via modified-Arrhenius type function (mol m^-2 s^-1) "jmax" Maximum rate of electron transport (mol m^-2 s^-1) "jmax25" Jmax normalized to 25°C via modified-Arrhenius type function (mol m^-2 s^-1) "Rd" Dark respiration rate (mol m^-2 s^-1)

source
ClimaLand.Canopy.set_historical_cache!Function
function set_historical_cache!(p, Y0, model::PModel, canopy)

The P-model requires a cache of optimal Vcmax25, Jmax25, and ξ that represent past acclimated values. Before the simulation, we need to have some physically meaningful initial values for these variables, which live in p.canopy.photosynthesis.OptVars.

This method assumes that the acclimation is to the initial conditions of the simulation. Note that if the initial condition is e.g., nighttime, then initially the optimal Vcmax and Jmax are zero, so it will take ~1 month (two e-folding timescales for α corresponding to 2 week acclimation timescale) for the model to reach a physically meaningful state.

An alternative to this approach is to initialize the initial optimal values to some reasonable values based on a spun-up simulation.

source
set_historical_cache!(p, Y0, m::AbstractPhotosynthesisModel, canopy)

For some canopy components (namely the P-model), we need values at t-1 to compute new values, so this function sets the historical cache values for the photosynthesis model. However, for other photosynthesis models this is not needed, so do nothing by default.

source
ClimaLand.Canopy.update_optimal_EMAFunction
update_optimal_EMA(
    parameters::PModelParameters{FT},
    constants::PModelConstants{FT},
    OptVars::NamedTuple{(:ξ_opt, :Vcmax25_opt, :Jmax25_opt), Tuple{FT, FT, FT}},
    T_canopy::FT,
    P_air::FT,
    VPD::FT,
    ca::FT,
    βm::FT,
    APAR::FT,
    local_noon_mask::FT,
) where {FT}

This function updates the optimal photosynthetic capacities Vcmax25, Jmax25 and sensitivity of stomatal conductance to dryness (ξ) using an exponential moving average (EMA) that computes new optimal values at local noon, following Mengoli et al. (2022).

Args:

  • parameters: PModelParameters object containing the model parameters.
  • constants: PModelConstants object containing the model constants.
  • OptVars: NamedTuple containing the current optimal values of ξ, Vcmax25, and Jmax25.
  • T_canopy: Canopy temperature (K).
  • P_air: Ambient air pressure (Pa).
  • VPD: Vapor pressure deficit (Pa).
  • ca: Ambient CO2 concentration (mol/mol).
  • βm: Soil moisture stress factor (unitless).
  • APAR: Absorbed photosynthetically active radiation (mol photons m^-2 s^-1).
  • local_noon_mask: A mask (0 or 1) indicating whether the current time is within the local noon window.

Returns:

  • NamedTuple with updated optimal values of ξ, Vcmax25, and Jmax25

Reference: Mengoli, G., Agustí-Panareda, A., Boussetta, S., Harrison, S. P., Trotta, C., & Prentice, I. C. (2022). Ecosystem photosynthesis in land-surface models: A first-principles approach incorporating acclimation. Journal of Advances in Modeling Earth Systems, 14, e2021MS002767. https://doi.org/10.1029/2021MS002767

source
ClimaLand.Canopy.make_PModel_callbackFunction
function make_PModel_callback(
    ::Type{FT},
    start_date::Dates.DateTime,
    dt::Union{AbstractFloat, Dates.Period},
    canopy,
    longitude = nothing
) where {FT <: AbstractFloat}

This constructs a FrequencyBasedCallback for the P-model that updates the optimal photosynthetic capacities using an exponential moving average at local noon.

We check for local noon using the provided longitude (once passing in lat/lon for point/column domains, this can be automatically extracted from the domain axes) every dt. The time of local noon is expressed in seconds UTC and neglects the effects of obliquity and eccentricity, so it is constant throughout the year. This presents an error of up to ~20 minutes, but it is sufficient for our application here (since meteorological drivers are often updated at coarser time intervals anyway).

Args

  • FT: The floating-point type used in the model (e.g., Float32, Float64).
  • start_date: datetime object for the start of the simulation (UTC).
  • dt: timestep
  • canopy: the canopy object containing the P-model parameters and constants.
  • longitude: optional longitude in degrees for local noon calculation (default is nothing). If we are on a ClimaLand.Domains.Point, this will need to be supplied explicitly. Otherwise, if we are on a field, then the longitude at each point can be automatically extracted from the field axes (THIS STILL NEEDS TO BE TESTED)
source
ClimaLand.Canopy.intrinsic_quantum_yieldFunction
intrinsic_quantum_yield(
    T::FT,
    c::FT,
    ϕa0::FT,
    ϕa1::FT,
    ϕa2::FT
) where {FT}

Computes the intrinsic quantum yield of photosynthesis ϕ (mol/mol) as a function of temperature T (K) and a calibratable parameter c (unitless). The functional form given in Bernacchi et al (2003) and used in Stocker et al. (2020) is a second order polynomial in T (deg C) with coefficients ϕa0, ϕa1, and ϕa2.

source
ClimaLand.Canopy.compute_KmmFunction
compute_Kmm(
    T::FT,
    p::FT,
    Kc25::FT,
    Ko25::FT,
    ΔHkc::FT,
    ΔHko::FT,
    To::FT,
    R::FT,
    oi::FT
) where {FT}

Computes the effective Michaelis-Menten coefficient for Rubisco-limited photosynthesis (Kmm), in units Pa, as a function of temperature T (K), atmospheric pressure p (Pa), and constants: Kc25 (Michaelis-Menten coefficient for CO2 at 25 °C), Ko25 (Michaelis-Menten coefficient for O2 at 25 °C), ΔHkc (effective enthalpy of activation for Kc), ΔHko (effective enthalpy of activation for Ko), To (reference temperature, typically 298.15 K), R (universal gas constant), and oi (O2 mixing ratio, typically 0.209).

source
ClimaLand.Canopy.optimal_co2_ratio_c3Function
optimal_co2_ratio_c3(
    Kmm::FT,
    Γstar::FT,
    ηstar::FT,
    ca::FT,
    VPD::FT,
    β::FT,
    Drel::FT
) where {FT}

The p-model assumptions, that 1) plants optimize the relative costs of transpiration per unit carbon assimlated and costs of maintaining carboxylation capacity per unit carbon assimilated;

  1. coordination hypothesis (assimilation is limited simultaneously by both light and Rubisco)

are applied to compute the optimal ratio of intercellular to ambient CO2 concentration (χ) and auxiliary variables ξ, mj, and mc. mj and mc represent capacities for light and Rubisco- limited photosynthesis, respectively.

Parameters: Kmm (effective Michaelis-Menten coefficient for Rubisco-limited photosynthesis, Pa), Γstar (CO2 compensation point, Pa), ηstar (viscosity ratio), ca_pp (ambient CO2 partial pressure, Pa), VPD (vapor pressure deficit, Pa), β (moisture stress factor, unitless), Drel = 1.6 (relative diffusivity of water vapor with respect to CO2, unitless).

source
ClimaLand.Canopy.intercellular_co2_pmodelFunction
intercellular_co2_pmodel(ξ::FT, ca_pp::FT, Γstar::FT, VPD::FT) where {FT}

Computes the intercellular co2 concentration (ci) as a function of the optimal ξ (sensitivity to dryness), ca_pp (ambient CO2 partial pressure), Γstar (CO2 compensation point), and VPD (vapor pressure deficit).

source
ClimaLand.Canopy.gs_co2_pmodelFunction
gs_co2_pmodel(
    χ::FT,
    ca::FT,
    A::FT
) where {FT}

Computes the stomatal conductance of CO2 (gs_co2), in units of mol CO2/m^2/s via Fick's law. Parameters are the ratio of intercellular to ambient CO2 concentration (χ), the ambient CO2 concentration (ca, in mol/mol), and the assimilation rate (A, mol m^-2 s^-1). This is related to the conductance of water by a factor Drel (default value = 1.6).

source
ClimaLand.Canopy.gs_h2o_pmodelFunction
gs_h2o_pmodel(
    χ::FT,
    ca::FT,
    A::FT,
    Drel::FT
) where {FT}

Computes the stomatal conductance of H2O (gs_h2o), in units of mol H2O/m^2/s via Fick's law. Parameters are the ratio of intercellular to ambient CO2 concentration (χ), the ambient CO2 concentration (ca, in mol/mol), the assimilation rate (A, mol m^-2 s^-1), and the relative conductivity ratio Drel (unitless).

source
ClimaLand.Canopy.vcmax_pmodelFunction
vcmax_pmodel(
    ϕ0::FT,
    APAR::FT,
    mprime::FT,
    mc::FT
    βm::FT
) where {FT}

Computes the maximum rate of carboxylation assuming optimality and Aj = Ac using the intrinsic quantum yield (ϕ0), absorbed photosynthetically active radiation (APAR), Jmax-adjusted capacity (mprime), a Rubisco-limited capacity (mc), and empirical soil moisture stress factor (βm). See Eqns 16 and 6 in Stocker et al. (2020).

source
ClimaLand.Canopy.compute_LUEFunction
compute_LUE(
    ϕ0::FT,
    β::FT,
    mprime::FT,
    Mc::FT
) where {FT}

Computes light use efficiency (LUE) in kg C/mol from intrinsic quantum yield (ϕ0), moisture stress factor (β), and a Jmax modified capacity (mprime); see Eqn 17 and 19 in Stocker et al. (2020). Mc is the molar mass of carbon (kg/mol) = 0.0120107 kg/mol.

source
ClimaLand.Canopy.compute_mj_with_jmax_limitationFunction
compute_mj_with_jmax_limitation(
    mj::FT,
    cstar::FT
) where {FT}

Computes m' such that Aj = ϕ0 APAR * m' (a LUE model) by assuming that dA/dJmax = c is constant. cstar is defined as 4c, a free parameter. Wang etal (2017) derive cstar = 0.412 at STP and using Vcmax/Jmax = 1.88.

source
ClimaLand.Canopy.co2_compensation_pmodelFunction
co2_compensation_pmode(
    T::FT,
    To::FT,
    p::FT,
    R::FT,
    ΔHΓstar::FT
    Γstar25::FT
) where {FT}

Computes the CO2 compensation point (Γstar), in units Pa, as a function of temperature T (K) and pressure p (Pa). See Equation B5 of Stocker et al. (2020).

source
ClimaLand.Canopy.quadratic_soil_moisture_stressFunction
quadratic_soil_moisture_stress(
    θ::FT,
    meanalpha::FT = FT(1.0),
    a_hat::FT = FT(0.0),
    b_hat::FT = FT(0.685),
    θ0::FT = FT(0.0),
    θ1::FT = FT(0.6)
) where {FT}

Computes an empirical soil moisture stress factor (β) according to the quadratic functional form of Stocker et al. (2020), Eq 21, using the soil moisture (θ), AET/PET (meanalpha), and parameters a_hat, b_hat, θ0, and θ1. Default values are from the original paper.

source
ClimaLand.Canopy.compute_APAR_canopy_molesFunction
compute_APAR_canopy_moles(
    f_abs::FT,
    par_d::FT,
    λ_γ_PAR::FT,
    lightspeed::FT,
    planck_h::FT,
    N_a::FT
) where {FT}

Computes the absorbed photosynthetically active radiation for the canopy in mol photons m^-2 s^-1, given the fraction of absorbed PAR (f_abs), the PAR downwelling flux (par_d, in W m^-2), and the wavelength of PAR (λ_γ_PAR, in m), and the physical constants necessary to compute the energy per mol PAR photons.

source
ClimaLand.Canopy.compute_APAR_leaf_molesFunction
compute_APAR_leaf_moles(
    f_abs::FT,
    par_d::FT,
    λ_γ_PAR::FT,
    lightspeed::FT,
    planck_h::FT,
    N_a::FT,
    LAI::FT
) where {FT}

Computes the absorbed photosynthetically active radiation for the leaf in mol photons m^-2 s^-1, given the fraction of absorbed PAR (f_abs), the PAR downwelling flux (par_d, in W m^-2), and the wavelength of PAR (λ_γ_PAR, in m), the leaf area index LAI, and the physical constants necessary to compute the energy per mol PAR photons.

source