API Reference

Terrarium.AbstractClosureRelationType
abstract type AbstractClosureRelation

Base type for prognostic variable closure relations for differential equations of the form:

\[\frac{\partial g(u)}{\partial t} = F(u)\]

where F represents the RHS tendency as a function of the state variable u, and g(u) is a closure or constitutive relation that maps u to the physical units matching the tendency. Common examples in soil hydrothermal modeling are temperature-enthalpy and saturation-pressure relations.

source
Terrarium.AuxiliaryVariableType
AuxiliaryVariable{VD<:VarDims} <: AbstractVariable

Represents an auxiliary (sometimes called "diagnostic") variable with the given name and dims on the spatial grid.

source
Terrarium.ColumnGridType
ColumnGrid{NF,RectGrid<:Grids.RectilinearGrid} <: AbstractLandGrid

Represents a set of laterally independent vertical columns with dimensions (x, y, z) where x is the column dimension, y=1 is constant, and z is the vertical axis.

source
Terrarium.ConstantSoilCarbonDenistyType
struct ConstantSoilCarbonDenisty{NF} <: Terrarium.AbstractSoilBiogeochemistry{NF}

Naive implementation of soil biogeochemistry that just assumes there to be a constant organic content in all soil layers.

Properties:

  • ρ_soc::Any: Soil organic carbon density [kg/m^3]

  • ρ_org::Any: Pure organic matter density [kg/m^3]

  • por_org::Any: Natural porosity of organic material

source
Terrarium.FieldInitializersType
struct FieldInitializers{VarInits<:NamedTuple} <: Terrarium.AbstractInitializer

Initializer type that allows for direct specification of Field initializer functions.

Example:

initializer = FieldInitializers(:temperature => )
source
Terrarium.FieldInitializersMethod
FieldInitializers(; vars...)

Creates a FieldInitializers for the given variables. The first value of the pair is the name of the state variable while the second value should be a Function or callable struct of the form f(x,y,z)::Real where z is the vertical coordinate, decreasing with depth.

Example:

initializer = FieldInitializers(
    temperature = (x,y,z) -> 0.01*abs(z) + exp(z)*sin(2π*z)
)
source
Terrarium.GlobalRingGridType
struct GlobalRingGrid{NF, RingGrid<:SpeedyWeather.RingGrids.AbstractGrid, RectGrid<:Oceananigans.Grids.AbstractGrid} <: Terrarium.AbstractLandGrid{NF}

Convenience wrapper around ColumnGrid that defines the columns as points on a global RingGrid.

source
Terrarium.HomogeneousSoilType
struct HomogeneousSoil{NF} <: Terrarium.AbstractStratigraphy{NF}

Represents a soil stratigraphy of well mixed material with homogeneous soil texture.

Properties:

  • texture::CryoGrid.Soils.SoilTexture{NF, NF, NF} where NF: Material composition of mineral soil componnet
source
Terrarium.InverseQuadraticType

The inverse quadratic (or "quadratic parallel") bulk thermal conductivity formula (Cosenza et al. 2003):

\[k = [\sum_{i=1}^N θᵢ\sqrt{kᵢ}]^2\]

Cosenza, P., Guérin, R., and Tabbagh, A.: Relationship between thermal conductivity and water content of soils using numerical modelling, European Journal of Soil Science, 54, 581–588, https://doi.org/10.1046/j.1365-2389.2003.00539.x, 2003.

source
Terrarium.LUEPhotosynthesisType
struct LUEPhotosynthesis{NF} <: Terrarium.AbstractPhotosynthesis

Photosynthesis implementation from PALADYN (Willeit 2016) for C3 PFTs following the general light use efficiency model described in Haxeltine and Prentice 1996.

Authors: Maha Badri and Matteo Willeit

Properties:

  • τ25::Any: Value of τ at 25°C

  • Kc25::Any: Value of Kc at 25°C

  • Ko25::Any: Value of Ko at 25°C

  • q10_τ::Any: q10 for temperature-sensitive parameter τ

  • q10_Kc::Any: q10 for temperature-sensitive parameter Kc

  • q10_Ko::Any: q10 for temperature-sensitive parameter Ko

  • α_leaf::Any: Leaf albedo in PAR range [-]

  • cq::Any: Conversion factor for solar radiation at 550 nm from J/m² to mol/m² [mol/J]

  • k_ext::Any: Extinction coefficient for radiation through vegetation [-]

  • α_a::Any: Fraction of PAR assimilated at ecosystem level, relative to leaf level [-]

  • t_CO2_high::Any: Parameter, PFT specific [°C]

  • t_CO2_low::Any: Parameter, PFT specific [°C]

  • t_photos_high::Any: Parameter, PFT specific [°C]

  • t_photos_low::Any: Parameter, PFT specific [°C]

  • α_C3::Any: Intrinsic quantum efficiency of CO2 uptake in C3 plants [mol/mol]

  • C_mass::Any: Atomic mass of carbon [gC/mol]

  • θ_r::Any: Shape parameter [-]

  • day_length::Any: Day length [h/day]

  • sec_day::Any: Seconds per day [s/day]

