4#include "../../quadrature/quadrature.hpp"
18template <
typename ScalarT >
31 bool single_quadpoint_ =
false;
76 , treat_boundary_( treat_boundary )
77 , diagonal_( diagonal )
78 , operator_apply_mode_( operator_apply_mode )
79 , operator_communication_mode_( operator_communication_mode )
80 , operator_stored_matrix_mode_( operator_stored_matrix_mode )
82 , send_buffers_( domain )
83 , recv_buffers_( domain )
92 KOKKOS_INLINE_FUNCTION
94 const int local_subdomain_id,
100 return util::has_flag( mask_( local_subdomain_id, x_cell, y_cell, r_cell ), flag );
107 operator_apply_mode_ = operator_apply_mode;
108 operator_communication_mode_ = operator_communication_mode;
134 operator_stored_matrix_mode_ = operator_stored_matrix_mode;
140 domain_, operator_stored_matrix_mode_, level_range, GCAElements );
147 KOKKOS_INLINE_FUNCTION
149 const int local_subdomain_id,
158 local_matrix_storage_.
set_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge, mat );
164 KOKKOS_INLINE_FUNCTION
166 const int local_subdomain_id,
170 const int wedge )
const
175 if ( !local_matrix_storage_.
has_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge ) )
177 Kokkos::abort(
"No matrix found at that spatial index." );
179 return local_matrix_storage_.
get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge );
197 if ( src_.extent( 0 ) != dst_.extent( 0 ) || src_.extent( 1 ) != dst_.extent( 1 ) ||
198 src_.extent( 2 ) != dst_.extent( 2 ) || src_.extent( 3 ) != dst_.extent( 3 ) )
200 throw std::runtime_error(
"LaplaceSimple: src/dst mismatch" );
203 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
205 throw std::runtime_error(
"LaplaceSimple: src/dst mismatch" );
213 domain_, dst_, send_buffers_, recv_buffers_ );
218 KOKKOS_INLINE_FUNCTION
void
219 operator()(
const int local_subdomain_id,
const int x_cell,
const int y_cell,
const int r_cell )
const
221 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
222 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
223 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
232 A[0] = local_matrix_storage_.
get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 0 );
233 A[1] = local_matrix_storage_.
get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 1 );
237 if ( local_matrix_storage_.
has_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 0 ) &&
238 local_matrix_storage_.
has_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 1 ) )
240 A[0] = local_matrix_storage_.
get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 0 );
241 A[1] = local_matrix_storage_.
get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 1 );
257 if ( treat_boundary_ )
262 boundary_mask.
fill( 1.0 );
266 for (
int i = 0; i < 6; i++ )
268 for (
int j = 0; j < 6; j++ )
270 if ( i != j && ( i < 3 || j < 3 ) )
272 boundary_mask( i, j ) = 0.0;
278 if ( r_cell + 1 == radii_.extent( 1 ) - 1 )
281 for (
int i = 0; i < 6; i++ )
283 for (
int j = 0; j < 6; j++ )
285 if ( i != j && ( i >= 3 || j >= 3 ) )
287 boundary_mask( i, j ) = 0.0;
302 dst[0] = A[0] * src[0];
303 dst[1] = A[1] * src[1];
316 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
317 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
329 for (
int i = 0; i < 8; i++ )
331 src_local_hex[i] = src_(
332 local_subdomain_id, x_cell + hex_offset_x[i], y_cell + hex_offset_y[i], r_cell + hex_offset_r[i] );
337 for (
int q = 0; q < num_quad_points; q++ )
339 const auto w = single_quadpoint_ ? quad_weights_1x1_[q] : quad_weights_3x2_[q];
340 const auto qp = single_quadpoint_ ? quad_points_1x1_[q] : quad_points_3x2_[q];
348 wedge, qp, w, r_1, r_2, wedge_phy_surf, k_local_hex, grad, jdet_keval_quadweight );
355 k_eval +=
shape( k, qp ) * k_local_hex[wedge]( k );
358 jdet_keval_quadweight *= k_eval;
360 fused_local_mv( src_local_hex, dst_local_hex, wedge, jdet_keval_quadweight, grad, r_cell );
364 for (
int i = 0; i < 8; i++ )
369 x_cell + hex_offset_x[i],
370 y_cell + hex_offset_y[i],
371 r_cell + hex_offset_r[i] ),
389 const auto det = J.det();
390 const auto abs_det = Kokkos::abs( det );
395 grad[k] = J_inv_transposed *
grad_shape( k, quad_point );
397 jdet_quadweight = quad_weight * abs_det;
402 KOKKOS_INLINE_FUNCTION
404 const int local_subdomain_id,
408 const int wedge )
const
416 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
417 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
428 for (
int q = 0; q < num_quad_points; q++ )
430 const auto w = single_quadpoint_ ? quad_weights_1x1_[q] : quad_weights_3x2_[q];
431 const auto qp = single_quadpoint_ ? quad_points_1x1_[q] : quad_points_3x2_[q];
436 wedge, qp, w, r_1, r_2, wedge_phy_surf, k_local_hex, grad, jdet_keval_quadweight );
443 k_eval +=
shape( k, qp ) * k_local_hex[wedge]( k );
446 jdet_keval_quadweight *= k_eval;
453 A( i, j ) += jdet_keval_quadweight * grad[j].
dot( grad[i] );
459 if ( treat_boundary_ )
462 boundary_mask.
fill( 1.0 );
467 for (
int i = 0; i < 6; i++ )
469 for (
int j = 0; j < 6; j++ )
472 if ( i != j && ( i < 3 || j < 3 ) )
474 boundary_mask( i, j ) = 0.0;
480 if ( r_cell + 1 == radii_.extent( 1 ) - 1 )
483 for (
int i = 0; i < 6; i++ )
485 for (
int j = 0; j < 6; j++ )
488 if ( i != j && ( i >= 3 || j >= 3 ) )
490 boundary_mask( i, j ) = 0.0;
510 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
511 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
512 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
517 const bool at_cmb = r_cell == 0;
518 const bool at_surface = r_cell + 1 == radii_.extent( 1 ) - 1;
520 int surface_shift = 0;
525 if ( treat_boundary_ && at_cmb )
530 else if ( treat_boundary_ && at_surface )
541 grad[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
550 dst_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] +=
551 jdet_keval_quadweight * ( grad[j] ).dot( grad_u );
557 if ( diagonal_ || ( ( treat_boundary_ && ( at_cmb || at_surface ) ) ) )
562 const auto grad_u_diag =
563 grad[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
565 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
566 jdet_keval_quadweight * ( grad[i] ).dot( grad_u_diag );
Definition div_k_grad.hpp:20
void set_operator_apply_and_communication_modes(const linalg::OperatorApplyMode operator_apply_mode, const linalg::OperatorCommunicationMode operator_communication_mode)
Definition div_k_grad.hpp:103
grid::Grid3DDataVec< ScalarT, 3 > get_grid() const
Getter for grid member.
Definition div_k_grad.hpp:124
static constexpr int LocalMatrixDim
Definition div_k_grad.hpp:25
void set_stored_matrix_mode(linalg::OperatorStoredMatrixMode operator_stored_matrix_mode, int level_range, grid::Grid4DDataScalar< ScalarType > GCAElements)
Definition div_k_grad.hpp:129
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition div_k_grad.hpp:219
dense::Mat< ScalarT, LocalMatrixDim, LocalMatrixDim > assemble_local_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge) const
assemble the local matrix and return it for a given element, wedge, and vectorial component (determin...
Definition div_k_grad.hpp:403
const grid::Grid4DDataScalar< ScalarType > & k_grid_data()
Getter for coefficient.
Definition div_k_grad.hpp:115
grid::Grid2DDataScalar< ScalarT > get_radii() const
Getter for radii member.
Definition div_k_grad.hpp:121
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition div_k_grad.hpp:187
void set_local_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge, const dense::Mat< ScalarT, LocalMatrixDim, LocalMatrixDim > &mat) const
Set the local matrix stored in the operator.
Definition div_k_grad.hpp:148
void fused_local_mv(ScalarType src_local_hex[8], ScalarType dst_local_hex[8], const int wedge, const ScalarType jdet_keval_quadweight, dense::Vec< ScalarType, 3 > *grad, int r_cell) const
Definition div_k_grad.hpp:502
void assemble_trial_test_vecs(const int wedge, const dense::Vec< ScalarType, 3 > &quad_point, const ScalarType quad_weight, const ScalarT r_1, const ScalarT r_2, dense::Vec< ScalarT, 3 >(*wedge_phy_surf)[3], const dense::Vec< ScalarT, 6 > *k_local_hex, dense::Vec< ScalarType, 3 > *grad, ScalarType &jdet_quadweight) const
Definition div_k_grad.hpp:377
linalg::OperatorStoredMatrixMode get_stored_matrix_mode()
Definition div_k_grad.hpp:144
void set_single_quadpoint(bool v)
S/Getter for quadpoint member.
Definition div_k_grad.hpp:127
ScalarT ScalarType
Definition div_k_grad.hpp:24
void set_diagonal(bool v)
S/Getter for diagonal member.
Definition div_k_grad.hpp:112
dense::Mat< ScalarT, LocalMatrixDim, LocalMatrixDim > get_local_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge) const
Retrives the local matrix if there is stored local matrices, the desired local matrix is loaded and r...
Definition div_k_grad.hpp:165
const grid::shell::DistributedDomain & get_domain() const
Getter for domain member.
Definition div_k_grad.hpp:118
bool has_flag(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, grid::shell::ShellBoundaryFlag flag) const
Getter for mask member.
Definition div_k_grad.hpp:93
DivKGrad(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &mask, const grid::Grid4DDataScalar< ScalarType > &k, bool treat_boundary, bool diagonal, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively, linalg::OperatorStoredMatrixMode operator_stored_matrix_mode=linalg::OperatorStoredMatrixMode::Off)
Definition div_k_grad.hpp:59
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
bool has_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge) const
Checks for presence of a local matrix for a certain element.
Definition local_matrix_storage.hpp:223
dense::Mat< ScalarT, LocalMatrixDim, LocalMatrixDim > get_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge) const
Retrives the local matrix if there is stored local matrices, the desired local matrix is loaded and r...
Definition local_matrix_storage.hpp:175
void set_matrix(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int wedge, dense::Mat< ScalarT, LocalMatrixDim, LocalMatrixDim > mat) const
Set the local matrix stored in the operator.
Definition local_matrix_storage.hpp:118
Concept for types that can be used as Galerkin coarse-grid operators in a multigrid hierarchy....
Definition operator.hpp:81
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 void quad_felippa_1x1_quad_points(dense::Vec< T, 3 >(&quad_points)[quad_felippa_1x1_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:36
constexpr void quad_felippa_1x1_quad_weights(T(&quad_weights)[quad_felippa_1x1_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:43
constexpr int quad_felippa_1x1_num_quad_points
Definition wedge/quadrature/quadrature.hpp:32
constexpr int num_nodes_per_wedge_surface
Definition kernel_helpers.hpp:6
void atomically_add_local_wedge_scalar_coefficients(const grid::Grid4DDataScalar< T > &global_coefficients, const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const dense::Vec< T, 6 >(&local_coefficients)[2])
Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coeff...
Definition kernel_helpers.hpp:407
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
ShellBoundaryFlag
FlagLike that indicates boundary types for the thick spherical shell.
Definition shell/bit_masks.hpp:12
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.
OperatorStoredMatrixMode
Modes for applying stored matrices.
Definition operator.hpp:47
@ Selective
Use stored matrices on selected, marked elements only, assemble on all others.
@ Full
Use stored matrices on all elements.
@ Off
Do not use stored matrices.
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
void fill(const T value)
Definition mat.hpp:201
constexpr Mat diagonal() const
Definition mat.hpp:377
Mat & hadamard_product(const Mat &mat)
Definition mat.hpp:213
T dot(const Vec &other) const
Definition vec.hpp:39