Domains, Variables, Fields, and Data arrays.

A Model contains Domains, each of which contain Variables defining Fields which contain Data arrays, and Reactions with ReactionMethods that operate on the Fields to calculate model time evolution.

Model

Domains

PALEOboxes.DomainType
Domain

A model region containing Variables and Reactions that act on them.

Domain spatial size is defined by grid, which may be nothing to define a scalar Domain, or an AbstractMesh to define a spatially-resolved Domain with multiple cells.

Named data_dims may be set by set_data_dimension! to allow Variables with additional non-spatial dimensions, eg to represent quantities on a wavelength grid.

source

Grids

PALEOboxes.internal_sizeFunction
internal_size(::Type{<:AbstractSpace}, mesh::AbstractMesh; [subdomain=""] [space=:cell]) -> NTuple{ndims, Int}

Array size to use for model Variables.

All AbstractMesh concrete subtypes (UnstructuredVectorGrid, CartesianLinearGrid, ...) should implement this method.

Optional Keyword Arguments

  • subdomain::String="": a named subdomain
source
PALEOboxes.cartesian_sizeFunction
cartesian_size(mesh::AbstractMesh) -> NTuple{ndims, Int}

Optional (only regular Cartesian grids should implement this method): Array size of Cartesian Domain.

NB: this may be different from internal_size if the mesh implements a mapping eg to a Vector for internal model Variables.

source
PALEOboxes.Grids.set_subdomain!Function
set_subdomain!(grid::PB.AbstractMesh, subdomainname::AbstractString, subdom::PB.AbstractSubdomain, allowcreate::Bool=false)

Set Subdomain

source
PALEOboxes.Grids.UnstructuredColumnGridType
UnstructuredColumnGrid <: PB.AbstractMesh

Minimal Grid for a Vector Domain composed of columns (not necessarily forming a 2-D array).

Fields

  • ncells::Int total number of cells in this Domain
  • Icolumns::Vector{Vector{Int}}: Icolumns[n] should be the indices of column n, in order from surface to floor, where n is also the index of

any associated boundary surface.

  • z_coords::Vector{FixedCoord}: z coordinates of cell mid point, lower surface, upper surface
  • columnnames::Vector{Symbol}: optional column names
source
PALEOboxes.Grids.CartesianLinearGridType
CartesianLinearGrid <: PB.AbstractMesh

nD grid with netcdf CF1.0 coordinates, using Vectors for PALEO internal representation of Variables, with a mapping linear indices <–> some subset of grid indices.

The linear indices mapping should be set with set_linear_index. Conversion of Field values between the PALEO internal representation (a Vector with a linear index) and a Cartesian Array with multiple dimensions (for import and export of model output) is then implemented by cartesian_to_internal and internal_to_cartesian.

Fields

  • ncells::Int64: number of cells in Domain = length(linear_index) (may be subset of total points in prod(dims))
  • dimnames::Vector{String}: names of dimensions (ordered list)
  • dims::Vector{Int}: sizes of dimensions (ordered list)
  • coords::Vector{Vector{Float64}}: attached cell-centre coordinates for each dimension (ordered list)
  • coords_edges::Vector{Vector{Float64}}: attached cell-edge coordinates for each dimension (ordered list)
source
PALEOboxes.Grids.CartesianArrayGridType
CartesianArrayGrid <: PB.AbstractMesh

nD grid with netcdf CF1.0 coordinates, using n-dimensional Arrays for PALEO internal representation of Variables

Fields

  • ncells::Int64: number of cells in Domain = length(linear_index) (may be subset of total points in prod(dims))
  • dimnames::Vector{String}: names of dimensions (ordered list)
  • dims::Vector{Int}: sizes of dimensions (ordered list)
  • coords::Vector{Vector{Float64}}: attached cell-centre coordinates for each dimension (ordered list)
  • coords_edges::Vector{Vector{Float64}}: attached cell-edge coordinates for each dimension (ordered list)
source

Subdomains

PALEOboxes.Grids.BoundarySubdomainType
BoundarySubdomain <: PB.AbstractSubdomain

A 2D subdomain corresponding to the 2D boundary Domain associated with a 3D interior Domain:

  • indices[ibnd] is the index of the 3D interior cell corresponding to a 2D boundary cell ibnd.
source
PALEOboxes.Grids.InteriorSubdomainType
InteriorSubdomain <: PB.AbstractSubdomain

