4#include "../../quadrature/quadrature.hpp"
17template <
typename ScalarT >
59 , treat_boundary_( treat_boundary )
61 , operator_apply_mode_( operator_apply_mode )
62 , operator_communication_mode_( operator_communication_mode )
64 , send_buffers_( domain )
65 , recv_buffers_( domain )
80 if ( src_.extent( 0 ) != dst_.extent( 0 ) || src_.extent( 1 ) != dst_.extent( 1 ) ||
81 src_.extent( 2 ) != dst_.extent( 2 ) || src_.extent( 3 ) != dst_.extent( 3 ) )
83 throw std::runtime_error(
"LaplaceSimple: src/dst mismatch" );
86 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
88 throw std::runtime_error(
"LaplaceSimple: src/dst mismatch" );
100 domain_, dst_, send_buffers_, recv_buffers_ );
105 KOKKOS_INLINE_FUNCTION
void
106 operator()(
const int local_subdomain_id,
const int x_cell,
const int y_cell,
const int r_cell )
const
113 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
114 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
120 ScalarT quad_weights[num_quad_points];
130 for (
int i = 0; i < 8; i++ )
132 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
133 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
134 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
136 src_local_hex[i] = src_(
137 local_subdomain_id, x_cell + hex_offset_x[i], y_cell + hex_offset_y[i], r_cell + hex_offset_r[i] );
140 const bool at_bot_boundary =
147 for (
int q = 0; q < num_quad_points; q++ )
149 const auto quad_point = quad_points[q];
150 const auto quad_weight = quad_weights[q];
154 const auto J =
jac( wedge_phy_surf[wedge], r_1, r_2, quad_points[q] );
155 const auto det = J.det();
156 const auto abs_det = Kokkos::abs( det );
157 const auto J_inv_transposed = J.inv_transposed( det );
163 grad_phy[k] = J_inv_transposed *
grad_shape( k, quad_point );
168 diagonal( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
170 else if ( treat_boundary_ && at_bot_boundary )
173 dirichlet_bot( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
175 else if ( treat_boundary_ && at_top_boundary )
178 dirichlet_top( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
182 neumann( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
187 for (
int i = 0; i < 8; i++ )
189 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
190 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
191 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
195 local_subdomain_id, x_cell + hex_offset_x[i], y_cell + hex_offset_y[i], r_cell + hex_offset_r[i] ),
208 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
209 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
210 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
218 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
224 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
225 quad_weight * grad_phy[i].
dot( grad_u ) * abs_det;
237 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
238 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
239 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
247 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
253 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
254 quad_weight * grad_phy[i].
dot( grad_u ) * abs_det;
258 for (
int i = 0; i < 3; i++ )
260 const auto grad_u_diag =
261 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
263 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
264 quad_weight * grad_phy[i].
dot( grad_u_diag ) * abs_det;
276 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
277 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
278 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
283 for (
int j = 0; j < 3; j++ )
286 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
290 for (
int i = 0; i < 3; i++ )
292 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
293 quad_weight * grad_phy[i].
dot( grad_u ) * abs_det;
299 const auto grad_u_diag =
300 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
302 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
303 quad_weight * grad_phy[i].
dot( grad_u_diag ) * abs_det;
315 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
316 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
317 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
324 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
326 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
327 quad_weight * grad_phy[i].
dot( grad_u ) * abs_det;
Definition laplace.hpp:19
void neumann(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:200
void dirichlet_bot(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:229
ScalarT ScalarType
Definition laplace.hpp:23
void diagonal(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:307
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition laplace.hpp:68
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition laplace.hpp:106
Laplace(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &mask, bool treat_boundary, bool diagonal, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition laplace.hpp:45
void dirichlet_top(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:268
Parallel data structure organizing the thick spherical shell metadata for distributed (MPI parallel) ...
Definition spherical_shell.hpp:2498
Q1 scalar finite element vector on a distributed shell grid.
Definition vector_q1.hpp:21
const grid::Grid4DDataScalar< ScalarType > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:137
Timer supporting RAII scope or manual stop.
Definition timer.hpp:270
void stop()
Stop the timer and record elapsed time.
Definition timer.hpp:289
Concept for types that behave like linear operators.
Definition operator.hpp:57
void unpack_and_reduce_local_subdomain_boundaries(const grid::shell::DistributedDomain &domain, const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_recv_buffers, CommunicationReduction reduction=CommunicationReduction::SUM)
Unpacks and reduces local subdomain boundaries.
Definition communication.hpp:672
void pack_send_and_recv_local_subdomain_boundaries(const grid::shell::DistributedDomain &domain, const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_send_buffers, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_recv_buffers)
Packs, sends and recvs local subdomain boundaries using two sets of buffers.
Definition communication.hpp:242
Definition boundary_mass.hpp:14
constexpr void quad_felippa_3x2_quad_weights(T(&quad_weights)[quad_felippa_3x2_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:93
constexpr int quad_felippa_3x2_num_quad_points
Definition wedge/quadrature/quadrature.hpp:66
constexpr void quad_felippa_3x2_quad_points(dense::Vec< T, 3 >(&quad_points)[quad_felippa_3x2_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:70
constexpr int num_nodes_per_wedge_surface
Definition kernel_helpers.hpp:6
void wedge_surface_physical_coords(dense::Vec< T, 3 >(&wedge_surf_phy_coords)[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const grid::Grid3DDataVec< T, 3 > &lateral_grid, const int local_subdomain_id, const int x_cell, const int y_cell)
Extracts the (unit sphere) surface vertex coords of the two wedges of a hex cell.
Definition kernel_helpers.hpp:26
constexpr int num_wedges_per_hex_cell
Definition kernel_helpers.hpp:5
constexpr int num_nodes_per_wedge
Definition kernel_helpers.hpp:7
constexpr dense::Vec< T, 3 > grad_shape(const int node_idx, const T xi, const T eta, const T zeta)
Gradient of the full shape function:
Definition integrands.hpp:228
constexpr dense::Mat< T, 3, 3 > jac(const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &p3_phy, const T r_1, const T r_2, const T xi, const T eta, const T zeta)
Definition integrands.hpp:643
Kokkos::MDRangePolicy< Kokkos::Rank< 4 > > local_domain_md_range_policy_cells(const DistributedDomain &distributed_domain)
Definition spherical_shell.hpp:2668
Kokkos::View< ScalarType ***[VecDim], Layout > Grid3DDataVec
Definition grid_types.hpp:40
Kokkos::View< ScalarType ****, Layout > Grid4DDataScalar
Definition grid_types.hpp:25
Kokkos::View< ScalarType **, Layout > Grid2DDataScalar
Definition grid_types.hpp:19
OperatorApplyMode
Modes for applying an operator to a vector.
Definition operator.hpp:30
@ Replace
Overwrite the destination vector.
OperatorCommunicationMode
Modes for communication during operator application.
Definition operator.hpp:40
@ CommunicateAdditively
Communicate and add results.
constexpr bool has_flag(E mask_value, E flag) noexcept
Checks if a bitmask value contains a specific flag.
Definition bit_masking.hpp:43
T dot(const Vec &other) const
Definition vec.hpp:39