source
Terrarium.MedlynStomatalConductanceType
struct MedlynStomatalConductance{NF} <: Terrarium.AbstractStomatalConductance

Stomatal conductance implementation from PALADYN (Willeit 2016) following the optimal stomatal conductance model (Medlyn et al. 2011).

Authors: Maha Badri and Matteo Willeit

Properties:

  • g1::Any: Parameter in optimal stomatal conductance formulation, Lin et al. 2015 [-], PFT specific
source
Terrarium.NamespaceType
Namespace <: AbstractVariable

Represents a new variable namespace, typically from a subcomponent of the model. It is (currently) assumed that tha name of the namesapce corresponds to a property defined on the model type.

source
Terrarium.PALADYNAutotrophicRespirationType
struct PALADYNAutotrophicRespiration{NF} <: Terrarium.AbstractAutotrophicRespiration

Autotrophic respiration implementation from PALADYN (Willeit 2016).

Authors: Maha Badri and Matteo Willeit

Properties:

  • cn_sapwood::Any: Sapwood parameter

  • cn_root::Any: Root parameter

  • aws::Any: Ratio of total to respiring stem carbon, Cox 2001, PFT specific [-]

source
Terrarium.PALADYNCarbonDynamicsType
struct PALADYNCarbonDynamics{NF} <: Terrarium.AbstractVegetationCarbonDynamics

Vegetation carbon dynamics implementation following PALADYN (Willeit 2016) but considering only the sum of the vegetation carbon pools. The subsequent splitting into Cleaf, Cstem, C_root is not implemented for now.

Authors: Maha Badri

Properties:

  • SLA::Any: Specific leaf area (Kattge et al. 2011) [m²/kgC], PFT specific

  • awl::Any: Allometric coefficient, modified from Cox 2001 to account for bwl=1 [kgC/m²], PFT specific

  • LAI_min::Any: Minimum Leaf Area Index modified from Clark et al. 2011 [m²/m²], PFT specific

  • LAI_max::Any: Maximum Leaf Area Index modified from Clark et al. 2011 [m²/m²], PFT specific

  • γL::Any: Leaf turnover rate (Kattge et al. 2011) [1/year], PFT specific

  • γR::Any: Root turnover rate [1/year], PFT specific

  • γS::Any: Stem turnover rate modified from Clark et al. 2011 [1/year], PFT specific

source
Terrarium.PALADYNPhenologyType
struct PALADYNPhenology{NF} <: Terrarium.AbstractPhenology

Vegetation phenology implementation from PALADYN (Willeit 2016).

Authors: Maha Badri and Matteo Willeit

Properties:

source
Terrarium.PALADYNVegetationDynamicsType
struct PALADYNVegetationDynamics{NF} <: Terrarium.AbstractVegetationDynamics

Vegetation dynamics implementation following PALADYN (Willeit 2016) for a single PFT based on the Lotka–Volterra approach.

Authors: Maha Badri

Properties:

  • ν_seed::Any: Vegetation seed fraction [-]

  • γv_min::Any: Minimum vegetation disturbance rate [1/year]

source
Terrarium.PhysicalConstantsType
struct PhysicalConstants{NF}

A collection of general physical constants that do not (usually) need to be varied in parameter calibration.

source
Terrarium.PrescribedHydraulicsType
PrescribedHydraulics{NF} <: AbstractSoilHydraulicProperties

Represents a simple case where soil hydraulic properties are directly prescribed to constant values. This is mostly provided just for testing, although it may be useful in certain cases where direct measurements of hydraulic properites are available.

Properties:

  • cond_sat::Any: Hydraulic conductivity at saturation [m/s]

  • porosity::Any: Prescribed soil porosity [-]

  • field_capacity::Any: Prescribed field capacity [-]

  • wilting_point::Any: Prescribed wilting point [-]