A 3D subdomain corresponding to the 3D interior Domain associated with a 2D boundary Domain: indices[iint] is either:

  • missing if iint is the index of an interior cell in the 3D Domain, or
  • the index of the 2D boundary Domain cell corresponding to the boundary-adjacent 3D interior Domain cell iint
source
PALEOboxes.Grids.subdomain_viewFunction
subdomain_view(values, subdomain::BoundarySubdomain) -> view

Create a view on values in an interior Domain to access cells corresponding to indices in a boundary Domain.

source
subdomain_view(values, subdomain::InteriorSubdomain) -> var

Return unmodified values from a boundary Domain, indices to access from interior supplied by subdomain_indices

source
subdomain_view(values, subdomain::Nothing) -> values

Fallback when subdomain == nothing

source
PALEOboxes.Grids.subdomain_indicesFunction
subdomain_indices(subdomain::BoundarySubdomain) -> nothing

No additional indices required to access Variables in an interior Domain from a boundary Domain (view created by subdomain_view is sufficient).

source
subdomain_indices(subdomain::InteriorSubdomain) -> indices

Return indices to access Variables in a boundary Domain from interior Domain (will have missing entries where interior cells do not correspond to boundary)

source
subdomain_indices(subdomain::Nothing) -> nothing

fallback when subdomain == nothing

source

Regions, Dimensions and Coordinates

PALEOboxes.Grids.get_regionFunction
get_region(grid::Union{PB.AbstractMesh, Nothing}, values; selectargs...) -> values_subset, (dim_subset::NamedDimension, ...)

Return the subset of values given by selectargs (Grid-specific keywords eg cell=, column=, ...) and corresponding dimensions (with attached coordinates).

source
get_region(grid::Nothing, values) -> values[]

Fallback for Domain with no grid, assumed 1 cell

source
get_region(grid::UnstructuredVectorGrid, values; cell) -> 
    values_subset, (dim_subset::NamedDimension, ...)

Keywords for region selection:

  • cell::Union{Int, Symbol}: an Int, or a Symbol to look up in cellnames
source
get_region(grid::UnstructuredColumnGrid, values; column, [cell=nothing]) -> 
    values_subset, (dim_subset::NamedDimension, ...)

Keywords for region selection:

  • column::Union{Int, Symbol}: (may be an Int, or a Symbol to look up in columnames)
  • cell::Int: optional cell index within column, highest cell is cell 1
source
get_region(grid::Union{CartesianLinearGrid{2}, CartesianArrayGrid{2}} , internalvalues; [i=i_idx], [j=j_idx]) ->
    arrayvalues_subset, (dim_subset::NamedDimension, ...)

Keywords for region selection:

  • i::Int: optional, slice along first dimension
  • j::Int: optional, slice along second dimension

internalvalues are transformed if needed from internal Field representation as a Vector length ncells, to an Array (2D if neither i, j arguments present, 1D if i or j present, 0D ie one cell if both present)

source
get_region(grid::Union{CartesianLinearGrid{3}, CartesianArrayGrid{3}}, internalvalues; [i=i_idx], [j=j_idx]) ->
    arrayvalues_subset, (dim_subset::NamedDimension, ...)

Keywords for region selection:

  • i::Int: optional, slice along first dimension
  • j::Int: optional, slice along second dimension
  • k::Int: optional, slice along third dimension

internalvalues are transformed if needed from internal Field representation as a Vector length ncells, to an Array (3D if neither i, j, k arguments present, 2D if one of i, j or k present, 1D if two present, 0D ie one cell if i, j, k all specified).

source

Grids may define name dimensions and attach coordinates for the convenience of output visualisation. Any coordinate information required by Reactions should be supplied as Variables.

PALEOboxes.NamedDimensionType
NamedDimension

A named dimension, with optional attached fixed coordinates coords

PALEO convention is that where possible coords contains three elements, for cell midpoints, lower edges, upper edges, in that order.

source

Variables

PALEO Variables exist in Domains and represent biogeochemical or other quantities. They are defined by Reactions as VariableReactions which are then linked to create VariableDomains, and contain Fields which represent data (a subtype of AbstractData) in a function space (a subtype of AbstractSpace) with dimensions defined by the Domain grid (a subtype of AbstractMesh)

Dataflow dependency between Variables is represented by two pairings:

