Classes | |
| struct | BaseCorners |
| struct | BoundaryConditionMapping |
| class | DistributedDomain |
| Parallel data structure organizing the thick spherical shell metadata for distributed (MPI parallel) simulations. More... | |
| class | DomainInfo |
| Information about the thick spherical shell mesh. More... | |
| struct | SubdomainDistribution |
| class | SubdomainInfo |
| (Sortable) Globally unique identifier for a single subdomain of a diamond. More... | |
| class | SubdomainNeighborhood |
| Neighborhood information of a single subdomain. More... | |
Typedefs | |
| using | BoundaryConditions = BoundaryConditionMapping[2] |
| template<std::floating_point T> | |
| using | MemoizationCache = std::map< std::pair< int, int >, dense::Vec< T, 3 > > |
| using | SubdomainToRankDistributionFunction = std::function< mpi::MPIRank(const SubdomainInfo &, const int, const int) > |
Enumerations | |
| enum class | ShellBoundaryFlag : uint8_t { NO_FLAG = 0 , INNER = 1 << 0 , BOUNDARY = 1 << 1 , CMB = BOUNDARY | ( 1 << 2 ) , SURFACE = BOUNDARY | ( 1 << 3 ) , ALL = INNER | BOUNDARY | CMB | SURFACE } |
| FlagLike that indicates boundary types for the thick spherical shell. More... | |
| enum class | BoundaryConditionFlag : uint8_t { NEUMANN = 0 , DIRICHLET = 1 , FREESLIP = 2 } |
| FlagLike that indicates the type of boundary condition More... | |
Functions | |
| BoundaryConditionFlag | get_boundary_condition_flag (const BoundaryConditions bcs, ShellBoundaryFlag sbf) |
| Retrieve the boundary condition flag that is associated with a location in the shell e.g. SURFACE -> DIRICHLET. | |
| void | set_boundary_condition_flag (BoundaryConditions &bcs, ShellBoundaryFlag sbf, BoundaryConditionFlag bcf) |
| Set the boundary condition flag that is associated with a location in the shell e.g. SURFACE -> DIRICHLET. | |
| ShellBoundaryFlag | get_shell_boundary_flag (const BoundaryConditions bcs, BoundaryConditionFlag bcf) |
| Retrieve the ShellBoundary flag associated with a certain boundary condition type/flag. | |
| Grid4DDataScalar< ShellBoundaryFlag > | setup_boundary_mask_data (const DistributedDomain &domain) |
| Set up mask data for a distributed shell domain. The mask encodes boundary information for each grid node. | |
| template<std::floating_point T> | |
| std::vector< T > | uniform_shell_radii (T r_min, T r_max, int num_shells) |
| Computes the radial shell radii for a uniform grid. | |
| template<std::floating_point T> | |
| std::function< T(T) > | make_tanh_boundary_cluster (T k) |
| Map to be used in mapped_shell_radii. | |
| template<std::floating_point T> | |
| std::vector< T > | mapped_shell_radii (T r_min, T r_max, int num_shells, const std::function< T(T) > &map) |
| Computes the radial shell radii for a non-uniformly distributed grid. | |
| template<std::floating_point T> | |
| T | min_radial_h (const std::vector< T > &shell_radii) |
| Computes the min absolute distance of two entries in the passed vector of shell radii. | |
| template<std::floating_point T> | |
| T | max_radial_h (const std::vector< T > &shell_radii) |
| Computes the max absolute distance of two entries in the passed vector of shell radii. | |
| template<std::floating_point T> | |
| dense::Vec< T, 3 > | compute_node_recursive (int i, int j, const BaseCorners< T > &corners, MemoizationCache< T > &cache) |
| Computes the coordinates for a specific node (i, j) in the final refined grid. Uses recursion and memoization, sourcing base points from the BaseCorners struct. | |
| template<std::floating_point T> | |
| void | compute_subdomain (const typename Grid3DDataVec< T, 3 >::HostMirror &subdomain_coords_host, int subdomain_idx, const BaseCorners< T > &corners, int i_start_incl, int i_end_incl, int j_start_incl, int j_end_incl) |
| Generates coordinates for a rectangular subdomain of the refined spherical grid. | |
| template<std::floating_point T> | |
| void | unit_sphere_single_shell_subdomain_coords (const typename Grid3DDataVec< T, 3 >::HostMirror &subdomain_coords_host, int subdomain_idx, int diamond_id, int ntan, int i_start_incl, int i_end_incl, int j_start_incl, int j_end_incl) |
| template<std::floating_point T> | |
| void | unit_sphere_single_shell_subdomain_coords (const typename Grid3DDataVec< T, 3 >::HostMirror &subdomain_coords_host, int subdomain_idx, int diamond_id, int global_refinements, int num_subdomains_per_side, int subdomain_i, int subdomain_j) |
| std::ostream & | operator<< (std::ostream &os, const SubdomainInfo &si) |
| mpi::MPIRank | subdomain_to_rank_all_root (const SubdomainInfo &subdomain_info, const int num_subdomains_per_diamond_laterally, const int num_subdomains_per_diamond_radially) |
| Assigns all subdomains to root (rank 0). | |
| mpi::MPIRank | subdomain_to_rank_iterate_diamond_subdomains (const SubdomainInfo &subdomain_info, const int num_subdomains_per_diamond_laterally, const int num_subdomains_per_diamond_radially) |
| Distributes subdomains to ranks as evenly as possible. | |
| SubdomainDistribution | subdomain_distribution (const DistributedDomain &domain) |
| template<typename ValueType > | |
| Grid4DDataScalar< ValueType > | allocate_scalar_grid (const std::string label, const DistributedDomain &distributed_domain) |
| template<typename ValueType > | |
| Grid4DDataScalar< ValueType > | allocate_scalar_grid (const std::string label, const DistributedDomain &distributed_domain, const int level) |
| Kokkos::MDRangePolicy< Kokkos::Rank< 4 > > | local_domain_md_range_policy_nodes (const DistributedDomain &distributed_domain) |
| Kokkos::MDRangePolicy< Kokkos::Rank< 3 > > | local_domain_md_range_policy_cells_lateral (const DistributedDomain &distributed_domain) |
| Kokkos::MDRangePolicy< Kokkos::Rank< 4 > > | local_domain_md_range_policy_cells (const DistributedDomain &distributed_domain) |
| Kokkos::RangePolicy | local_domain_md_range_policy_cells_linearized (const DistributedDomain &distributed_domain) |
| template<std::floating_point T> | |
| Grid3DDataVec< T, 3 > | subdomain_unit_sphere_single_shell_coords (const DistributedDomain &domain) |
| Returns an initialized grid with the coordinates of all subdomains' nodes projected to the unit sphere. | |
| template<std::floating_point T> | |
| Grid2DDataScalar< T > | subdomain_shell_radii (const DistributedDomain &domain) |
| Returns an initialized grid with the radii of all subdomain nodes. | |
| Grid2DDataScalar< int > | subdomain_shell_idx (const DistributedDomain &domain) |
| Returns an initialized grid with the shell index of all subdomain nodes. | |
| template<typename CoordsShellType , typename CoordsRadiiType > | |
| dense::Vec< typename CoordsShellType::value_type, 3 > | coords (const int subdomain, const int x, const int y, const int r, const CoordsShellType &coords_shell, const CoordsRadiiType &coords_radii) |
| template<typename CoordsShellType , typename CoordsRadiiType > | |
| dense::Vec< typename CoordsShellType::value_type, 3 > | coords (const dense::Vec< int, 4 > subdomain_x_y_r, const CoordsShellType &coords_shell, const CoordsRadiiType &coords_radii) |
| using terra::grid::shell::BoundaryConditions = typedef BoundaryConditionMapping[2] |
| using terra::grid::shell::MemoizationCache = typedef std::map< std::pair< int, int >, dense::Vec< T, 3 > > |
| using terra::grid::shell::SubdomainToRankDistributionFunction = typedef std::function< mpi::MPIRank( const SubdomainInfo&, const int, const int ) > |
|
strong |
|
strong |
|
inline |
|
inline |
| dense::Vec< T, 3 > terra::grid::shell::compute_node_recursive | ( | int | i, |
| int | j, | ||
| const BaseCorners< T > & | corners, | ||
| MemoizationCache< T > & | cache | ||
| ) |
Computes the coordinates for a specific node (i, j) in the final refined grid. Uses recursion and memoization, sourcing base points from the BaseCorners struct.
| i | Row index (0 to corners.N). |
| j | Column index (0 to corners.N). |
| corners | Struct containing base corner coordinates and N = ntan - 1. |
| cache | Cache to store/retrieve already computed nodes. |
| void terra::grid::shell::compute_subdomain | ( | const typename Grid3DDataVec< T, 3 >::HostMirror & | subdomain_coords_host, |
| int | subdomain_idx, | ||
| const BaseCorners< T > & | corners, | ||
| int | i_start_incl, | ||
| int | i_end_incl, | ||
| int | j_start_incl, | ||
| int | j_end_incl | ||
| ) |
Generates coordinates for a rectangular subdomain of the refined spherical grid.
| subdomain_coords_host | a properly sized host-allocated view that is filled with the coordinates of the points |
| corners | Struct containing the base corner points and N = ntan - 1. |
| i_start_incl | Starting row index (inclusive) of the subdomain (global index). |
| i_end_incl | Ending row index (inclusive) of the subdomain (global index). |
| j_start_incl | Starting column index (inclusive) of the subdomain (global index). |
| j_end_incl | Ending column index (inclusive) of the subdomain (global index). |
| dense::Vec< typename CoordsShellType::value_type, 3 > terra::grid::shell::coords | ( | const dense::Vec< int, 4 > | subdomain_x_y_r, |
| const CoordsShellType & | coords_shell, | ||
| const CoordsRadiiType & | coords_radii | ||
| ) |
| dense::Vec< typename CoordsShellType::value_type, 3 > terra::grid::shell::coords | ( | const int | subdomain, |
| const int | x, | ||
| const int | y, | ||
| const int | r, | ||
| const CoordsShellType & | coords_shell, | ||
| const CoordsRadiiType & | coords_radii | ||
| ) |
| BoundaryConditionFlag terra::grid::shell::get_boundary_condition_flag | ( | const BoundaryConditions | bcs, |
| ShellBoundaryFlag | sbf | ||
| ) |
Retrieve the boundary condition flag that is associated with a location in the shell e.g. SURFACE -> DIRICHLET.
| ShellBoundaryFlag terra::grid::shell::get_shell_boundary_flag | ( | const BoundaryConditions | bcs, |
| BoundaryConditionFlag | bcf | ||
| ) |
Retrieve the ShellBoundary flag associated with a certain boundary condition type/flag.
|
inline |
|
inline |
|
inline |
|
inline |
| std::function< T(T) > terra::grid::shell::make_tanh_boundary_cluster | ( | T | k | ) |
Map to be used in mapped_shell_radii.
Can be used to have smaller elements near the shell boundaries.
Returns a function
\[ f(s) := \begin{cases} s & k \leq 0 \\ \frac{1}{2} \frac{\tanh( k( 2s-1 ) )}{\tanh(k) + 1} & k > 0 \end{cases} \]
| k | \(k = 0\) results in a uniform distribution (linear function), \(k > 0\) refines at the boundary (roughly: \(k \approx 1\): mild clustering, \(k \approx 2\): strong clustering) |
| std::vector< T > terra::grid::shell::mapped_shell_radii | ( | T | r_min, |
| T | r_max, | ||
| int | num_shells, | ||
| const std::function< T(T) > & | map | ||
| ) |
Computes the radial shell radii for a non-uniformly distributed grid.
Note that a shell is a 2D manifold in 3D space. A layer is a 3D volume in 3D space – it is sandwiched by two shells (one on each side).
The shell radii are generated from a uniform parameter ( \(N\) is the number of shells)
\[ s_i = \frac{i}{N-1}, \; i = 0,\dots,N-1, \]
which is mapped to a redistributed parameter
\[ t_i = f(s_i) \]
using a user-supplied function
\[ f : [0,1] \rightarrow [0,1]. \]
The physical radii are then computed as
\[ r_i = r_{\min} + (r_{\max} - r_{\min}) \, t_i . \]
The function \(f(s) = s\) therefore maps the shells uniformly. In areas where the gradient of \(f\) is small, shells are closer together. In areas where the gradient of \(f\) is large, shells are further apart.
Try make_tanh_boundary_cluster which returns a parameterized function that can be passed to this function.
| T | Floating-point type. |
| r_min | Radius of the innermost shell. |
| r_max | Radius of the outermost shell. |
| num_shells | Number of shells \( N \). |
| map | Mapping function \( f : [0,1] \rightarrow [0,1] \) controlling shell distribution. It should be monotone increasing and satisfy \( f(0)=0 \) and \( f(1)=1 \). |
| T terra::grid::shell::max_radial_h | ( | const std::vector< T > & | shell_radii | ) |
Computes the max absolute distance of two entries in the passed vector of shell radii.
| T terra::grid::shell::min_radial_h | ( | const std::vector< T > & | shell_radii | ) |
Computes the min absolute distance of two entries in the passed vector of shell radii.
|
inline |
| void terra::grid::shell::set_boundary_condition_flag | ( | BoundaryConditions & | bcs, |
| ShellBoundaryFlag | sbf, | ||
| BoundaryConditionFlag | bcf | ||
| ) |
Set the boundary condition flag that is associated with a location in the shell e.g. SURFACE -> DIRICHLET.
|
inline |
Set up mask data for a distributed shell domain. The mask encodes boundary information for each grid node.
| domain | Distributed shell domain. |
|
inline |
|
inline |
Returns an initialized grid with the shell index of all subdomain nodes.
The layout is
grid( local_subdomain_id, r_idx ) = global_shell_idx
| Grid2DDataScalar< T > terra::grid::shell::subdomain_shell_radii | ( | const DistributedDomain & | domain | ) |
Returns an initialized grid with the radii of all subdomain nodes.
The layout is
grid( local_subdomain_id, r_idx ) = radius
|
inline |
Assigns all subdomains to root (rank 0).
|
inline |
Distributes subdomains to ranks as evenly as possible.
Not really sophisticated, just iterates over subdomain indices. Tries to concentrate ranks in as few diamonds as possible. But does not optimize surface to volume ratio. For that we need something like a z- or Hilbert curve.
| Grid3DDataVec< T, 3 > terra::grid::shell::subdomain_unit_sphere_single_shell_coords | ( | const DistributedDomain & | domain | ) |
Returns an initialized grid with the coordinates of all subdomains' nodes projected to the unit sphere.
The layout is
grid( local_subdomain_id, x_idx, y_idx, node_coord )
where node_coord is in {0, 1, 2} and refers to the cartesian coordinate of the point.
| std::vector< T > terra::grid::shell::uniform_shell_radii | ( | T | r_min, |
| T | r_max, | ||
| int | num_shells | ||
| ) |
Computes the radial shell radii for a uniform grid.
Note that a shell is a 2D manifold in 3D space. A layer is a 3D volume in 3D space - it is sandwiched by two shells (one on each side).
| r_min | Radius of the innermost shell. |
| r_max | Radius of the outermost shell. |
| num_shells | Number of shells. |
| void terra::grid::shell::unit_sphere_single_shell_subdomain_coords | ( | const typename Grid3DDataVec< T, 3 >::HostMirror & | subdomain_coords_host, |
| int | subdomain_idx, | ||
| int | diamond_id, | ||
| int | global_refinements, | ||
| int | num_subdomains_per_side, | ||
| int | subdomain_i, | ||
| int | subdomain_j | ||
| ) |
| void terra::grid::shell::unit_sphere_single_shell_subdomain_coords | ( | const typename Grid3DDataVec< T, 3 >::HostMirror & | subdomain_coords_host, |
| int | subdomain_idx, | ||
| int | diamond_id, | ||
| int | ntan, | ||
| int | i_start_incl, | ||
| int | i_end_incl, | ||
| int | j_start_incl, | ||
| int | j_end_incl | ||
| ) |