source
Terrarium.SURFEXHydraulicsType
struct SURFEXHydraulics{NF} <: Terrarium.AbstractSoilHydraulicProperties{NF}

SURFEX parameterization of mineral soil porosity (Masson et al. 2013).

source
Terrarium.SimulationType
struct Simulation{Model<:Terrarium.AbstractModel, TimeStepperCache<:Terrarium.AbstractTimeStepperCache, StateVars<:Terrarium.AbstractStateVariables} <: Terrarium.AbstractSimulation

Represents the state of a "simulation" for a given model. A simulation is here defined as a particular choice of model in conjunction with values for all state variables (including the Clock) and any necessary state caches for the time-stepper.

source
Terrarium.SoilEnergyBalanceType
struct SoilEnergyBalance{NF, HeatOperator<:Terrarium.AbstractHeatOperator, ThermalProps<:(SoilThermalProperties{NF})} <: Terrarium.AbstractSoilEnergyBalance{NF}

Standard implementation of the soil energy balance accounting for freezing and thawing of pore water/ice.

Properties:

  • operator::Terrarium.AbstractHeatOperator: Heat transport operator

  • thermal_properties::SoilThermalProperties: Soil thermal properties

source
Terrarium.SoilHydrologyType
struct SoilHydrology{NF, SoilWaterFluxes<:Terrarium.AbstractSoilWaterFluxes, SoilHydraulicProperties<:Terrarium.AbstractSoilHydraulicProperties{NF}, FC<:FreezeCurves.FreezeCurve} <: Terrarium.AbstractSoilHydrology{NF}

Properties:

  • fluxes::Terrarium.AbstractSoilWaterFluxes: Soil water flux scheme

  • hydraulic_properties::Terrarium.AbstractSoilHydraulicProperties: Soil hydraulic properties parameterization

  • freezecurve::FreezeCurves.FreezeCurve: Soil freezing and water retention curve(s) from FreezeCurves.jl

source
Terrarium.SoilModelType
struct SoilModel{NF, GridType<:Terrarium.AbstractLandGrid{NF}, Stratigraphy<:Terrarium.AbstractStratigraphy, SoilEnergy<:Terrarium.AbstractSoilEnergyBalance, SoilHydrology<:Terrarium.AbstractSoilHydrology, Biogeochemistry<:Terrarium.AbstractSoilBiogeochemistry, Constants<:PhysicalConstants{NF}, BoundaryConditions<:Terrarium.AbstractBoundaryConditions, Initializer<:Terrarium.AbstractInitializer, TimeStepper<:Terrarium.AbstractTimeStepper} <: Terrarium.AbstractSoilModel{NF}

General implementation of a 1D column model of soil energy, water, and carbon transport.

Properties:

  • grid::Terrarium.AbstractLandGrid: Spatial grid type

  • strat::Terrarium.AbstractStratigraphy: Stratigraphy of the soil

  • energy::Terrarium.AbstractSoilEnergyBalance: Soil energy balance

  • hydrology::Terrarium.AbstractSoilHydrology: Soil hydrology/water balance

  • biogeochem::Terrarium.AbstractSoilBiogeochemistry: Soil biogeochemistry

  • constants::PhysicalConstants: Physical constants

  • boundary_conditions::Terrarium.AbstractBoundaryConditions: Boundary conditions

  • initializer::Terrarium.AbstractInitializer: State variable initializer

  • time_stepping::Terrarium.AbstractTimeStepper: Timestepping scheme

source
Terrarium.SoilThermalPropertiesType
struct SoilThermalProperties{NF, CondWeighting}

Properties:

  • cond::SoilThermalConductivities: Thermal conductivities for all constituents

  • cond_bulk::Any: Thermal conductivity mixing scheme

  • heatcap::SoilHeatCapacities: Thermal conductivities for all constituents

source
Terrarium.StateVariablesType
struct StateVariables{prognames, tendnames, auxnames, subnames, closurenames, ProgVars, TendVars, AuxVars, SubVars, Closures, ClockType} <: Terrarium.AbstractStateVariables

Container type for all Fields corresponding to state variables defined by a model. StateVariables partitions the fields into three categories: prognostic, tendencies, and auxiliary. Prognostic variables are those which characterize the state of the system and are assigned tendencies to be integrated by the timestepper. Auxiliary fields are additional state variables derived from the prognostic state variables but which are conditionally independent of their values at the previous time step given the current prognostic state. It is worth noting that tendencies are also treated internally as auxiliary variables; however, they are assigned their own category here since they need to be handled separately by the timestepping scheme.

