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
PALEOboxes.Model — TypeModelA biogeochemical model consisting of Domains, created from a YAML configuration file using create_model_from_config.
Domains
PALEOboxes.Domain — TypeDomainA 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.
Dimensions and Coordinates
PALEO approximately follows the Common Data Model used by NetCDF and other geoscience data formats.
Domains provide named dimensions, to which coordinate variables may be attached (these are just normal PALEO variables).
- Domain spatial dimensions are provided by the Domain grid (see
Grids) - Additional Domain data dimensions (eg a wavelength grid) may be set by Reactions (see
set_data_dimension!).
PALEOboxes.NamedDimension — TypeNamedDimension(name, size)A named dimension
PALEOboxes.get_dimensions — Functionfunction get_dimensions(obj) -> Vector{NamedDimension}Get all dimensions for PALEO object obj
PALEOboxes.get_dimension — Functionfunction get_dimension(obj, dimname) -> NamedDimensionGet all dimension dimname for PALEO object obj
PALEOboxes.get_coordinates — Functionfunction get_coordinates(obj, dimname) -> coordinates::VectorGet coordinates (if any) attached to dimname for PALEO object obj
PALEO convention is that where possible coordinates contains:
- three variable names, for cell midpoints, lower edges, upper edges, in that order.
- two variable names, for cell midpoints and a bounds variable (with a bounds dimension as the first dimensions)
PALEOboxes.set_coordinates! — Functionfunction set_coordinates!(obj, dimname, coordinates::Vector{String})Set coordinates (variable names) attached to dimname for PALEO object obj
PALEO convention is that where possible coordinates contains:
- three variable names, for cell midpoints, lower edges, upper edges, in that order.
- two variable names, for cell midpoints and a bounds variable (with a bounds dimension as the first dimensions)
Grids
PALEOboxes.AbstractMesh — TypeAbstractMeshDefines additional geometric and topological information for PB.Domain
Concrete subtypes should implement methods:
available_spaces tuple of PB.AbstractSpace types supported.
PB.has_internal_cartesian true if subtype uses an internal spatial representation that differs from the external (cartesian) representation.
PB.internal_size, optionally PB.cartesian_size
PB.Grids.set_subdomain!, PB.Grids.get_subdomain
PB.Grids.create_default_cellrange
Internal and cartesian representations
If a subtype uses an internal representation that differs from the external (cartesian) representation, it should define PB.has_internal_cartesian true and implement cartesian_to_internal, internal_to_cartesian, and PB.cartesian_size.
PALEOboxes.Grids.available_spaces — Functionavailable_spaces(grid::PB.AbstractMesh) -> NTuple{PB.Space}Tuple of Spaces supported by this grid
PALEOboxes.has_internal_cartesian — Functionhas_internal_cartesian(mesh::AbstractMesh, Space::Type{<:AbstractSpace})::Booltrue if Mesh uses different internal and external (cartesian) spatial array representations for Space
PALEOboxes.internal_size — Functioninternal_size(Space::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
PALEOboxes.cartesian_size — Functioncartesian_size(Space::Type{<:AbstractSpace}, 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.
PALEOboxes.Grids.cartesian_to_internal — Functioncartesian_to_internal(grid::PB.AbstractMeshOrNothing, griddata::AbstractArray)Map a cartesian (external) array to internal representation
cartesiantointernal(grid::CartesianLinearGrid, griddata::AbstractArray) -> lindata::Vector
Convert Cartesian Array griddata to Vector on grid.linear_index
PALEOboxes.Grids.internal_to_cartesian — Functioncartesian_to_internal(grid::PB.AbstractMeshOrNothing, griddata::AbstractArray)Map an internal array to cartesian (external) representation
internaltocartesian(grid::CartesianLinearGrid, internaldata::AbstractArray [,missing_value=missing]) -> griddata::Array
Convert Array internaldata (first index on grid.linear_index, up to 2 additional indices for additional data dimensions) to Cartesian Array griddata (with missing_value where no data)
PALEOboxes.Grids.set_subdomain! — Functionset_subdomain!(grid::PB.AbstractMesh, subdomainname::AbstractString, subdom::PB.AbstractSubdomain, allowcreate::Bool=false)Set Subdomain
PALEOboxes.Grids.get_subdomain — Functionget_subdomain(grid::PB.AbstractMesh, subdomainname::AbstractString) -> PB.AbstractSubdomainGet Subdomain
PALEOboxes.Grids.UnstructuredVectorGrid — TypeUnstructuredVectorGrid <: PB.AbstractMeshMinimal Grid for a Vector Domain, defines only some named cells for plotting
PALEOboxes.Grids.UnstructuredColumnGrid — TypeUnstructuredColumnGrid <: PB.AbstractMeshMinimal Grid for a Vector Domain composed of columns (not necessarily forming a 2-D array).
Fields
ncells::Inttotal number of cells in this DomainIcolumns::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.
columnnames::Vector{Symbol}: optional column namescoordinates::Dict{String, Vector{String}}: optional attached coordinates
PALEOboxes.Grids.CartesianLinearGrid — TypeCartesianLinearGrid <: PB.AbstractMeshnD 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 inprod(dims))dimensions::Vector{PB.NamedDimensions}: names and sizes of cartesian spatial dimensions (ordered list)coordinates::Dict{String, Vector{String}}: optional attached coordinates
PALEOboxes.Grids.CartesianArrayGrid — TypeCartesianArrayGrid <: PB.AbstractMeshnD grid with netcdf CF1.0 coordinates, using n-dimensional Arrays for PALEO internal representation of Variables
Fields
ncells::Int64: number of cells in Domain = product of cartesian dimension sizesdimensions::Vector{PB.NamedDimensions}: names and sizes of cartesian spatial dimensions (ordered list)coordinates::Dict{String, Vector{String}}: optional attached coordinates
Subdomains
PALEOboxes.AbstractSubdomain — TypeAbstractSubdomainDefines the relationship between two PB.Domains by mapping indices in one Domain to related indices in the other, eg interior cells adjacent to a boundary.
Concrete subtypes should implement:
PALEOboxes.Grids.BoundarySubdomain — TypeBoundarySubdomain <: PB.AbstractSubdomainA 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 cellibnd.
PALEOboxes.Grids.InteriorSubdomain — TypeInteriorSubdomain <: PB.AbstractSubdomainA 3D subdomain corresponding to the 3D interior Domain associated with a 2D boundary Domain: indices[iint] is either:
missingifiintis 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
PALEOboxes.Grids.subdomain_view — Functionsubdomain_view(values, subdomain::BoundarySubdomain) -> viewCreate a view on values in an interior Domain to access cells corresponding to indices in a boundary Domain.
subdomain_view(values, subdomain::InteriorSubdomain) -> varReturn unmodified values from a boundary Domain, indices to access from interior supplied by subdomain_indices
subdomain_view(values, subdomain::Nothing) -> valuesFallback when subdomain == nothing
PALEOboxes.Grids.subdomain_indices — Functionsubdomain_indices(subdomain::BoundarySubdomain) -> nothingNo additional indices required to access Variables in an interior Domain from a boundary Domain (view created by subdomain_view is sufficient).
subdomain_indices(subdomain::InteriorSubdomain) -> indicesReturn indices to access Variables in a boundary Domain from interior Domain (will have missing entries where interior cells do not correspond to boundary)
subdomain_indices(subdomain::Nothing) -> nothingfallback when subdomain == nothing
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:
VariableReactions ofVariableTypeVT_ReactPropertyandVT_ReactDependencywhich are linked to createVariableDomPropDepsVariableReactions ofVariableTypeVT_ReactTargetandVT_ReactContributorwhich are linked to createVariableDomContribTargets.
PALEOboxes.VariableBase — TypeVariableBaseA Model biogeochemical Variable. Reactions access Variables using derived Types VariableReaction which are links to VariableDomains.
PALEOboxes.VariableType — Type@enum VariableTypeEnumeration of VariableBase subtypes. Allowed values:
VariableReaction:VT_ReactProperty,VT_ReactDependency,VT_ReactContributor,VT_ReactTargetVariableDomain:VT_DomPropDep,VT_DomContribTarget
PALEOboxes.VariableDomain — TypeVariableDomain <: VariableBaseAbstract base Type for a (Domain) model variable.
Defines a named Variable and corresponding data Fields that are linked to by VariableReactions.
PALEOboxes.VariableDomPropDep — TypeVariableDomPropDep <: VariableDomainModel (Domain) VariableDomain linking a single VariableReaction{VT_ReactProperty} to multiple VariableReaction{VT_ReactDependency}.
PALEOboxes.VariableDomContribTarget — TypeVariableDomContribTarget <: VariableDomainModel (Domain) VariableDomain linking a single VariableReaction{VT_ReactTarget} to multiple VariableReaction{VT_ReactContributor} and VariableReaction{VT_ReactDependency}.
Variable Attributes
PALEOboxes.Attribute — TypeAttribute{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)
PALEOboxes.StandardAttributes — ConstantStandardAttributesList of standard Variable attributes.
Some of these follow netCDF COARDS/CF conventions:
COARDS:https://ferret.pmel.noaa.gov/Ferret/documentation/coards-netcdf-conventions
- units (where possible should follow the Unidata udunits package https://docs.unidata.ucar.edu/udunits/current/)
- long_name
CF conventions: https://cfconventions.org/cf-conventions/cf-conventions.html#_description_of_the_data
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 cellsPALEOboxes.get_attribute — Functionget_attribute(var::VariableBase, name::Symbol, missing_value=missing) -> valueGet Variable attribute.
PALEOboxes.VariableFunction — Type@enum VariableFunctionAllowed 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.
PALEOboxes.VariablePhase — Type@enum VariablePhaseAllowed values of :vphase Variable Attribute, defining the component phase this Variable belongs to for multiphase cells.
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
:field_dataand optionally:data_dimsdefine the subtype ofAbstractData:spacedefines the subtype ofAbstractSpace
Examples:
:field_data=ScalarData,:space=ScalarSpacedefines a Domain scalar (0D) quantity.:field_data=IsotopeLinear,:space=ScalarSpacedefines a Domain scalar (0D) isotope quantity.:field_data=ScalarData,:space=CellSpacedefines a per-cell quantity in a spatial Domain.:field_data=ArrayScalarData,:data_dims =("wgrid",):space=ScalarSpacedefines a Domain scalar (0D) quantity on a wavelength grid, a Vector of length given by the value of the Domain "wgrid" data dimension (seeset_data_dimension!)
PALEOboxes.Field — TypeField{FieldData <: AbstractData, Space <: AbstractSpace, V, N, Mesh <: AbstractMeshOrNothing}
Field(FieldData, data_dims::NTuple{N, NamedDimensions}, data_type, Space, mesh::AbstractMeshOrNothing; allocatenans)
Field(existing_values, FieldData, data_dims::NTuple{N, NamedDimension}, data_type::Union{DataType, Missing}, Space, mesh::AbstractMeshOrNothing)A Field of values::V of type FieldData <: AbstractData defined on function space Space over grid of type Mesh and (optionally) with N data_dims::NTuple{N, NamedDimensions}.
values::V is usually Array or similar of elements with Type data_type and dimensions defined by Space, Mesh and data_dims, and is either created by allocate_values or supplied by existing_values.
PALEOboxes.get_field — Functionget_field(obj, ...) -> FieldGet Field from PALEO object obj
Spaces
PALEOboxes.AbstractSpace — TypeAbstractSpaceDefines a function Space within a Domain, on a mesh defined by a Grid
PALEOboxes.ScalarSpace — TypeScalarSpace <: AbstractSpaceA Domain position-independent quantity
PALEOboxes.CellSpace — TypeCellSpace <: AbstractSpaceA per-cell quantity. Use as Variable attribute :space to create a Variable with data array dimensions from Grid cells.
PALEOboxes.ColumnSpace — TypeColumnSpace <: AbstractSpaceA per-column quantity. Use as Variable attribute :space to create a Variable with data array dimensions from Grid columns.
Data types
PALEOboxes.AbstractData — TypeAbstractDataDefines 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
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
If the subtype needs to provide a thread-safe atomic addition operation eg to provide scalar accumulator variables for Domain totals with a tiled model, it should implement atomic_add! for the field values.
PALEOboxes.UndefinedData — TypeUndefinedData <: AbstractDataUndefined data type (no methods implemented). Used to indicate that a Variable can link to any data type.
PALEOboxes.ScalarData — TypeScalarData <: AbstractDataA Field scalar (eg a biogeochemical concentration)
PALEOboxes.ArrayScalarData — TypeArrayScalarData <: AbstractDataAn 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.
PALEOboxes.allocate_values — Functionallocate_values(
FieldData::Type{<:AbstractData}, data_dims::Tuple{Vararg{NamedDimension}}, data_type, Space::Type{<:AbstractSpace}, spatial_size::Tuple;
allocatenans::Bool,
) -> values::Vallocate Field.values::V (eg an Array of elements with Type data_type) for FieldData with dimensions defined by spatial_size and data_dims
PALEOboxes.check_values — Functioncheck_values(
existing_values,
FieldData::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.
PALEOboxes.zero_values! — Functionzero_values!(values, FieldData::Type{<:AbstractData}, data_dims::Tuple{Vararg{NamedDimension}}, Space::Type{<:AbstractSpace}, cellrange)Set values over spatial region cellrange to zero at start of main loop
PALEOboxes.dof_values — Functiondof_values(
FieldData::Type{<:AbstractData},
data_dims::Tuple{Vararg{NamedDimension}},
Space::Type{<:AbstractSpace},
mesh,
cellrange
) -> dof::IntReturn degrees-of-freedom for FieldData over spatial region cellrange.
PALEOboxes.init_values! — Functioninit_values!(
values, FieldData::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 FieldData type is used for a model (state) Variable that requires initialisation.
Arguments:
values: data to be zeroedinit_value::Symbol: one of :initialvalue, :normvalue, requesting type of initial value requiredattribv::VariableBase: Variable with attributes to use for initialisationconvertfn::Function: apply multiplierconvertfn(convertvalues, i)to initialisation value for cell i. Typically this is used to convert units eg concentration to mol.convertvalues: parameters (if any) required byconvertfn, eg a volume measure.cellrange: range of cells to initialiseinfo::::NTuple{3, String}: Tuple (varinfo, convertinfo, trsfrinfo) of identifier strings to use for log messages
PALEOboxes.copyfieldto! — Functioncopyfieldto!(
dest,
doff,
values,
FieldData::Type{<:AbstractData},
data_dims::Tuple{Vararg{NamedDimension}},
Space::Type{<:AbstractSpace},
cellrange
) -> num_copied::IntCopy 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 FieldData type needs to provide values for a numerical solver.
PALEOboxes.copytofield! — Functioncopytofield!(
values,
FieldData::Type{<:AbstractData},
data_dims::Tuple{Vararg{NamedDimension}},
Space::Type{<:AbstractSpace},
cellrange,
src,
soff
) -> num_copied::IntCopy 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 FieldData type needs to provide values for a numerical solver.
PALEOboxes.add_field! — Functionadd_field!(obj, f::Field ...)Add Field or Field to PALEO object obj
PALEOboxes.add_field_vec! — Functionadd_field_vec!(
dest,
FieldData::Type{<:AbstractData},
data_dims::Tuple{Vararg{NamedDimension}},
Space::Type{<:AbstractSpace},
a,
cellrange,
src,
soff
) -> num_added::IntImplement 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.
PALEOboxes.num_components — Methodnum_components(FieldData::Type{<:AbstractData}) -> Intget number of components (optional - implement if FieldData has a representation as components)
PALEOboxes.get_components — Functionget_components(values, FieldData::Type{<:AbstractData}) -> VectorConvert Field values to a Vector of components (optional - implement if FieldData has a representation as components)
Isotopes
PALEOboxes.AbstractIsotopeScalar — TypeAbstractIsotopeScalar <: AbstractDataAn 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
AbstractDatainterface
and optional isotope-specific functions, eg for a single isotope:
- isotope_totaldelta(::Type{<: IsotopeScalar}, total, delta) -> IsotopeScalar()
- get_delta(is::IsotopeScalar) -> delta
PALEOboxes.IsotopeLinear — TypeIsotopeLinear <: AbstractIsotopeScalarLinearized representation of isotopic composition, where moldelta = total * delta.
PALEOboxes.isotope_totaldelta — Functionisotope_totaldelta(::Type{IsotopeLinear}, total, delta) -> IsotopeLinearCreate 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.0PALEOboxes.get_total — Functiongeneric get_total for non-isotope variable
get_total(is::IsotopeLinear)Get isotope total
PALEOboxes.get_delta — Functionget_delta(is::IsotopeLinear)Get isotope delta (per mil)