PALEOboxes.VariableTypeType
@enum VariableType

Enumeration of VariableBase subtypes. Allowed values:

  • VariableReaction: VT_ReactProperty, VT_ReactDependency, VT_ReactContributor, VT_ReactTarget
  • VariableDomain : VT_DomPropDep, VT_DomContribTarget
source
PALEOboxes.VariableDomContribTargetType
VariableDomContribTarget <: VariableDomain

Model (Domain) VariableDomain linking a single VariableReaction{VT_ReactTarget} to multiple VariableReaction{VT_ReactContributor} and VariableReaction{VT_ReactDependency}.

source

Variable Attributes

PALEOboxes.AttributeType
Attribute{T}

Definition for Variable attribute name. Defines a data type T, a default_value, required (true if always present ie added with default_value when Variable is created), units, and an optional description.

Note that Variable attributes are stored as a per-Variable Dict(name => value), these definitions are only used to provide defaults, check types, and provide descriptive metadata.

ParseFromString should usually be Nothing: a value of Type T is then required when calling set_attribute!. If ParseFromString is true, then set_attribute! will accept an AbstractString and call Base.parse(T, strvalue) to convert to T. This allows eg an enum-valued Attribute to be defined by Attribute{EnumType, true} and implementing parse(EnumType, rawvalue::AbstractString)

source
PALEOboxes.StandardAttributesConstant
StandardAttributes

List of standard Variable attributes.

Some of these follow netCDF COARDS/CF conventions:

COARDS:https://ferret.pmel.noaa.gov/Ferret/documentation/coards-netcdf-conventions

CF conventions: https://cfconventions.org/cf-conventions/cf-conventions.html#_description_of_the_data

source
32×5 DataFrame
 name                       default_value  required  units      description
 Symbol                     Any            Bool      String     String
─────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
 data_dims                  ()                 true             Tuple of variable data dimension names, or empty for a scalar
 description                                   true             human-readable description
 field_data                 UndefinedData      true             AbstractData type Variable contains
 initial_delta              0.0                true  per mil    initial value for isotope variable to be applied eg to state or constant variable
 initial_value              0.0                true             initial value to be applied eg to state or constant variable
 is_constant                false              true             true if variable is not changed after initialisation
 norm_value                 1.0                true             normalisation value to be passed through to solver and optionally provided as model variable
 space                      CellSpace          true             function space Variable is defined on
 units                                         true             where possible should follow netcdf conventions
 vfunction                  VF_Undefined       true             host function (to label state variables etc)
 advect                     false             false             true to apply advective transport to this tracer
 advect_zmin                0.0               false  m          minimum height for transport
 charge                     missing           false             species charge
 check_length               true              false             true to check length matches length of linked VariableDomain
 deposition_velocity        missing           false  cm s-1     surface deposition velocity for atmospheric tracer
 diffusivity                missing           false  cm^2/s     species diffusivity
 diffusivity_speciesname                      false             species name to define diffusivity and charge
 gamma                      missing           false             bioirrigation scaling factor for sediment solute
 initialize_to_zero         false             false             request initialize to zero at start of each timestep.
 long_name                                    false             netcdf long descriptive name
 mesh                       default           false             grid mesh on which Variable is defined (empty for Domain spatial scalar)
 operatorID                 Int64[]           false             Reaction operatorIDs that modify this Variable
 rainout                    missing           false             normalized rainout rate for atmospheric tracer
 rate_processname                             false             process name for this reaction rate variable
 rate_species               String[]          false             Vector of reactant and product species for this reaction rate variable
 rate_stoichiometry         Float64[]         false             Vector of reactant and product stoichiometries for this reaction rate variable
 safe_name                                    false             optional short or escaped name for compatibility with other software
 specific_light_extinction  0.0               false  m^2 mol-1  wavelength-independent specific extinction in water column
 standard_name                                false             netcdf CF conventions standard name
 totalnames                 String[]          false             Vector of total Variable names for this species
 vertical_movement          0.0               false  m d-1      vertical advective transport (+ve upwards) in column
 vphase                     VP_Undefined      false             phase for concentrations in multiphase cells
PALEOboxes.VariableFunctionType
@enum VariableFunction

Allowed values of :vfunction Variable Attribute, defining the Variable function to the host ODE or DAE solver.

Explicit ODE problems with dS/dt = F(S) consist of pairs of S::VFStateExplicit, F::VFDeriv Variables.

