Loading...
Searching...
No Matches
terra::fv::hex Namespace Reference

Namespaces

namespace  operators
 

Classes

struct  DirichletBCs
 Constant Dirichlet boundary conditions for the inner (CMB) and outer (surface) radial boundaries of the spherical shell. More...
 

Concepts

concept  FVProjectionFunctor
 

Functions

template<typename ScalarType , typename GridScalarType , FVProjectionFunctor< ScalarType, GridScalarType > Functor>
void l2_project_analytical_to_fv (linalg::VectorFVScalar< ScalarType > &dst, const Functor &src, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii)
 L2 projection of an analytical function into a finite volume function.
 
template<typename ScalarType , typename GridScalarType >
void l2_project_fv_to_fe (linalg::VectorQ1Scalar< ScalarType > &dst, const linalg::VectorFVScalar< ScalarType > &src, const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii, std::vector< linalg::VectorQ1Scalar< ScalarType > > &tmps)
 L2 projection from a finite volume function into a Q1 (wedge) finite element function.
 
template<typename ScalarType , typename GridScalarType >
void l2_project_fv_to_fe_lumped (linalg::VectorQ1Scalar< ScalarType > &dst, const linalg::VectorFVScalar< ScalarType > &src, const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii, std::vector< linalg::VectorQ1Scalar< ScalarType > > &tmps)
 Bound-preserving FV-to-FE transfer via a lumped mass matrix.
 
template<typename ScalarType , typename GridScalarType >
void l2_project_fe_to_fv (linalg::VectorFVScalar< ScalarType > &dst, const linalg::VectorQ1Scalar< ScalarType > &src, const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii)
 L2 projection from a Q1 finite element function into a finite volume function.
 
template<typename ScalarT >
void apply_dirichlet_bcs (grid::Grid4DDataScalar< ScalarT > data, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &boundary_mask, const DirichletBCs< ScalarT > &bcs, const grid::shell::DistributedDomain &domain)
 Strongly enforce Dirichlet boundary conditions on the radial shell boundaries.
 
template<typename ScalarT >
void apply_dirichlet_bcs (linalg::VectorFVScalar< ScalarT > &T, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &boundary_mask, const DirichletBCs< ScalarT > &bcs, const grid::shell::DistributedDomain &domain)
 
template<typename ScalarType , typename GridScalarType >
void compute_cell_centers (linalg::VectorFVVec< ScalarType, 3 > &dst, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii)
 Computes cell centers and writes to a vector valued finite volume function.
 
template<typename ScalarType , typename GridScalarType >
void initialize_cell_centers (linalg::VectorFVVec< ScalarType, 3 > &dst, const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< GridScalarType, 3 > &coords_shell, const grid::Grid2DDataScalar< GridScalarType > &coords_radii)
 Computes cell center positions once and populates ghost layers via MPI communication.
 

Function Documentation

◆ apply_dirichlet_bcs() [1/2]

template<typename ScalarT >
void terra::fv::hex::apply_dirichlet_bcs ( grid::Grid4DDataScalar< ScalarT >  data,
const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &  boundary_mask,
const DirichletBCs< ScalarT > &  bcs,
const grid::shell::DistributedDomain domain 
)

Strongly enforce Dirichlet boundary conditions on the radial shell boundaries.

For each FV cell that lies on the CMB ( \(r = 1\), innermost real cell of the innermost radial subdomain) or on the outer surface ( \(r = N_r - 1\), outermost real cell of the outermost radial subdomain), the scalar field is overwritten with the prescribed value.

Boundary cells are identified via the boundary_mask:

  • CMB cells: r_cell == 1 and boundary_mask(id, 0, 0, 0) == ShellBoundaryFlag::CMB.
  • Surface cells: r_cell == last and boundary_mask(id, 0, 0, last_q1) == ShellBoundaryFlag::SURFACE.

This strong enforcement is first-order in time (the boundary value is fixed after each step). For cells adjacent to the boundary the scheme is spatially consistent: the next time step's stencil will see the correct Dirichlet value at the boundary cell.

