3#include "../../quadrature/quadrature.hpp"
16template <
typename ScalarT,
int VecDim = 3 >
58 , treat_boundary_( treat_boundary )
60 , operator_apply_mode_( operator_apply_mode )
61 , operator_communication_mode_( operator_communication_mode )
63 , send_buffers_( domain )
64 , recv_buffers_( domain )
73 operator_apply_mode_ = operator_apply_mode;
74 operator_communication_mode_ = operator_communication_mode;
89 if ( src_.extent( 0 ) != dst_.extent( 0 ) || src_.extent( 1 ) != dst_.extent( 1 ) ||
90 src_.extent( 2 ) != dst_.extent( 2 ) || src_.extent( 3 ) != dst_.extent( 3 ) )
92 throw std::runtime_error(
"VectorLaplace: src/dst mismatch" );
95 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
97 throw std::runtime_error(
"VectorLaplace: src/dst mismatch" );
100 util::Timer timer_kernel(
"vector_laplace_kernel" );
110 domain_, dst_, send_buffers_, recv_buffers_ );
115 KOKKOS_INLINE_FUNCTION
void
116 operator()(
const int local_subdomain_id,
const int x_cell,
const int y_cell,
const int r_cell )
const
123 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
124 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
130 ScalarT quad_weights[num_quad_points];
140 for (
int dimi = 0; dimi < 3; ++dimi )
142 for (
int dimj = 0; dimj < 3; ++dimj )
144 if ( diagonal_ and dimi != dimj )
150 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
151 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
152 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
153 for (
int i = 0; i < 8; i++ )
155 src_local_hex[i] = src_(
157 x_cell + hex_offset_x[i],
158 y_cell + hex_offset_y[i],
159 r_cell + hex_offset_r[i],
164 for (
int q = 0; q < num_quad_points; q++ )
166 const auto quad_weight = quad_weights[q];
171 jac( wedge_phy_surf[wedge], r_1, r_2, quad_points[q] );
172 const auto det = J.det();
173 const auto abs_det = Kokkos::abs( det );
178 k_eval +=
shape( k, quad_points[q] ) * k_local_hex[wedge]( k );
188 sym_grad_i[k] = ( grad_i + grad_i.
transposed() ) * 0.5;
193 sym_grad_j[k] = ( grad_j + grad_j.
transposed() ) * 0.5;
210 else if ( treat_boundary_ && r_cell == 0 )
225 else if ( treat_boundary_ && r_cell + 1 == radii_.extent( 1 ) - 1 )
257 for (
int i = 0; i < 8; i++ )
259 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
260 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
261 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
266 x_cell + hex_offset_x[i],
267 y_cell + hex_offset_y[i],
268 r_cell + hex_offset_r[i],
286 const int dimj )
const
288 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
289 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
290 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
298 grad_u = grad_u + sym_grad_i[j] *
299 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]];
305 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
306 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u ) * abs_det;
320 const int dimj )
const
322 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
323 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
324 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
331 grad_u = grad_u + sym_grad_i[j] *
332 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]];
338 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
339 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u ) * abs_det;
345 for (
int i = 0; i < 3; i++ )
347 const auto grad_u_diag =
348 sym_grad_i[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
350 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
351 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u_diag ) * abs_det;
366 const int dimj )
const
368 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
369 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
370 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
376 for (
int j = 0; j < 3; j++ )
378 grad_u = grad_u + sym_grad_i[j] *
379 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]];
383 for (
int i = 0; i < 3; i++ )
385 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
386 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u ) * abs_det;
394 const auto grad_u_diag =
395 sym_grad_i[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
397 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
398 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u_diag ) * abs_det;
413 const int dimj )
const
415 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
416 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
417 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
422 const auto grad_u_diag =
423 sym_grad_i[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
425 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
426 2 * quad_weight * k_eval * ( sym_grad_j[i] ).double_contract( grad_u_diag ) * abs_det;
Send and receive buffers for all process-local subdomain boundaries.
Definition communication.hpp:56
Definition epsilon.hpp:18
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition epsilon.hpp:77
void diagonal(ScalarType src_local_hex[8], ScalarType dst_local_hex[8], const ScalarType k_eval, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, dense::Mat< ScalarType, 3, 3 > *sym_grad_i, dense::Mat< ScalarType, 3, 3 > *sym_grad_j, const int dimi, const int dimj) const
Definition epsilon.hpp:403
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition epsilon.hpp:116
void set_operator_apply_and_communication_modes(const linalg::OperatorApplyMode operator_apply_mode, const linalg::OperatorCommunicationMode operator_communication_mode)
Definition epsilon.hpp:69
void dirichlet_cmb(ScalarType src_local_hex[8], ScalarType dst_local_hex[8], const ScalarType k_eval, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, dense::Mat< ScalarType, 3, 3 > *sym_grad_i, dense::Mat< ScalarType, 3, 3 > *sym_grad_j, const int dimi, const int dimj) const
Definition epsilon.hpp:310
const grid::Grid4DDataScalar< ScalarType > & k_grid_data()
Definition epsilon.hpp:67
ScalarT ScalarType
Definition epsilon.hpp:22
void neumann(ScalarType src_local_hex[8], ScalarType dst_local_hex[8], const ScalarType k_eval, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, dense::Mat< ScalarType, 3, 3 > *sym_grad_i, dense::Mat< ScalarType, 3, 3 > *sym_grad_j, const int dimi, const int dimj) const
Definition epsilon.hpp:276
Epsilon(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const grid::Grid4DDataScalar< ScalarT > &k, bool treat_boundary, bool diagonal, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition epsilon.hpp:44
void dirichlet_surface(ScalarType src_local_hex[8], ScalarType dst_local_hex[8], const ScalarType k_eval, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, dense::Mat< ScalarType, 3, 3 > *sym_grad_i, dense::Mat< ScalarType, 3, 3 > *sym_grad_j, const int dimi, const int dimj) const
Definition epsilon.hpp:356
Parallel data structure organizing the thick spherical shell metadata for distributed (MPI parallel) ...
Definition spherical_shell.hpp:2498
Static assertion: VectorQ1Scalar satisfies VectorLike concept.
Definition vector_q1.hpp:162
const grid::Grid4DDataVec< ScalarType, VecDim > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:280
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
void extract_local_wedge_scalar_coefficients(dense::Vec< T, 6 >(&local_coefficients)[2], const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const grid::Grid4DDataScalar< T > &global_coefficients)
Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient v...
Definition kernel_helpers.hpp:306
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 T shape(const int node_idx, const T xi, const T eta, const T zeta)
(Tensor-product) Shape function.
Definition integrands.hpp:146
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 ****[VecDim], Layout > Grid4DDataVec
Definition grid_types.hpp:43
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 Mat< T, Cols, Rows > transposed() const
Definition mat.hpp:187
static constexpr Mat from_single_col_vec(const Vec< T, Cols > &col, const int d)
Definition mat.hpp:66