CFDomains

Documentation for CFDomains.

General

CFDomains.HVLayoutType
struct HVLayout{rank} end
layout = HVLayout(rank)

Singleton type describing a multi-layer data layout where horizontal layers are contiguous in memory. rank is the number of horizontal indices (1 or 2).

source
CFDomains.HyperDiffusionMethod
filter = HyperDiffusion(fieldtype::Symbol, domain, niter::Int, nu)                  # user-friendly
filter = HyperDiffusion{fieldtype, D, F, X}(domain::D, niter::Int, nu::F, extra::X) # internal

Return a filter that applies Laplacian diffusion iterated niter times with hyperdiffusive coefficient nu ≥ 0 on fields of type fieldtype. Supported field types are :scalar, :vector, :vector_curl, :vector_div.

The filter is to be used as:

# out-of-place, allocates `scratch`, allocates and returns `filtered_data`
filtered_data = filter(void, in_data, void)

# mutating, non-allocating
out_data = filter(out_data, in_data, scratch)

# in-place, non-allocating
inout_data = filter(inout_data, inout_data, scratch)

filter(out_data, in_data, scratch) returns hyperdiffusion!(filter.domain, filter, out_data, in_data, scratch). If scratch::Void, then

scratch = scratch_space(filter, in_data)

The user-friendly HyperDiffusion may be specialized for specific types of domain. The specialized method should call the internal constructor which takes an extra argument extra, stored as filter.extra for later use by hyperdiffusion!. The default user-friendly constructor sets extra=nothing.

Similarly, scratch_space(filter, in_data) calls scratch_hyperdiff(f.domain, Val(fieldtype), in_data) which by default returns nothing. scratch_hyperdiff may be specialized for specific domain types and fieldtype.

source
CFDomains.ShellType
multi_layer_domain = Shell(nz::Int, layer::AbstractDomain, layout)

Return a multi-layer domain made of nz layers with data layout specified by layout. Unless you know what you are doing, it is recommended to use rather:

multi_layer_domain = shell(nz::Int, layer::AbstractDomain)

which gets the data layout from data_layout(layer). Otherwise, multi_layer_domain may be non-optimal or non-usable.

source
CFDomains.SigmaCoordinateType
sigma = SigmaCoordinate(N, ptop) <: PressureCoordinate{N}

Pressure based sigma-coordinate for N levels with top pressure ptop. Pressure levels are linear in vertical coordinate k : k/N = (ps-p)/(ps-ptop) where k ranges from 0 (ground) to N (model top).

source
CFDomains.VHLayoutType
struct VHLayout{rank} end
layout = VHLayout(rank)

Singleton type describing a multi-layer data layout where vertical columns are contiguous. rank is the number of horizontal indices (1 or 2).

source
CFDomains.allocate_fieldFunction
field = allocate_field(kind::Symbol, domain::AbstractDomain, precision::Type)

Allocate a field of the given kind and precision over the given domain. Typical values for precision are Float32, Float64 or ForwardDiff.Dual. Depending on the domain, valid values for kind may include :scalar, :vector, :scalar_spec, scalar_spat (spectral/spatial representation of a scalar field), :vector_spec, vector_spat (spectral/spatial representation of a vector field).

Internally, allocate_field(kind::Symbol, domain, F) returns allocate_field(Val(kind), domain, F). To specialize allocate_field for MyDomain <: AbstractDomain, one must provide methods for:

allocate_field(::Val{kind}, domain::MyDomain, F)

where symbol kind::Symbol is one of the valid field kinds for that domain.

source
CFDomains.allocate_fieldsFunction
fields = allocate_fields(kinds::Tuple, domain::AbstractDomain, F::Type)
fields = allocate_fields(kinds::NamedTuple, domain::AbstractDomain, F::Type)

Allocate a (named) tuple of fields according to the provided kinds. For instance:

fields = allocate_fields((:vector, :scalar), domain, F)
fields = allocate_fields((a=:vector, b=:scalar), domain, F)

are equivalent to, respectively:

fields = (allocate_field(:vector, domain, F), allocate_field(:scalar, domain, F))
fields = (a=allocate_field(:vector, domain, F), b=allocate_field(:scalar, domain, F))
source
CFDomains.data_layoutMethod
layout = data_layout(shell::Shell)

Return layout describing the data layout of multi-layer domain shell.

layout = data_layout(domain::Domain)

Return layout describing the preferred data layout for a shell made of layers of type Domain.

Typical values for layout are the singletons HVLayout() (layers are contiguous in memory) and VHLayout (columns are contiguous in memory).

source
CFDomains.hyperdiff_shell!Function
new_data = hyperdiff_shell!(layer, layout, filter, out_data, in_data, scratch)

Apply hyperdiffusive filtertoindata. The result is written intooutdataand returned. Ifdomain::SpectralDomain, the filter applies to spectral coefficients. Ifout_data::Void, it is adequately allocated. SeeHyperDiffusion.