Parameters
T[in/out] Transported scalar field (boundary cells overwritten).
boundary_maskNode-based boundary flag array (Q1 layout, CMB at r=0, SURFACE at r=last).
bcsPrescribed BC values and flags.
domainDistributed domain (used for the loop range policy).

Apply Dirichlet BCs directly to a raw FV grid view.

Sets both the real boundary cell and the adjacent ghost cell outside the physical domain. The ghost cell is never filled by update_fv_ghost_layers (no subdomain neighbour exists beyond a physical boundary) so it must be set here to give the correct inflow/diffusion BC when FCT predictor stencils read across the boundary face.

This overload is the low-level implementation; the VectorFVScalar overload delegates to it. It can also be called directly on intermediate FCT buffers (e.g. T_L) that are raw views.

◆ apply_dirichlet_bcs() [2/2]

template<typename ScalarT >
void terra::fv::hex::apply_dirichlet_bcs ( linalg::VectorFVScalar< ScalarT > &  T,
const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &  boundary_mask,
const DirichletBCs< ScalarT > &  bcs,
const grid::shell::DistributedDomain domain 
)

◆ compute_cell_centers()

template<typename ScalarType , typename GridScalarType >
void terra::fv::hex::compute_cell_centers ( linalg::VectorFVVec< ScalarType, 3 > &  dst,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii 
)

Computes cell centers and writes to a vector valued finite volume function.

This function computes for every finite volume cell \(K\) the value

\[ u_K = \frac{1}{|K|} \int_K x \ dx \]

with

\[ |K| = \int_K 1 \ dx \]

and writes \(u_K\) into the respective finite volume cell.

Parameters
dst[out] finite volume function that is being written to
coords_shellcoords of the unit shell
coords_radiiradii of all shells

◆ initialize_cell_centers()

template<typename ScalarType , typename GridScalarType >
void terra::fv::hex::initialize_cell_centers ( linalg::VectorFVVec< ScalarType, 3 > &  dst,
const grid::shell::DistributedDomain domain,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii 
)

Computes cell center positions once and populates ghost layers via MPI communication.

This is the recommended one-time initialization for cell centers at application startup. After this call, dst contains valid cell centers in both real cells and all ghost layers, so kernels can look up neighbour cell centers without recomputing geometry.

Parameters
dst[out] FV vector field receiving cell centers (real cells + ghost layers).
domainDistributed domain (used for ghost layer communication).
coords_shellLateral node coordinates of the unit sphere surface.
coords_radiiRadial shell radii.

◆ l2_project_analytical_to_fv()

template<typename ScalarType , typename GridScalarType , FVProjectionFunctor< ScalarType, GridScalarType > Functor>
void terra::fv::hex::l2_project_analytical_to_fv ( linalg::VectorFVScalar< ScalarType > &  dst,
const Functor &  src,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii 
)

L2 projection of an analytical function into a finite volume function.

Use this function if you want to represent an analytical function with a finite volume function.

Given some function \(u = u(x)\) (as a FVProjectionFunctor), this function computes for every finite volume cell \(K\) the value

\[ u_K = \frac{1}{|K|} \int_K u(x) \ dx \]

with

\[ |K| = \int_K 1 \ dx \]

and writes \(u_K\) into the respective finite volume cell.

Note
The operator() method of the functor must be annotated with KOKKOS_INLINE_FUNCTION. Otherwise, this might not work on GPUs.
Parameters
dst[out] finite volume function that is being written to
srcfunctor that implements the FVProjectionFunctor concept
coords_shellcoords of the unit shell
coords_radiiradii of all shells

◆ l2_project_fe_to_fv()

template<typename ScalarType , typename GridScalarType >
void terra::fv::hex::l2_project_fe_to_fv ( linalg::VectorFVScalar< ScalarType > &  dst,
const linalg::VectorQ1Scalar< ScalarType > &  src,
const grid::shell::DistributedDomain domain,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii 
)

L2 projection from a Q1 finite element function into a finite volume function.

For each FV cell \(K\) computes the exact volume-weighted cell average of the Q1 field:

\[ u_K^{FV} = \frac{1}{|K|} \int_K u_h^{FE}(x) \, dx = \frac{1}{|K|} \int_K \sum_i u_i \, \phi_i(\xi(x)) \, dx \]