source
Terrarium.TemperatureEnergyClosureType
TemperatureEnergyClosure

Defines the constitutive relationship between temperature and the internal energy, U, of the system, i.e:

\[U(T) = T\times C(T) - L_f \theta_{wi} (1 - F(T))\]

where T is temperature, C(T) is the temperature-dependent heat capacity, L_f is the volumetric latent heat of fusion, and F(T) is the constitutive relation between temperature and the unfrozen fraction of pore water. Note that this formulation implies that the zero

source
Terrarium.VegetationModelType
struct VegetationModel{NF, Photosynthesis<:Terrarium.AbstractPhotosynthesis, StomatalConducatance<:Terrarium.AbstractStomatalConductance, AutotrophicRespiration<:Terrarium.AbstractAutotrophicRespiration, CarbonDynamics<:Terrarium.AbstractVegetationCarbonDynamics, VegetationDynamics<:Terrarium.AbstractVegetationDynamics, Phenology<:Terrarium.AbstractPhenology, GridType<:Terrarium.AbstractLandGrid{NF}, Constants<:PhysicalConstants{NF}, BoundaryConditions<:Terrarium.AbstractBoundaryConditions, Initializer<:Terrarium.AbstractInitializer, TimeStepper<:Terrarium.AbstractTimeStepper} <: Terrarium.AbstractVegetationModel{NF}

Model for natural (unmanaged) vegetation processes for a single plant functional type (PFT). Multiple PFTs can be later handled with a TiledVegetationModel type that composes multiple VegetationModels with different parameters for each PFT.

Properties:

  • grid::Terrarium.AbstractLandGrid: Spatial grid type

  • photosynthesis::Terrarium.AbstractPhotosynthesis: Photosynthesis scheme

  • stomatal_conductance::Terrarium.AbstractStomatalConductance: Stomatal conducantance scheme

  • autotrophic_respiration::Terrarium.AbstractAutotrophicRespiration: Autotrophic respiration scheme

  • phenology::Terrarium.AbstractPhenology: Phenology scheme

  • carbon_dynamics::Terrarium.AbstractVegetationCarbonDynamics: Vegetation carbon pool dynamics

  • vegetation_dynamics::Terrarium.AbstractVegetationDynamics: Vegetation population density or coverage fraction dynamics

  • constants::PhysicalConstants: Physical constants

  • boundary_conditions::Terrarium.AbstractBoundaryConditions: Boundary conditions

  • initializer::Terrarium.AbstractInitializer: State variable initializer

  • time_stepping::Terrarium.AbstractTimeStepper: Timestepping scheme

source
Terrarium.XYType
XY <: VarDims

Indicator type for variables that should be assigned a 2D (lateral only) field on their associated grid.

source
Terrarium.XYZType
XYZ <: VarDims

Indicator type for variables that should be assigned a 3D field on their associated grid.

source
Terrarium.compute_APARMethod
compute_APAR(photo, swdown, LAI)

Computes absorbed PAR limited by the fraction of PAR assimilated at ecosystem level APAR [mol/m²/day], Eq. 62, PALADYN (Willeit 2016).

source
Terrarium.compute_AgMethod
compute_Ag(photo, c_1, c_2, APAR, Vc_max, β)

Computes the daily gross photosynthesis Ag [gC/m²/day], Eqn 2, Haxeltine & Prentice 1996

source
Terrarium.compute_AndMethod
compute_And(photo, c_1, c_2, APAR, Vc_max, β, Rd)

Computes the total daytime net photosynthesis And [gC/m²/day], Eqn 19, Haxeltine & Prentice 1996 + Eq. 65, PALADYN (Willeit 2016).

source
Terrarium.compute_C_veg_tendMethod
compute_C_veg_tend(vegcarbon_dynamics, LAI_b, NPP)

Computes the C_veg tendency based on NPP and the balanced Leaf Area Index LAI_b, Eq. 72, PALADYN (Willeit 2016)

source
Terrarium.compute_JE_JCMethod
compute_JE_JC(photo, c_1, c_2, APAR, Vc_max)

Computes the PAR-limited and the rubisco-activity-limited photosynthesis rates JE and JC [gC/m²/day], Eqn 3+5, Haxeltine & Prentice 1996.

source
Terrarium.compute_NPPMethod
compute_NPP(autoresp, GPP, Ra)

Computes Net Primary Productivity NPP as the difference between Gross Primary Production GPP and autotrophic respiration Ra in [kgC/m²/day].