source
CFDomains.hyperdiffusion!Method
new_data = hyperdiffusion!(domain, filter, out_data, in_data, scratch)

Apply hyperdiffusive filtertoindata. The result is written intooutdataand returned. Ifdomain::SpectralDomain, the filter applies to spectral coefficients. Ifout_data::Void, it is adequately allocated. SeeHyperDiffusion.

source
CFDomains.mass_coordinateFunction
mcoord = mass_coordinate(pcoord::PressureCoordinate, metric_cov=1)

Return the mass-based coordinate deduced from pcoord and the covariant metric factor metric_cov. metric_cov can be a scalar or a vector. In the latter case metric_cov[ij]) is the metric factor at horizontal positionij. The objectmcoord` can then be used with:

m = mass_level(k, ij, masstot, mcoord)

With metric_cov==1, masstot should be in Pa ( kg/m²⋅(m/s²) ). With metric_cov in m², masstot should be in kg⋅(m/s²). m has the same unit as mass_tot. If masstot is covariant (integral over a cell) then metric_cov should include the cell area. See also mass_level.

source
CFDomains.mass_coordinateMethod
mcoord = mass_coordinate(mcoord::MassCoordinate)

Return mcoord itself, unchanged. Interim function for backwards compatibility.

source
CFDomains.mass_levelFunction
m = mass_level(k, ij, masstot, mcoord::MassCoordinate{N})

Return mass m in level k/2 and at horizontal position ij, as prescribed by vertical coordinate mcoord and total mass masstot.

So-called full levels correspond to odd values k=1,2...2N-1 while interfaces between full levels (so-called half-levels) correspond to even values k=0,2...2N

masstot may be:

  • per unit surface, with unit Pa
  • per unit non-dimensional surface (e.g. on the unit sphere), with unit kg.(m/s²))

Which convention is appropriate depends on the metric_factor provided when constructing mcoord. See also mass_coordinate.

source
CFDomains.periodize!Method
periodize!(data, box::AbstractBox, mgr)

Enforce horizontally-periodic boundary conditions on array data representing grid point values in box. data may also be a collection, in which case periodize! is applied to each element of the collection. Call periodize! on data obtained by computations involving horizontal averaging/differencing.

source
CFDomains.pressure_levelFunction
p = pressure_level(k, ps, vcoord::PressureCoordinate{N})

Returns pressure p corresponding to level k/2 as prescribed by vertical coordinate vcoord and surface pressure ps.

So-called full levels correspond to odd values k=1,2...2N-1 while interfaces between full levels (so-called half-levels) correspond to even values k=0,2...2N

source
CFDomains.shellMethod
multi_layer_domain = shell(nz::Int, layer::AbstractDomain)

Return a multi-layer domain made of nz layers with data layout specified by data_layout(layer).

source
CFDomains.transpose!Method
x_ji = transpose!(x_ji, mgr, x_ij)
y_ji = transpose!(void, mgr, y_ij)

Transposes x_ij and writes the result into x_ji, which may be ::Void, in which case it is allocated.

When working with shells it is sometimes useful to transpose fields for performance. CFDomains.transpose! can be specialized for specific managers, for instance:

import CFDomains: transpose!, Void
using Strided: @strided
function transpose!(x, ::MultiThread, y)
   @strided permutedims!(x, y, (2,1))
   return x # otherwise returns a StridedView
end
transpose!(::Void, ::MultiThread, y) = permutedims(y, (2,1)) # for non-ambiguity
source

Vertical interpolation

CFDomains.VerticalInterpolation.interpolate!Function
interpolated = interpolate!(mgr::LoopManager, domain::Shell, field, coord, refs, increasing)

Interpolate field defined on 3D domain to reference values refs of coord. If coord increases with level number, increasing==true and vice-versa. refs must be sorted according to increasing. out may be void. Example:

T_ref = interpolate!(mgr, void, domain, [850., 500.], temperature, pressure, false)

With these arguments, interpolate! is not meant to be specialized for a specific domain type. Instead, specialize: interpolated = interpolate!(mgr, out, layout::Layout, refs, field, coord, increasing::Val)

source

Lazy expressions and operators

CFDomains.LazyExpressions.pdvFunction
fa = pdv(fun, a)
fa, fb = pdv(fun, a, b)
fa, fb, fc = pdv(fun, a, b, c)

Return the partial derivatives of scalar function fun evaluated at input a, .... This function is implemented only when the package ForwardDiff is loaded either directly from the main program or via some dependency.

source

Virtual zero-filled arrays

CFDomains.ZeroArrays.zero_arrayMethod
z = zero_array(a::AbstractArray)

Return z which behaves like a read-only array with the same axes as a, filled with zeros.

z[i...] does not check bounds and always returns the special value Zero(), which behaves like 0 when passed to +, * and muladd.

source