Skip to content
Jorn Bruggeman edited this page Feb 26, 2024 · 77 revisions

FABM comes with a "pyfabm" package for the Python programming language. This package enables you to access FABM's biogeochemical models directly from Python. For instance, you can enumerate a model's parameters and variables, or obtain temporal derivatives and diagnostics for any given environment and model state. In combination with a Python-based time integration scheme (e.g., scipy.integrate.solve_ivp), this also allows you to perform model simulations.

You can try pyfabm online without installing anything:

Binder

Installing

Getting pyfabm from Anaconda

Pre-built versions of pyfabm are available for Windows, Linux and Mac. These are based on the latest stable FABM release. Unless you want to modify the FABM source code or the underlying biogeochemical models, we recommend using such a pre-built version:

  1. Ensure you have Anaconda:

    • Linux/Mac: execute conda --version in a terminal
    • Windows: look for "Anaconda prompt" in the start menu

    If you do not have Anaconda, install Miniconda on LinuxWindows, or Mac.

  2. Open a terminal (Windows: "Anaconda prompt" in the start menu) and create an isolated environment with pyfabm:

    conda create -n fabm -c conda-forge pyfabm
    

    If you experience problems, we recommend you first execute conda update conda to ensure your conda is up to date. Should this fail because of lack of permissions, we recommend you install Miniconda as described under step 1. After you have an up-to-date conda, retry the conda create ... command.

  3. Activate the fabm environment in your terminal whenever you want to use pyfabm:

    conda activate fabm
    

Manual build and install

If you want to modify the FABM source code, for instance, to use new or modified biogeochemical models, you will need to build pyfabm yourself.

First, get the FABM source code.

To build pyfabm, you will need:

To install, open a command line terminal and execute:

python -m pip install <FABMDIR>

where <FABMDIR> is the path to the FABM source code.

If needed, you can customize the installation by creating a setup.cfg file in the FABM source directory. This file could contain, for instance,

[build_ext]
cmake_opts=-DFABM_EXTRA_INSTITUTES=foo -DFABM_FOO_BASE="<FOODIR>"
force=1
debug=1

In this example, cmake_opts contains any additional options that will be passed to cmake. force=1 specifies that the build directory needs to be be emptied before each build. debug=1 specifies that FABM should be built in debug mode instead of release mode.

After the pyfabm package has been built and installed, you can access it from python by typing import pyfabm. Please try this before continuing.

Examples

In addition to examples below, pyfabm comes with several Jupyter Notebooks. These can be found at <FABMDIR>/testcases/python. You can also try them online: Binder

Enumerating variables

import pyfabm
model = pyfabm.Model("fabm.yaml")
for variable in model.state_variables:
   print(f"  {variable.name} = {variable.long_name} ({variable.units})")

Run a simulation

import scipy.integrate
from matplotlib import pyplot
import pyfabm

# Create model (loads fabm.yaml)
model = pyfabm.Model()

# Configure the environment
# Note: the set of environmental dependencies depends on the loaded biogeochemical model.
model.cell_thickness = 10.0
model.dependencies["surface_downwelling_photosynthetic_radiative_flux"].value = 50.0
model.dependencies["downwelling_photosynthetic_radiative_flux"].value = 25.0

# Verify the model is ready to be used
assert model.start(), "One or more model dependencies have not been fulfilled."

# Time derivative
def dy(t, y):
    model.state[:] = y
    return model.getRates()

# Time-integrate over one year (note: FABM's internal time unit is seconds!)
result = scipy.integrate.solve_ivp(dy, [0.0, 365.0 * 86400], model.state)

# Plot results
fig, ax = pyplot.subplots()
ax.plot(result.t / 86400, result.y.T)
ax.legend([variable.path for variable in model.state_variables])
pyplot.show()

Included utilities

pyfabm comes with a number of utilities:

  • fabm_complete_yaml takes a FABM configuration (fabm.yaml) and adds descriptions of all parameters and variables (e.g., long names, units) in comments at the end of each line.
  • fabm_describe_model takes a FABM configuration (fabm.yaml) and writes out a list of all variables contained in the resulting model.
  • fabm_evaluate evaluates a FABM configuration (fabm.yaml) for a particular set of inputs (state variables, environmental dependencies). The values of these inputs can be provided in a YAML file with key: value pairs, in a NetCDF file or on the command line with -v/--values. This script can be used to diagnose model behaviour under circumstances that are known to cause problems, e.g., they crash 3D simulations.
  • fabm_stress_test takes a FABM configuration (fabm.yaml) and verifies that the resulting biogeochemical model returns valid derivatives under a wide variety of inputs (state variables, environmental dependencies). Different tests can be run, using either random values or extremes for inputs, and running randomized or exhaustive tests.

Each of the utilities can be run on the command line with -h as argument; this will list information about the script and its command line arguments.

Should you wan to use these scripts as example of how to use pyfabm: they can be found at <FABMDIR>/src/pyfabm/utils.