An implicit ODE problem with dU/dt = F(U) where Total variables U are functions U(S) of State variables S will consist of pairs of U::VFTotal and F::VFDeriv Variables, and also the same number of S::VF_StateTotal (in no particular order).

Algebraic constraints C(S) = 0 include variables C::VFConstraint and the same number of S::VFState, with no corresponding VF_Deriv.

Not all solvers support all combinations.

source

TODO standard Variable names.

Use NetCDF CF standard names https://cfconventions.org/Data/cf-standard-names/76/build/cf-standard-name-table.html ?

Fields

Field, data and function space are defined by Variable Attributes in combination with an AbstractMesh

Examples:

PALEOboxes.FieldType
Field{D <: AbstractData, S <: AbstractSpace, V, N, M}

A Field of values::V of data type D defined on function space S over mesh::M and (optionally) with N data_dims::NTuple{N, NamedDimensions}.

source

Spaces

PALEOboxes.CellSpaceType
CellSpace <: AbstractSpace

A per-cell quantity. Use as Variable attribute :space to create a Variable with data array dimensions from Grid cells.

source
PALEOboxes.ColumnSpaceType
ColumnSpace <: AbstractSpace

A per-column quantity. Use as Variable attribute :space to create a Variable with data array dimensions from Grid columns.

source

Data types

PALEOboxes.AbstractDataType
AbstractData

Defines a Data type that can be composed with an AbstractSpace to form a Field

Concrete subtypes should implement:

allocate_values, check_values, zero_values!, dof_values, get_values_output

If the subtype needs to provide values for a numerical solver (eg as a state variable), it also needs to implement:

init_values!, copyfieldto!, copytofield!, add_field!, add_field_vec!

If the subtype has a representation as components, it should implement: num_components, get_components

source
PALEOboxes.UndefinedDataType
UndefinedData <: AbstractData

Undefined data type (no methods implemented). Used to indicate that a Variable can link to any data type.

source
PALEOboxes.ArrayScalarDataType
ArrayScalarData <: AbstractData

An Array of Field scalars (eg a intensity on a wavelength grid).

NB: only implements minimal methods sufficient to allow use as a model internal variable, not as a state variable.

source
PALEOboxes.allocate_valuesFunction
allocate_values(
    field_data::Type{<:AbstractData}, data_dims::Tuple{Vararg{NamedDimension}}, data_type, space::Type{<:AbstractSpace}, spatial_size::Tuple;
    thread_safe::Bool, allocatenans::Bool,
) -> values

allocate Field.values (eg an Array) for field_data with dimensions defined by spatial_size and data_dims

source
PALEOboxes.check_valuesFunction
check_values(
    existing_values, 
    field_data::Type{<:AbstractData},
    data_dims::Tuple{Vararg{NamedDimension}}, 
    data_type, 
    space::Type{<:AbstractSpace}, 
    spatial_size::Tuple{Integer, Vararg{Integer}}
)

Check existing_values is of suitable type, size etc for use as Field.values, throw exception if not.

source
PALEOboxes.zero_values!Function
zero_values!(values, field_data::Type{<:AbstractData}, data_dims::Tuple{Vararg{NamedDimension}}, space::Type{<:AbstractSpace}, cellrange)

Set values over spatial region cellrange to zero at start of main loop

source
PALEOboxes.dof_valuesFunction
dof_values(
    field_data::Type{<:AbstractData}, 
    data_dims::Tuple{Vararg{NamedDimension}}, 
    space::Type{<:AbstractSpace}, 
    mesh, 
    cellrange
) -> dof::Int

Return degrees-of-freedom for field_data over spatial region cellrange.

source
PALEOboxes.get_values_outputFunction

Optional: sanitize values for storing as model output. Default implementation is usually OK - only implement eg for Atomic types that should be converted to standard types for storage

source

sanitized version of values, suitable for storing as output

source
PALEOboxes.init_values!Function
init_values!(
    values, field_data::Type{<:AbstractData}, data_dims::Tuple{Vararg{NamedDimension}}, space::Type{<:AbstractSpace},
    init_value::Symbol, attribv::VariableBase, convertfn, convertvalues, cellrange, info::NTuple{3, String}
)

Initialize values at model start to init_value over region cellrange using information from Variable attribv attributes, scaled by convertfn and convertvalues.

