16template <
typename ScalarT >
31 bool lumped_diagonal_;
47 const bool diagonal =
false,
48 const bool lumped_diagonal =
false,
55 , diagonal_( diagonal )
56 , lumped_diagonal_( lumped_diagonal )
57 , operator_apply_mode_( operator_apply_mode )
58 , operator_communication_mode_( operator_communication_mode )
60 , send_buffers_( domain )
61 , recv_buffers_( domain )
79 Kokkos::abort(
"Operator apply mode not implemented." );
91 std::vector< std::unique_ptr< std::array< int, 11 > > > expected_recvs_metadata;
92 std::vector< std::unique_ptr< MPI_Request > > expected_recvs_requests;
95 domain_, dst_, send_buffers_, recv_buffers_ );
101 Kokkos::abort(
"Communication mode not implemented." );
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
119 ScalarT quad_weights[num_quad_points];
129 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
130 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
139 for (
int q = 0; q < num_quad_points; q++ )
147 const ScalarT shape_i =
shape_lat( i, quad_points[q] ) *
shape_rad( i, quad_points[q] );
148 const ScalarT shape_j =
shape_lat( j, quad_points[q] ) *
shape_rad( j, quad_points[q] );
151 quad_weights[q] * ( shape_i * shape_j * r * r * grad_r * det_jac_lat[wedge][q] );
162 else if ( lumped_diagonal_ )
175 dst[0] = A[0] * src[0];
176 dst[1] = A[1] * src[1];
ScalarT ScalarType
Definition mass.hpp:22
void set_lumped_diagonal(bool v)
S/Getter for lumped diagonal member.
Definition mass.hpp:68
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition mass.hpp:70
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition mass.hpp:106
void set_diagonal(bool v)
S/Getter for diagonal member.
Definition mass.hpp:65
Mass(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const bool diagonal=false, const bool lumped_diagonal=false, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition mass.hpp:43
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
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 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
constexpr T grad_forward_map_rad(const T r_1, const T r_2)
Definition integrands.hpp:596
constexpr T shape_rad(const int node_idx, const T zeta)
Radial shape function.
Definition integrands.hpp:86
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 T shape_lat(const int node_idx, const T xi, const T eta)
Lateral shape function.
Definition integrands.hpp:118
constexpr T forward_map_rad(const T r_1, const T r_2, const T zeta)
Definition integrands.hpp:579
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 void jacobian_lat_determinant(T(&det_jac_lat)[num_wedges_per_hex_cell][NumQuadPoints], const dense::Vec< T, 3 > wedge_surf_phy_coords[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const dense::Vec< T, 3 > quad_points[NumQuadPoints])
Like jacobian_lat_inverse_transposed_and_determinant() but only computes the determinant (cheaper if ...
Definition kernel_helpers.hpp:158
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.
@ Add
Add to the destination vector.
OperatorCommunicationMode
Modes for communication during operator application.
Definition operator.hpp:40
@ SkipCommunication
Do not communicate.
@ CommunicateAdditively
Communicate and add results.
constexpr Mat diagonal() const
Definition mat.hpp:377
static constexpr Mat diagonal_from_vec(const Vec< T, Rows > &diagonal)
Definition mat.hpp:79
void fill(const T value)
Definition vec.hpp:30