Fields
ClimaCore.Fields.Field — TypeField(values, space)A set of values defined at each point of a space.
ClimaCore.Fields.coordinate_field — Functioncoordinate_field(space::AbstractSpace)Return a pointer to the input space's coordinates Field.
ClimaCore.Fields.local_geometry_field — Functionlocal_geometry_field(space::AbstractSpace)Return a pointer to the input space's LocalGeometry Field.
Base.zeros — Methodzeros(space::AbstractSpace)Create a new field on space that is zero everywhere.
Base.ones — Methodones(space::AbstractSpace)Create a new field on space that is one everywhere.
Base.sum — Methodsum([f=identity,] v::Field)Approximate integration of v or f.(v) over the domain. In an AbstractSpectralElementSpace, an integral over the entire space is computed by summation over the elements of the integrand multiplied by the Jacobian determinants and the quadrature weights at each node within an element. Hence, sum is computed by summation of the field values multiplied by the Jacobian determinants and quadrature weights:
\[\sum_i f(v_i) W_i J_i \approx \int_\Omega f(v) \, d \Omega\]
where $v_i$ is the value at each node, and $f$ is the identity function if not specified.
If v is a distributed field, this uses a ClimaComms.allreduce operation.
ClimaCore.Fields.local_sum — FunctionFields.local_sum(v::Field)Compute the approximate integral of v over the domain local to the current context.
See sum for the integral over the full domain.
Statistics.mean — Methodmean([f=identity, ]v::Field)The mean value of field or f.(field) over the domain, weighted by area. Similar to sum, in an AbstractSpectralElementSpace, this is computed by summation of the field values multiplied by the Jacobian determinants and quadrature weights:
\[\frac{\sum_i f(v_i) W_i J_i}{\sum_i W_i J_i} \approx \frac{\int_\Omega f(v) \, d \Omega}{\int_\Omega \, d \Omega}\]
where $v_i$ is the Field value at each node, and $f$ is the identity function if not specified.
If v is a distributed field, this uses a ClimaComms.allreduce operation.
LinearAlgebra.norm — Methodnorm(v::Field, p=2; normalize=true)The approximate $L^p$ norm of v, where $L^p$ represents the space of measurable functions for which the p-th power of the absolute value is Lebesgue integrable, that is:
\[\| v \|_p = \left( \int_\Omega |v|^p d \Omega \right)^{1/p}\]
where $|v|$ is defined to be the absolute value if $v$ is a scalar-valued Field, or the 2-norm if it is a vector-valued Field or composite Field (see LinearAlgebra.norm). Similar to sum and mean, in an AbstractSpectralElementSpace, this is computed by summation of the field values multiplied by the Jacobian determinants and quadrature weights. If normalize=true (the default), then internally the discrete norm is divided by the sum of the Jacobian determinants and quadrature weights:
\[\left(\frac{\sum_i |v_i|^p W_i J_i}{\sum_i W_i J_i}\right)^{1/p} \approx \left(\frac{\int_\Omega |v|^p \, d \Omega}{\int_\Omega \, d \Omega}\right)^{1/p}\]
If p=Inf, then the norm is the maximum of the absolute values
\[\max_i |v_i| \approx \sup_{\Omega} |v|\]
Consequently all norms should have the same units for all $p$ (being the same as calling norm on a single value).
If normalize=false, then the denominator term is omitted, and so the result will be the norm as described above multiplied by the length/area/volume of $\Omega$ to the power of $1/p$.
ClimaCore.Fields.set! — Functionset!(f::Function, field::Field, args = ())Apply function f to populate values in field field. f must have a function signature with signature f(::LocalGeometry[, args...]). Additional arguments may be passed to f with args.
Example
using ClimaCore.Fields
using ClimaCore.CommonSpaces
ᶜspace = ExtrudedCubedSphereSpace(Float64;
z_elem = 10,
z_min = 0,
z_max = 1,
radius = 10,
h_elem = 10,
n_quad_points = 4,
staggering = CellCenter(),
)
x = Fields.Field(Float64, ᶜspace)
Fields.set!(x) do lg
sin(lg.coordinates.z)
endClimaCore.Grids.ColumnIndex — TypeColumnIndex(ij,h)An index into a column of a field. This can be used as an argument to getindex of a Field, to return a field on that column.
Example
colidx = ColumnIndex((1,1),1)
field[colidx]ClimaCore.Fields.bycolumn — FunctionFields.bycolumn(fn, space)Call fn(colidx) to every ColumnIndex colidx of space. This can be used to apply multiple column-wise operations in a single pass, making use of multiple threads.
On GPUs this will simply evaluate f once with colidx=: (i.e. it doesn't perform evaluation by columns). This may change in future.
Example
∇ = GradientF2C()
div = DivergenceC2F()
bycolumn(axes(f)) do colidx
@. ∇f[colidx] = ∇(f[colidx])
@. df[colidx] = div(∇f[colidx])
endClimaCore.Fields.Δz_field — FunctionΔz_field(field::Field)
Δz_field(space::AbstractSpace)Return a pointer to the input space's Field containing the Δz values on the same space as the given field.