Loading...
Searching...
No Matches
vector_mass.hpp
Go to the documentation of this file.
1
2#pragma once
3
5#include "dense/vec.hpp"
8#include "linalg/operator.hpp"
9#include "linalg/vector.hpp"
10#include "linalg/vector_q1.hpp"
11
13
14template < typename ScalarT, int VecDim = 3 >
16{
17 public:
20 using ScalarType = ScalarT;
21
22 private:
24
27
28 bool diagonal_;
29
30 linalg::OperatorApplyMode operator_apply_mode_;
31 linalg::OperatorCommunicationMode operator_communication_mode_;
32
35
38
39 public:
44 const bool diagonal,
46 linalg::OperatorCommunicationMode operator_communication_mode =
48 : domain_( domain )
49 , grid_( grid )
50 , radii_( radii )
51 , diagonal_( diagonal )
52 , operator_apply_mode_( operator_apply_mode )
53 , operator_communication_mode_( operator_communication_mode )
54 // TODO: we can reuse the send and recv buffers and pass in from the outside somehow
55 , send_buffers_( domain )
56 , recv_buffers_( domain )
57 {}
58
59 void apply_impl( const SrcVectorType& src, DstVectorType& dst )
60 {
61 if ( operator_apply_mode_ == linalg::OperatorApplyMode::Replace )
62 {
63 assign( dst, 0 );
64 }
65
66 src_ = src.grid_data();
67 dst_ = dst.grid_data();
68
69 Kokkos::parallel_for( "matvec", grid::shell::local_domain_md_range_policy_cells( domain_ ), *this );
70
71 if ( operator_communication_mode_ == linalg::OperatorCommunicationMode::CommunicateAdditively )
72 {
74 domain_, dst_, send_buffers_, recv_buffers_ );
76 }
77 }
78
79 KOKKOS_INLINE_FUNCTION void
80 operator()( const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell ) const
81 {
82 // First all the r-independent stuff.
83 // Gather surface points for each wedge.
84
86 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
87
88 // Compute lateral part of Jacobian.
89
90 constexpr auto num_quad_points = quadrature::quad_felippa_3x2_num_quad_points;
91
92 dense::Vec< ScalarT, 3 > quad_points[num_quad_points];
93 ScalarT quad_weights[num_quad_points];
94
97
98 ScalarT det_jac_lat[num_wedges_per_hex_cell][num_quad_points] = {};
99
100 jacobian_lat_determinant( det_jac_lat, wedge_phy_surf, quad_points );
101
102 // Only now we introduce radially dependent terms.
103 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
104 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
105
106 // For now, compute the local element matrix. We'll improve that later.
108
109 const ScalarT grad_r = grad_forward_map_rad( r_1, r_2 );
110
111 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
112 {
113 for ( int q = 0; q < num_quad_points; q++ )
114 {
115 const ScalarT r = forward_map_rad( r_1, r_2, quad_points[q]( 2 ) );
116
117 for ( int i = 0; i < num_nodes_per_wedge; i++ )
118 {
119 for ( int j = 0; j < num_nodes_per_wedge; j++ )
120 {
121 const ScalarT shape_i = shape( i, quad_points[q] );
122 const ScalarT shape_j = shape( j, quad_points[q] );
123
124 A[wedge]( i, j ) +=
125 quad_weights[q] * ( shape_i * shape_j * r * r * grad_r * det_jac_lat[wedge][q] );
126 }
127 }
128 }
129 }
130
131 if ( diagonal_ )
132 {
133 A[0] = A[0].diagonal();
134 A[1] = A[1].diagonal();
135 }
136
137 for ( int d = 0; d < VecDim; d++ )
138 {
140 extract_local_wedge_vector_coefficients( src, local_subdomain_id, x_cell, y_cell, r_cell, d, src_ );
141
143
144 dst[0] = A[0] * src[0];
145 dst[1] = A[1] * src[1];
146
147 atomically_add_local_wedge_vector_coefficients( dst_, local_subdomain_id, x_cell, y_cell, r_cell, d, dst );
148 }
149 }
150};
151
153
154} // namespace terra::fe::wedge::operators::shell
ScalarT ScalarType
Definition vector_mass.hpp:20
VectorMass(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const bool diagonal, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition vector_mass.hpp:40
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition vector_mass.hpp:59
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition vector_mass.hpp:80
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
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
constexpr T grad_forward_map_rad(const T r_1, const T r_2)
Definition integrands.hpp:596
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
void atomically_add_local_wedge_vector_coefficients(const grid::Grid4DDataVec< T, VecDim > &global_coefficients, const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int d, 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:465
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_vector_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 int d, const grid::Grid4DDataVec< T, VecDim > &global_coefficients)
Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient v...
Definition kernel_helpers.hpp:356
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
constexpr T shape(const int node_idx, const T xi, const T eta, const T zeta)
(Tensor-product) Shape function.
Definition integrands.hpp:146
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 > 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.
Definition mat.hpp:10
constexpr Mat diagonal() const
Definition mat.hpp:377