Optional: only required if this field_data type is used for a model (state) Variable that requires initialisation.

Arguments:

  • values: data to be zeroed
  • init_value::Symbol: one of :initialvalue, :normvalue, requesting type of initial value required
  • attribv::VariableBase: Variable with attributes to use for initialisation
  • convertfn::Function: apply multiplier convertfn(convertvalues, i) to initialisation value for cell i. Typically this is used to convert units eg concentration to mol.
  • convertvalues: parameters (if any) required by convertfn, eg a volume measure.
  • cellrange: range of cells to initialise
  • info::::NTuple{3, String}: Tuple (varinfo, convertinfo, trsfrinfo) of identifier strings to use for log messages
source
PALEOboxes.copyfieldto!Function
copyfieldto!(
    dest,
    doff, 
    values, 
    field_data::Type{<:AbstractData}, 
    data_dims::Tuple{Vararg{NamedDimension}}, 
    space::Type{<:AbstractSpace}, 
    cellrange
) -> num_copied::Int

Copy Field.values values from spatial region defined by cellrange, to dest Array starting at index doff.

Number of values over whole Domain should equal degrees-of-freedom returned by dof_values

Required if this field_data type needs to provide values for a numerical solver.

source
PALEOboxes.copytofield!Function
copytofield!(
    values,
    field_data::Type{<:AbstractData},
    data_dims::Tuple{Vararg{NamedDimension}},
    space::Type{<:AbstractSpace},
    cellrange, 
    src, 
    soff
) -> num_copied::Int

Copy from src Array starting at index soff to Field.values values for spatial region defined by cellrange.

Number of values over whole Domain should equal degrees-of-freedom returned by dof_values

Required if this field_data type needs to provide values for a numerical solver.

source
PALEOboxes.add_field_vec!Function
add_field_vec!(
    dest, 
    field_data::Type{<:AbstractData}, 
    data_dims::Tuple{Vararg{NamedDimension}}, 
    space::Type{<:AbstractSpace}, 
    a, 
    cellrange,
    src, 
    soff
) -> num_added::Int

Implement dest += a*src where dest is a Field.values, src is an Array, a is a number, over region defined by cellrange, starting at index soff in src.

Returns number of elements of src used.

See copytofield!, copyfieldto! for the relationship between Array src and Field values dest.

source
PALEOboxes.num_componentsMethod
num_components(field_data::Type{<:AbstractData}) -> Int

get number of components (optional - implement if field_data has a representation as components)

source
PALEOboxes.get_componentsFunction
get_components(values, field_data::Type{<:AbstractData}) -> Vector

Convert Field values to a Vector of components (optional - implement if field_data has a representation as components)

source

Isotopes

PALEOboxes.AbstractIsotopeScalarType
AbstractIsotopeScalar <: AbstractData

An IsotopeScalar represents a quantity or flux with isotopic composition. Can be added, subtracted, multiplied by a scalar, and decomposed into components with the same (bulk) transport properties.

Implementation

Each IsotopeScalar should be added to IsotopeTypes, and implement:

  • get_total(is::IsotopeScalar) -> total
  • arithmetic operations +/- (ie IsotopeScalars can be added and subtracted) and *number, /number (IsotopeScalars can be scaled by a real number).
  • methods of the AbstractData interface

and optional isotope-specific functions, eg for a single isotope:

  • isotope_totaldelta(::Type{<: IsotopeScalar}, total, delta) -> IsotopeScalar()
  • get_delta(is::IsotopeScalar) -> delta
source
PALEOboxes.IsotopeLinearType
IsotopeLinear <: AbstractIsotopeScalar

Linearized representation of isotopic composition, where moldelta = total * delta.

source
PALEOboxes.isotope_totaldeltaFunction
isotope_totaldelta(::Type{IsotopeLinear}, total, delta) -> IsotopeLinear

Create an IsotopeLinear from total and delta

Examples:

julia> a = PB.isotope_totaldelta(PB.IsotopeLinear, 10.0, -2.0)
(v=10.0, v_moldelta=-20.0, ‰=-2.0)

julia> b = a*2
(v=20.0, v_moldelta=-40.0, ‰=-2.0)

julia> c = a + b
(v=30.0, v_moldelta=-60.0, ‰=-2.0)

julia> PB.get_total(c)
30.0

julia> PB.get_delta(c)
-2.0
source