Spaces
A Space represents a discretized function space over some domain. Currently two main discretizations are supported: Spectral Element Discretization (both Continuous Galerkin and Discontinuous Galerkin types) and a staggered Finite Difference Discretization. Combination of these two in the horizontal/vertical directions, respectively, is what we call a hybrid space.
Sketch of a 2DX hybrid discretization:

ClimaCore.Spaces — Module
ClimaCore.Spaces.Δz_data — Function
Finite Difference Spaces
ClimaCore.jl supports staggered Finite Difference discretizations. Finite Differences discretize an interval domain by approximating the function by a value at either the center of each element (also referred to as cell) (CenterFiniteDifferenceSpace), or the interfaces (faces in 3D, edges in 2D or points in 1D) between elements (FaceFiniteDifferenceSpace).
ClimaCore.Spaces.FiniteDifferenceSpace — Type
FiniteDifferenceSpace(
grid::Grids.FiniteDifferenceGrid,
staggering::Staggering
)A 1D finite-difference space, that lives on either:
- cell centers (where
staggeringisGrids.CellCenter) or - cell faces (where
staggeringisGrids.CellFace)
Users should construct either the center or face space from the mesh, then construct the other space from the original one: this internally reuses the same data structures, and avoids allocating additional memory.
Spectral Element Spaces
ClimaCore.Spaces.SpectralElementSpace1D — Type
SpectralElementSpace1D(grid::SpectralElementGrid1D)
SpectralElementSpace1D(
topology::Topologies.IntervalTopology,
quadrature_style::Quadratures.QuadratureStyle;
kwargs...
)sourceClimaCore.Spaces.SpectralElementSpace2D — Type
SpectralElementSpace2D(grid::SpectralElementGrid1D)
SpectralElementSpace2D(
topology::Topologies.Topology2D,
quadrature_style::Quadratures.QuadratureStyle;
kwargs...,
)sourceClimaCore.Spaces.SpectralElementSpaceSlab — Type
SpectralElementSpaceSlab <: AbstractSpaceA view into a SpectralElementSpace2D for a single slab.
ClimaCore.Spaces.node_horizontal_length_scale — Function
Spaces.node_horizontal_length_scale(space::AbstractSpectralElementSpace)The approximate length scale of the distance between nodes. This is defined as the length scale of the mesh (see Meshes.element_horizontal_length_scale), divided by the number of unique quadrature points along each dimension.
Extruded Finite Difference Spaces
ClimaCore.Spaces.ExtrudedFiniteDifferenceSpace — Type
ExtrudedFiniteDifferenceSpace(grid, staggering)
ExtrudedFiniteDifferenceSpace(
horizontal_space::AbstractSpace,
vertical_space::FiniteDifferenceSpace,
hypsography::Grids.HypsographyAdaption = Grids.Flat();
deep::Bool = false,
)An extruded finite-difference space, where the extruded direction is staggered, containing grid information at either
- cell centers (where
staggeringisGrids.CellCenter) or - cell faces (where
staggeringisGrids.CellFace)
Utilities
ClimaCore.Spaces.area — Function
Spaces.area(space::Spaces.AbstractSpace)The length/area/volume of space. This is computed as the sum of the quadrature weights $W_i$ multiplied by the Jacobian determinants $J_i$:
\[\sum_i W_i J_i \approx \int_\Omega \, d \Omega\]
If space is distributed, this uses a ClimaComms.allreduce operation.
ClimaCore.Spaces.local_area — Function
Spaces.local_area(space::Spaces.AbstractSpace)The length/area/volume of space local to the current context. See Spaces.area