source
Terrarium.compute_RaMethod
compute_Ra(
    autoresp,
    vegcarbon_dynamics,
    T_air,
    Rd,
    phen,
    C_veg,
    GPP
)

Computes autotrophic respiration Ra as the sum of maintenance respiration Rm and growth respiration Rg in [kgC/m²/day].

source
Terrarium.compute_RdMethod
compute_Rd(photo, Vc_max, β)

Computes the daily leaf respiration Rd [gC/m²/day], Eqn 10, Haxeltine & Prentice 1996 and Eq. 10 PALADYN (Willeit 2016).

source
Terrarium.compute_RmMethod
compute_Rm(
    autoresp,
    vegcarbon_dynamics,
    T_air,
    Rd,
    phen,
    C_veg
)

Computes maintenance respiration Rm in [kgC/m²/day].

source
Terrarium.compute_Vc_maxMethod
compute_Vc_max(
    photo,
    c_1,
    APAR,
    Kc,
    Ko,
    Γ_star,
    pres_i,
    pres_O2
)

Computes the maximum daily rate of net photosynthesis Vc_max [gC/m²/day], following the coordination hypothesis (acclimation), see Harrison 2021 Box 2. Note: this is not the same formula in PALADYN paper, this implementaion is taken from the code

source
Terrarium.compute_auxiliary!Function
compute_auxiliary!(state, model::AbstractModel)

Computes updates to all auxiliary variables based on the current prognostic state of the model.

source
Terrarium.compute_c1_c2Method
compute_c1_c2(photo, T_air, Γ_star, Kc, Ko, pres_i, pres_O2)

Computes factor for light-limited assimilation c_1 and factor for RuBisCO-limited assimilation c_2, Eqs. C4+C5, PALADYN (Willeit 2016).

source
Terrarium.compute_photosynthesisMethod
compute_photosynthesis(
    photo,
    T_air,
    swdown,
    pres,
    co2,
    LAI,
    λc
)

Computes Gross Primary Production GPPin [kgC/m²/day] and leaf respiration Rd in [gC/m²/day]

source
Terrarium.compute_pres_iMethod
compute_pres_i(photo, λc, pres_a)

Computes intercellular CO2 partial pressure [Pa], Eq. 67, PALADYN (Willeit 2016).

source
Terrarium.compute_tendencies!Function
compute_tendencies!(state, model::AbstractModel)

Computes tendencies for all prognostic state variables for model stored in the given state.

source
Terrarium.compute_vpdMethod
compute_vpd(T_air, q_air, pres)

Computes the vapor pressure deficit from air temperature, specific humidity, and surface pressure.

source
Terrarium.compute_Λ_locMethod
compute_Λ_loc(vegcarbon_dynamics, LAI_b)

Computes the local litterfall rate Λ_loc based on the balanced Leaf Area Index LAI_b (assuming evergreen PFTs), Eq. 75, PALADYN (Willeit 2016).

source
Terrarium.compute_βMethod
compute_β(photo)

Computes the soil-moisture limiting factor β, Eq. 66, PALADYN (Willeit 2016).

source
Terrarium.compute_λ_NPPMethod
compute_λ_NPP(vegcarbon_dynamics, LAI_b)

Computes λ_NPP,a factor determining the partitioning of NPP between increase of vegetation carbon of the existing vegetated area and spreading of the given PFT based on the balanced Leaf Area Index LAI_b, Eq. 74, PALADYN (Willeit 2016).

source
Terrarium.compute_λcMethod
compute_λc(stomcond, vpd)

Computes the ratio of leaf-internal and air CO2 concentration λc, derived from the optimal stomatal conductance model (Medlyn et al. 2011), Eq. 71, PALADYN (Willeit 2016).

source
Terrarium.compute_ν_starMethod
compute_ν_star(veg_dynamics, ν)

Computes ν_star which is the maximum between the current vegetation fraction ν and the seed fraction ν_seed [-], to ensure that a PFT is always seeded.

source
Terrarium.compute_ν_tendMethod
compute_ν_tend(
    veg_dynamics,
    vegcarbon_dynamics,
    LAI_b,
    C_veg,
    ν
)

Computes the vegetation fraction tendency for a single PFT, Eq. 73, PALADYN (Willeit 2016).

source
Terrarium.create_fieldMethod
create_field(var, init, bcs, grid, args; kwargs...)

Initializes a Field on grid for the variable, var, using the given initializer, init, and boundary conditions, bcs. Additional arguments are passed direclty to the Field constructor. The location of the Field is determined by the specified VarDims on var.

