FABM Coupling¶
Framework for Aquatic Biogeochemical Models (FABM) coupling interface. Provides state containers, input-variable adapters, and the chunked time loop for linking GOTM with FABM-compatible biogeochemical models.
See Biogeochemistry and FABM Coupling for the full scientific description of the coupling architecture.
Configuration¶
Configuration helpers for the Python-level FABM bridge.
- class pygotm.fabm.config.FABMConfig(use, config_path, freshwater_impact, repair_state, shade_feedback, albedo_feedback, surface_drag_feedback)[source]¶
Bases:
objectResolved FABM settings from a GOTM YAML document.
- Parameters:
- pygotm.fabm.config.fabm_enabled(document)[source]¶
Return whether
fabm.useis enabled in a GOTM YAML document.
Engine¶
FABMEngine.output_variable_specs() exposes the model-discovered NetCDF
surface: interior state and output-enabled diagnostics are reported as
z-profiles, while surface, bottom, and horizontal variables are reported as
scalars. Each spec carries the normalized output name, units, and long name
from pyfabm.
Thin Python wrapper around pyfabm.Model.
- class pygotm.fabm.engine.FABMEngine(config_path, *, model_factory=None)[source]¶
Bases:
objectManage one pyfabm model while pyGOTM owns transport and timestepping.
- initialize(nlev=None, h_col=None, *, skip_start=False)[source]¶
Construct and start the pyfabm model for a 1-D water column.
nlev is the number of grid cells (pyGOTM
nlev). When provided, the model is created withshape=(nlev,)so pyfabm allocates 1-D arrays of the correct size. h_col (length nlev, bottom→surface) is the initial cell-thickness array required by pyfabm beforestart().
- set_dependency(name, value)[source]¶
Set one FABM dependency from a scalar or contiguous float64 array.
- set_dependency_if_present(name, value)[source]¶
Set one optional FABM dependency and report whether it exists.
- get_rates(*, surface=True, bottom=True, time=None)[source]¶
Return FABM source rates for each layer.
When surface=True (default for backwards compat) the returned rates include the air-sea surface exchange distributed over ALL layers by pyfabm — which is physically wrong for a 1D column. Pass
surface=False, bottom=Falseto get bulk-only reaction rates that can safely be applied to every layer, then handle boundary exchange explicitly on the top/bottom layers.
- get_vertical_movement()[source]¶
Return vertical movement velocities for all interior state variables.
Returns an array of shape
(n_state_vars, nlev)in m s⁻¹ (positive upward), orNoneif the underlying model does not expose this API. Must be called afterget_rates()has updated the model’s internal diagnostics.
- diagnostic(name, *, copy=True)[source]¶
Return one FABM diagnostic by name, or
Nonewhen unavailable.
- start()[source]¶
Call model.start() after dependencies are set; update state buffers.
Use this when
initialize()was called withskip_start=Trueto defer start until the caller has supplied initial dependency values.- Return type:
Coupling¶
Map pyGOTM-owned physical fields into FABM dependencies.
- class pygotm.fabm.coupling.FABMDependencyEngine(*args, **kwargs)[source]¶
Bases:
ProtocolSubset of
FABMEngineused by coupling helpers.
- pygotm.fabm.coupling.apply_fabm_dependencies(engine, *, temperature, practical_salinity, density, cell_thickness, downwelling_photosynthetic_radiative_flux)[source]¶
Set the core GOTM-to-FABM dependencies in bulk.
- pygotm.fabm.coupling.copy_bioshade_feedback(engine, target)[source]¶
Copy FABM light-attenuation diagnostics into
targetwhen present.- Return type:
- Parameters:
engine (FABMDependencyEngine)
target (ndarray)
FABM Loop¶
FABM time loop driven by stored GOTM hydrodynamic state.
- pygotm.fabm.fabm_loop.run_fabm_chunk(engine, chunk_params, output, hydro_T, hydro_S, hydro_rho, hydro_h, hydro_nuh, hydro_rad, cc_in, out_index_base, forcing_u10=None, forcing_v10=None, forcing_yearday=None, forcing_secondsofday=None, forcing_precip=None, hydro_taub=None, step_offset=0, is_first_chunk=True)[source]¶
Run pyfabm for one physics chunk and return updated (cc, out_index).
- Return type:
- Parameters:
engine (FABMEngine)
chunk_params (RuntimeParams)
output (RuntimeOutput)
hydro_T (np.ndarray)
hydro_S (np.ndarray)
hydro_rho (np.ndarray)
hydro_h (np.ndarray)
hydro_nuh (np.ndarray)
hydro_rad (np.ndarray)
cc_in (np.ndarray | None)
out_index_base (int)
forcing_u10 (np.ndarray | None)
forcing_v10 (np.ndarray | None)
forcing_yearday (np.ndarray | None)
forcing_secondsofday (np.ndarray | None)
forcing_precip (np.ndarray | None)
hydro_taub (np.ndarray | None)
step_offset (int)
is_first_chunk (bool)
Parameters¶
- cc_innp.ndarray or None
None on first chunk — engine.start() sets initial conditions. Shape
(n_vars, nlev)array carrying FABM state on subsequent chunks.- out_index_baseint
Output slot counter before this chunk starts.
- is_first_chunkbool
True on the first chunk — triggers engine.start() and records the IC diagnostic slot.
Returns¶
- ccnp.ndarray
Updated FABM state at end of chunk, shape
(n_vars, nlev).- out_indexint
Updated output slot counter after this chunk.
- pygotm.fabm.fabm_loop.run_fabm_loop(engine, params, output, hydro_T, hydro_S, hydro_rho, hydro_h, hydro_nuh, hydro_rad, forcing_u10=None, forcing_v10=None, forcing_yearday=None, forcing_secondsofday=None, forcing_precip=None, hydro_taub=None)[source]¶
Run pyfabm over every stored GOTM hydro step, then fill reference outputs.
The compiled Numba loop already completed. Hydrodynamic state buffers hold T, S, rho, h, nuh, rad at every step (shape
(nt+1, nlev+1)). This function walks those steps in Python, coupling pyfabm to stored GOTM state without entering any Numba kernel.FABM profile and scalar outputs are written into
output.fabm_z_profilesandoutput.fabm_scalarsat the same output-slot indices as the physics outputs.Optional forcing_u10 / forcing_v10 (shape
(nt+1,)) provide surface wind components for FABM models that needwind_speed. forcing_yearday and forcing_secondsofday provide the fractional GOTM calendar day for models needingnumber_of_days_since_start_of_the_year.- Return type:
- Parameters:
engine (FABMEngine)
params (RuntimeParams)
output (RuntimeOutput)
hydro_T (ndarray)
hydro_S (ndarray)
hydro_rho (ndarray)
hydro_h (ndarray)
hydro_nuh (ndarray)
hydro_rad (ndarray)
forcing_u10 (ndarray | None)
forcing_v10 (ndarray | None)
forcing_yearday (ndarray | None)
forcing_secondsofday (ndarray | None)
forcing_precip (ndarray | None)
hydro_taub (ndarray | None)
State Containers¶
GOTM–FABM coupling layer — translation of gotm_fabm.F90.
Provides the interface between the General Ocean Turbulence Model (GOTM) and
the Framework for Aquatic Biogeochemical Models (FABM). The Fortran source
delegates biogeochemistry to the external FABM library (use fabm,
use fabm_types, use fabm_config, use field_manager); this module
mirrors that boundary and keeps pyGOTM physics independent of any specific
pyfabm installation.
Fortran state arrays (cc, cc_diag, rhs) and observation records
are held in FabmState and FabmObservation respectively.
Public interface: configure_gotm_fabm(), gotm_fabm_create_model(),
init_gotm_fabm(), init_var_gotm_fabm(),
init_gotm_fabm_state(), start_gotm_fabm(),
set_env_gotm_fabm(), do_gotm_fabm(), clean_gotm_fabm(),
register_observation(), register_scalar_observation(),
register_bulk_observation(), register_horizontal_observation(),
register_field(), calculate_derived_input(),
calculate_conserved_quantities(), right_hand_side_rhs(),
right_hand_side_ppdd(), do_repair_state(), light(),
save_diagnostics(), calendar_date_interface(),
FabmState, FabmObservation.
Original FORTRAN authors: Jorn Bruggeman.
- class pygotm.fabm.gotm_fabm.FabmObservation(kind, variable_id, data, relax_tau=None)[source]¶
Bases:
objectRegistered observed FABM variable and optional relaxation data.
- class pygotm.fabm.gotm_fabm.FabmState(fabm_calc=False, model=None, freshwater_impact=True, repair_state=False, nlev=0, dt=0.0, cc=None, rhs=None, pp=None, dd=None, bioshade=None, observations=<factory>, registered_fields=<factory>, environment=<factory>, diagnostics=<factory>)[source]¶
Bases:
objectpyGOTM-owned FABM coupling state.
- Parameters:
-
observations:
list[FabmObservation]¶
- pygotm.fabm.gotm_fabm.calculate_conserved_quantities(state, nlev, h, total)[source]¶
Calculate vertically integrated conserved quantities.
- pygotm.fabm.gotm_fabm.calculate_derived_input(state, nlev)[source]¶
Ask FABM to calculate derived environmental inputs.
- pygotm.fabm.gotm_fabm.calendar_date_interface(julian)[source]¶
Convert a Julian day count to a calendar date.
- pygotm.fabm.gotm_fabm.center_depths_single(nlev, h, depth)[source]¶
Calculate local depth below surface from layer height.
Used internally to compute light field, and may be used by biogeochemical models as well.
- pygotm.fabm.gotm_fabm.configure_gotm_fabm(state, cfg=None)[source]¶
Configure the FABM coupling from a YAML-like mapping.
- pygotm.fabm.gotm_fabm.do_repair_state(state, nlev, location='')[source]¶
Repair invalid FABM state values.
- pygotm.fabm.gotm_fabm.gotm_driver_fatal_error(self, location, message)[source]¶
FABM driver fatal-error callback.
- pygotm.fabm.gotm_fabm.gotm_fabm_create_model(state, namlst=None, *, model_factory=None)[source]¶
Create the external FABM model object.
- pygotm.fabm.gotm_fabm.init_gotm_fabm(state, nlev, dt, field_manager=None)[source]¶
Initialise FABM coupling arrays.
- pygotm.fabm.gotm_fabm.init_gotm_fabm_state(state, nlev)[source]¶
Initialise FABM state variables from the external model when available.
- pygotm.fabm.gotm_fabm.light(state, nlev, attenuation=0.04)[source]¶
Compute a simple exponentially decaying light/shading profile.
- pygotm.fabm.gotm_fabm.par_from_background_single(nlev, h, rad, light_A, light_g2, depth, par_col)[source]¶
Calculate photosynthetically active radiation (PAR).
- pygotm.fabm.gotm_fabm.par_with_bioext_from_attenuation_single(nlev, attenuation, h, rad, light_A, light_g2, depth, par_col)[source]¶
Calculate PAR using background and biotic extinction.
- pygotm.fabm.gotm_fabm.register_bulk_observation(state, variable_id, data, relax_tau)[source]¶
Register an observed depth-varying FABM variable.
- pygotm.fabm.gotm_fabm.register_field(state, variable, prefix='', dimensions=(), data0d=None, data1d=None, part_of_state=False, used=True)[source]¶
Register a FABM field with a field-manager-like registry.
- pygotm.fabm.gotm_fabm.register_horizontal_observation(state, horizontal_id, data, relax_tau)[source]¶
Register an observed horizontal-slice FABM variable.
- pygotm.fabm.gotm_fabm.register_observation(state, variable_id, data, relax_tau=None)[source]¶
Dispatch observation registration by data shape.
- pygotm.fabm.gotm_fabm.register_scalar_observation(state, scalar_id, data)[source]¶
Register an observed scalar FABM variable.
- pygotm.fabm.gotm_fabm.right_hand_side_ppdd(state, first, numc, nlev, cc, pp, dd)[source]¶
Calculate split production and destruction terms.
- pygotm.fabm.gotm_fabm.right_hand_side_rhs(state, first, numc, nlev, cc, rhs)[source]¶
Calculate FABM right-hand-side source terms.
- pygotm.fabm.gotm_fabm.set_env_gotm_fabm(state, **environment)[source]¶
Store environmental arrays/scalars passed from GOTM.
- pygotm.fabm.gotm_fabm.start_gotm_fabm(state, nlev, field_manager=None)[source]¶
Start FABM after GOTM fields have been registered.
- pygotm.fabm.gotm_fabm.step_fabm_post_rates_single(nlev, dt, n_interior, n_surface, n_bottom, bulk_rates, surf_rates, bot_rates, cc)[source]¶
Add pelagic sink and source terms for all depth levels.
- pygotm.fabm.gotm_fabm.step_fabm_transport_single(nlev, dt, cnpar, precip, has_vert_move, n_interior, vert_move, h_step, nuh_step, cc, y, ws, adv_cu, au, bu, cu, du, ru, qu, l_sour, q_sour, tau_r, y_obs)[source]¶
Vertical advection and residual movement; Vertical diffusion.
- Return type:
- Parameters:
nlev (int)
dt (float)
cnpar (float)
precip (float)
has_vert_move (int)
n_interior (int)
vert_move (ndarray)
h_step (ndarray)
nuh_step (ndarray)
cc (ndarray)
y (ndarray)
ws (ndarray)
adv_cu (ndarray)
au (ndarray)
bu (ndarray)
cu (ndarray)
du (ndarray)
ru (ndarray)
qu (ndarray)
l_sour (ndarray)
q_sour (ndarray)
tau_r (ndarray)
y_obs (ndarray)
State-buffer utilities for pyGOTM-owned FABM tracer arrays.
- class pygotm.fabm.state.FABMStateBuffer(variable_names, state, rates, bioshade, diagnostics=<factory>)[source]¶
Bases:
objectpyGOTM-owned FABM state, rates, diagnostics, and light feedback arrays.
- Parameters:
Input Variables¶
FABM observed-variable input — translation of gotm_fabm_input.F90.
Initialises and reads one or more data files containing observed profiles for a subset of FABM state variables, then registers those observations with the GOTM–FABM coupling layer for relaxation.
Fortran types and list structures are mapped to Python dataclasses:
type_input_variable→InputVariablefirst_input_variable(linked list) →FabmInputState(variables list)
Fortran dependencies use input and use gotm_fabm are satisfied by
pygotm.input.input and pygotm.fabm.gotm_fabm respectively.
Public interface: configure_gotm_fabm_input(), init_gotm_fabm_input(),
fabm_input_create(), append_input(),
FabmInputState, InputVariable.
Original FORTRAN authors: Jorn Bruggeman.
- class pygotm.fabm.gotm_fabm_input.FabmInputState[source]¶
Bases:
objectLinked-list replacement for FABM input variables.
-
variables:
list[InputVariable]¶
-
variables:
- class pygotm.fabm.gotm_fabm_input.InputVariable(name, scalar_input=None, profile_input=None, interior_id=None, horizontal_id=None, scalar_id=None, ncid=-1, relax_tau=1000000000000000.0, relax_tau_bot=1000000000000000.0, relax_tau_surf=1000000000000000.0, h_bot=0.0, h_surf=0.0, relax_tau_1d=None)[source]¶
Bases:
objectInformation on an observed FABM variable.
- Parameters:
-
scalar_input:
ScalarInput|None= None¶
-
profile_input:
ProfileInput|None= None¶
- pygotm.fabm.gotm_fabm_input.append_input(state, input_variable)[source]¶
Append an input variable to the current FABM input list.
- Return type:
- Parameters:
state (FabmInputState)
input_variable (InputVariable)
- pygotm.fabm.gotm_fabm_input.configure_gotm_fabm_input(state, cfg=None)[source]¶
Configure FABM input descriptors from a mapping.
- Return type:
- Parameters:
state (FabmInputState)