where \(\phi_i\) are the wedge Q1 shape functions and \(u_i\) are the Q1 nodal values surrounding the cell. The integral is evaluated with the same Felippa quadrature rule used throughout the codebase, so the result is consistent with l2_project_fv_to_fe.

Ghost layers of dst are populated via MPI exchange on exit so the result is immediately usable as input to FCT kernels.

Parameters
dst[out] FV scalar field.
srcQ1 finite element scalar field (read-only).
domainDistributed domain (used for MPI ghost-layer exchange).
coords_shellLateral node coordinates of the unit-sphere surface.
coords_radiiRadial shell radii.

◆ l2_project_fv_to_fe()

template<typename ScalarType , typename GridScalarType >
void terra::fv::hex::l2_project_fv_to_fe ( linalg::VectorQ1Scalar< ScalarType > &  dst,
const linalg::VectorFVScalar< ScalarType > &  src,
const grid::shell::DistributedDomain domain,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii,
std::vector< linalg::VectorQ1Scalar< ScalarType > > &  tmps 
)

L2 projection from a finite volume function into a Q1 (wedge) finite element function.

Use this function if you want to convert from the finite volume to the finite element space.

Given some finite volume function \(u_h^{FV}\) this function solves

\[ \int_\Omega u_h^{FE} v \ d\Omega = \int_\Omega u_h^{FV} v \ d\Omega, \quad \forall v \in V_h \]

for \(u_h^{FE}\). In matrix form:

\[ M u^{FE} = b \]

where \(M\) is the finite element mass matrix and

\[ b_i = \int_\Omega u^{FV}(x) \, \phi_i(x) \ dx = \sum_K u^{FV}_K \int_\Omega \phi_i(x) \ dx. \]

Internally, this function

  1. sets up \(b\) locally
  2. communicates \(b\) to complete the assembly
  3. solves the mass matrix system with CG
Parameters
dst[out] finite element function that is the solution of the mass matrix system
srcfinite volume vector to convert
domaindistributed domain (required for communication)
coords_shellcoords of the unit shell
coords_radiiradii of all shells
tmpsat least 5 tmp vectors (required for \(b\) and the CG solver)

◆ l2_project_fv_to_fe_lumped()

template<typename ScalarType , typename GridScalarType >
void terra::fv::hex::l2_project_fv_to_fe_lumped ( linalg::VectorQ1Scalar< ScalarType > &  dst,
const linalg::VectorFVScalar< ScalarType > &  src,
const grid::shell::DistributedDomain domain,
const grid::Grid3DDataVec< GridScalarType, 3 > &  coords_shell,
const grid::Grid2DDataScalar< GridScalarType > &  coords_radii,
std::vector< linalg::VectorQ1Scalar< ScalarType > > &  tmps 
)

Bound-preserving FV-to-FE transfer via a lumped mass matrix.

Computes, for every Q1 node \(i\),

\[ u_i^{FE} = \frac{b_i}{D_i}, \quad b_i = \sum_K u_K^{FV} \int_K \phi_i \, dx, \quad D_i = \sum_K \int_K \phi_i \, dx \]

where \(\phi_i\) are the Q1 wedge shape functions. Because every shape function is non-negative and the shape functions form a partition of unity, \(u_i^{FE}\) is a convex combination of the surrounding FV cell values. Consequently this projection is monotone: if \(u^{FV} \in [a, b]\) then \(u^{FE} \in [a, b]\), with no Gibbs-type over- or undershoots.

Compare with l2_project_fv_to_fe which uses a consistent mass matrix and solves a global linear system. The lumped variant is cheaper (no solver, no MPI reduction for convergence) but slightly less L2-accurate in smooth cases. It is the preferred choice when bound-preservation matters (e.g. compositional fields, density anomalies that feed directly into buoyancy).

Parameters
dst[out] FE scalar field receiving the projected values.
srcFV scalar field to project.
domainDistributed domain (required for MPI ghost communication).
coords_shellLateral node coordinates of the unit-sphere surface.
coords_radiiRadial shell radii.
tmpsAt least 2 temporary Q1 vectors (used for \(b\) and \(D\)).