source
Terrarium.fastmapMethod
fastmap(f::F, iter::NamedTuple...) where {F}

Same as map for NamedTuples but with guaranteed type stability. fastmap is a @generated function which unrolls calls to f into a loop-free tuple construction expression. All named tuples must have the same keys but in no particular order. The returned NamedTuple

source
Terrarium.fastmapMethod
fastmap(f::F, iter::NTuple{N,Any}...) where {F,N}

Same as map for NTuples but with guaranteed type stability. fastmap is a @generated function which unrolls calls to f into a loop-free tuple construction expression.

source
Terrarium.get_dtFunction
get_dt(timestepper::AbstractTimeStepper)

Get the current timestep size for the time stepper.

source
Terrarium.get_field_boundary_conditionsMethod
get_field_boundary_conditions(bcs::AbstractBoundaryConditions, var::AbstractVariable)

Retrieve the Field boundary conditions for the corresponding variable. Defaults to returning nothing if no BC is defined.

source
Terrarium.get_gridMethod
get_grid(model::AbstractModel)::AbstractGrid

Return the spatial grid associated with this model.

source
Terrarium.get_spacingMethod
get_spacing(spacing)

Return a Vector of vertical layer thicknesses according to the given discretization.

source
Terrarium.heatcapacityMethod
heatcapacity(props, fracs)

Computes the bulk heat capacity of a finite volume from the given volumetric fractions.

source
Terrarium.initialize!Function
initialize!(state, model::AbstractModel)
initialize!(state, model::AbstractModel, initializer::AbstractInitializer)

Calls initialize! on the model and its corresponding initializer. This method only needs to be implemented if initialization routines are necessary in addition to direct field/variable initializers.

source
Terrarium.initialize!Method

Resets the simulation clock and calls initialize!(state, model) on the underlying model which should reset all state variables to their values as defiend by the model initializer.

source
Terrarium.initializeMethod

Creates and initializes a Simulation for the given model with the given clock state. This method allocates all necessary Fields for the state variables and calls initialize!(sim). Note that this method is not type stable and should not be called in an Enzyme autodiff call.

source
Terrarium.is_adaptiveFunction
is_adaptive(timestepper::AbstractTimeStepper)

Returns true if the given time stepper is adaptive, false otherwise.

source
Terrarium.merge_duplicatesMethod
merge_duplicates(values::Tuple)

Filters out duplicates from the given tuple. Note that this method is not type stable or allocation-free!

source
Terrarium.namespaceMethod
namespace(name)

Convenience constructor method for Namespace provided only for consistency.

source
Terrarium.organic_fractionMethod
organic_fraction(idx, state, bgc::ConstantSoilCarbonDenisty)

Calculate the soil organic carbon fraction at the given grid index. For ConstantSoilCarbonDenisty, this simply returns a constant parameter. Implementations of soil carbon dynamics would want to compute this based on the prognostic state of the carbon content/pool stored in the soil.

source
Terrarium.organic_fractionMethod
organic_fraction(bgc::ConstantSoilCarbonDenisty)

Calculate the organic solid fraction based on the prescribed SOC and natural porosity/density of the organic material.

source
Terrarium.run!Method
run!(sim; steps, period, dt)

Run the simulation by steps or a period with dt timestep size (in seconds or Dates.Period).

source
Terrarium.timestep!Function
timestep!(state, model::AbstractModel, timestepper::AbstractTimeStepper, [dt = nothing])

Advance the model state by one time step, or by dt units of time.

source
Terrarium.timestep!Method

Convenience dispatch for timestep! that forwards to timestep!(state, model, get_time_stepping(model), dt).

source
Terrarium.timestep!Method
timestep!(sim)

Advance the simulation forward by one timestep with optional timestep size dt.

source
Terrarium.vardimsMethod
vardims(::AbstractVariable)

Retrieves the grid dimensions on which this variable is defined.

source
Terrarium.variablesFunction
variables(model::AbstractModel)

Return a Tuple of AbstractVariables (i.e. PrognosticVariable, AuxiliaryVariable, etc.) defined by the model. For models that consist of one or more sub-models, variables may optionally be grouped into namespaces by returning a NamedTuple where the keys correspond to the group names.

source
Terrarium.varnameMethod
varname(::AbstractVariable)
varname(::AbstractClosureRelation)
varname(::Pair{Symbol})

Retrieves the name of the given variable or closure. For closure relations, varname should return the name of the variable returned by the closure relation.

source