Ghost layer communication for scalar finite volume (FV) fields on the shell grid. More...
#include <algorithm>#include <tuple>#include <vector>#include <mpi.h>#include "grid/grid_types.hpp"#include "grid/shell/spherical_shell.hpp"#include "mpi/mpi.hpp"Go to the source code of this file.
Classes | |
| class | terra::communication::shell::FVGhostLayerBuffers< ScalarType > |
| Pre-allocated send/recv buffers and sorted face-pair lists for FV ghost layer comm. More... | |
| struct | terra::communication::shell::FVGhostLayerBuffers< ScalarType >::FacePair |
| All metadata needed to drive communication for one subdomain face. More... | |
| class | terra::communication::shell::FVGhostLayerVecBuffers< ScalarType, VecDim > |
| Pre-allocated send/recv buffers for vector-valued FV ghost layer communication. More... | |
| struct | terra::communication::shell::FVGhostLayerVecBuffers< ScalarType, VecDim >::FacePair |
Namespaces | |
| namespace | terra |
| namespace | terra::communication |
| namespace | terra::communication::shell |
| namespace | terra::communication::shell::fv_detail |
Functions | |
| template<typename ScalarType > | |
| void | terra::communication::shell::fv_detail::pack_inner_cells (const grid::Grid2DDataScalar< ScalarType > &buffer, const grid::Grid4DDataScalar< ScalarType > &data, const int local_subdomain_id, const grid::BoundaryFace face) |
Packs the innermost real-cell layer adjacent to face into buffer. | |
| template<typename ScalarType > | |
| void | terra::communication::shell::fv_detail::unpack_to_ghost (const grid::Grid2DDataScalar< ScalarType > &buffer, const grid::Grid4DDataScalar< ScalarType > &data, const int local_subdomain_id, const grid::BoundaryFace face, const grid::BoundaryDirection dir0, const grid::BoundaryDirection dir1) |
Writes buffer into the ghost cell layer at face. | |
| template<typename ScalarType > | |
| void | terra::communication::shell::update_fv_ghost_layers (const grid::shell::DistributedDomain &domain, const grid::Grid4DDataScalar< ScalarType > &data, FVGhostLayerBuffers< ScalarType > &bufs) |
| Fills all ghost layers of a scalar FV field from neighbouring subdomains. | |
| template<typename ScalarType > | |
| void | terra::communication::shell::update_fv_ghost_layers (const grid::shell::DistributedDomain &domain, const grid::Grid4DDataScalar< ScalarType > &data) |
| Convenience overload — allocates temporary buffers internally. | |
| template<typename ScalarType , int VecDim> | |
| void | terra::communication::shell::fv_detail::pack_inner_cells_vec (const grid::Grid3DDataScalar< ScalarType > &buffer, const grid::Grid4DDataVec< ScalarType, VecDim > &data, const int local_subdomain_id, const grid::BoundaryFace face) |
Packs the innermost real-cell layer adjacent to face into buffer for a vector field. | |
| template<typename ScalarType , int VecDim> | |
| void | terra::communication::shell::fv_detail::unpack_to_ghost_vec (const grid::Grid3DDataScalar< ScalarType > &buffer, const grid::Grid4DDataVec< ScalarType, VecDim > &data, const int local_subdomain_id, const grid::BoundaryFace face, const grid::BoundaryDirection dir0, const grid::BoundaryDirection dir1) |
Writes buffer into the ghost cell layer at face for a vector field. | |
| template<typename ScalarType , int VecDim> | |
| void | terra::communication::shell::update_fv_ghost_layers (const grid::shell::DistributedDomain &domain, const grid::Grid4DDataVec< ScalarType, VecDim > &data, FVGhostLayerVecBuffers< ScalarType, VecDim > &bufs) |
| Fills all ghost layers of a vector-valued FV field from neighbouring subdomains. | |
| template<typename ScalarType , int VecDim> | |
| void | terra::communication::shell::update_fv_ghost_layers (const grid::shell::DistributedDomain &domain, const grid::Grid4DDataVec< ScalarType, VecDim > &data) |
| Convenience overload for vector fields — allocates temporary buffers internally. | |
Variables | |
| constexpr int | terra::communication::shell::MPI_TAG_FV_GHOST_LAYERS = 200 |
| MPI tag for FV ghost layer messages (distinct from MPI_TAG_BOUNDARY_DATA = 100). | |
| constexpr int | terra::communication::shell::MPI_TAG_FV_VEC_GHOST_LAYERS = 201 |
| MPI tag for vector-valued FV ghost layer messages. | |
Ghost layer communication for scalar finite volume (FV) fields on the shell grid.
FV cells are owned exclusively by one subdomain. Each cell array has a one-cell-wide ghost layer on every side so that operator stencils can read neighbour values without extra index checks. Unlike the FE communication in communication.hpp — which adds contributions from shared nodes — here we simply copy the innermost real cells of the neighbour into our ghost cells (no reduction required).
Layout reminder (Grid4DDataScalar, shape [n_subdomains, N_lat+1, N_lat+1, N_rad+1]):
Only face-to-face communication is implemented here. A 6-neighbour cell stencil (±x, ±y, ±r) never touches edge or vertex ghost cells, so face communication is sufficient.
Typical usage (operator that is applied repeatedly):