Loading...
Searching...
No Matches
epsilon_simple.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "../../quadrature/quadrature.hpp"
5#include "dense/vec.hpp"
9#include "linalg/operator.hpp"
10#include "linalg/vector.hpp"
11#include "linalg/vector_q1.hpp"
12#include "util/timer.hpp"
13
15
16template < typename ScalarT, int VecDim = 3 >
18{
19 public:
22 using ScalarType = ScalarT;
23
24 private:
26
30
31 bool treat_boundary_;
32 bool diagonal_;
33
34 linalg::OperatorApplyMode operator_apply_mode_;
35 linalg::OperatorCommunicationMode operator_communication_mode_;
36
39
42
43 public:
49 bool treat_boundary,
50 bool diagonal,
52 linalg::OperatorCommunicationMode operator_communication_mode =
54 : domain_( domain )
55 , grid_( grid )
56 , radii_( radii )
57 , k_( k )
58 , treat_boundary_( treat_boundary )
59 , diagonal_( diagonal )
60 , operator_apply_mode_( operator_apply_mode )
61 , operator_communication_mode_( operator_communication_mode )
62 // TODO: we can reuse the send and recv buffers and pass in from the outside somehow
63 , send_buffers_( domain )
64 , recv_buffers_( domain )
65 {}
66
68
70 const linalg::OperatorApplyMode operator_apply_mode,
71 const linalg::OperatorCommunicationMode operator_communication_mode )
72 {
73 operator_apply_mode_ = operator_apply_mode;
74 operator_communication_mode_ = operator_communication_mode;
75 }
76
77 void apply_impl( const SrcVectorType& src, DstVectorType& dst )
78 {
79 util::Timer timer_apply( "vector_laplace_apply" );
80
81 if ( operator_apply_mode_ == linalg::OperatorApplyMode::Replace )
82 {
83 assign( dst, 0 );
84 }
85
86 src_ = src.grid_data();
87 dst_ = dst.grid_data();
88
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 ) )
91 {
92 throw std::runtime_error( "VectorLaplace: src/dst mismatch" );
93 }
94
95 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
96 {
97 throw std::runtime_error( "VectorLaplace: src/dst mismatch" );
98 }
99
100 util::Timer timer_kernel( "vector_laplace_kernel" );
101 Kokkos::parallel_for( "matvec", grid::shell::local_domain_md_range_policy_cells( domain_ ), *this );
102 Kokkos::fence();
103 timer_kernel.stop();
104
105 if ( operator_communication_mode_ == linalg::OperatorCommunicationMode::CommunicateAdditively )
106 {
107 util::Timer timer_comm( "vector_laplace_comm" );
108
110 domain_, dst_, send_buffers_, recv_buffers_ );
112 }
113 }
114
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
117 {
118 // Gather surface points for each wedge.
120 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
121
122 // Gather wedge radii.
123 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
124 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
125
126 // Quadrature points.
127 constexpr auto num_quad_points = quadrature::quad_felippa_3x2_num_quad_points;
128
129 dense::Vec< ScalarT, 3 > quad_points[num_quad_points];
130 ScalarT quad_weights[num_quad_points];
131
134
136 extract_local_wedge_scalar_coefficients( k, local_subdomain_id, x_cell, y_cell, r_cell, k_ );
137
138 // Compute the local element matrix.
140
141 // FE dimensions: velocity coupling components of epsilon operator
142 for ( int dimi = 0; dimi < 3; ++dimi )
143 {
144 for ( int dimj = 0; dimj < 3; ++dimj )
145 {
146 if ( diagonal_ and dimi != dimj )
147 continue;
148
149 // spatial dimensions: quadrature points and wedge
150 for ( int q = 0; q < num_quad_points; q++ )
151 {
152 const auto w = quad_weights[q];
153
154 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
155 {
156 dense::Mat< ScalarT, 3, 3 > J = jac( wedge_phy_surf[wedge], r_1, r_2, quad_points[q] );
157 const auto det = J.det();
158 const auto abs_det = Kokkos::abs( det );
159 const auto J_inv_transposed = J.inv_transposed( det );
160 ScalarType k_eval = 0.0;
161 for ( int j = 0; j < num_nodes_per_wedge; j++ )
162 {
163 k_eval += shape( j, quad_points[q] ) * k[wedge]( j );
164 }
165 // FE dimensions: local DoFs/associated shape functions
166 for ( int i = 0; i < num_nodes_per_wedge; i++ )
167 {
168 // basis functions are vectors with VecDim components -> build tensorial gradients
171 grad_shape( i, quad_points[q] ), dimi );
172 dense::Mat< ScalarT, 3, 3 > sym_grad_i = ( grad_i + grad_i.transposed() );
173 //const auto grad_i = J_inv_transposed * grad_shape( i, quad_points[q] );
174
175 for ( int j = 0; j < num_nodes_per_wedge; j++ )
176 {
179 grad_shape( j, quad_points[q] ), dimj );
180
181 dense::Mat< ScalarT, 3, 3 > sym_grad_j = ( grad_j + grad_j.transposed() );
182
183 //const auto grad_j = J_inv_transposed * grad_shape( j, quad_points[q] );
184
185 A[wedge]( i + num_nodes_per_wedge * dimi, j + num_nodes_per_wedge * dimj ) +=
186 0.5 * w * k_eval * ( ( sym_grad_i ).double_contract( sym_grad_j ) * abs_det );
187 }
188 }
189 }
190 }
191 }
192 }
193
194 if ( treat_boundary_ )
195 {
196 for ( int dimi = 0; dimi < 3; ++dimi )
197 {
198 for ( int dimj = 0; dimj < 3; ++dimj )
199 {
200 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
201 {
202 dense::Mat< ScalarT, 18, 18 > boundary_mask;
203 boundary_mask.fill( 1.0 );
204
205 if ( r_cell == 0 )
206 {
207 // Inner boundary (CMB).
208 for ( int i = 0; i < 6; i++ )
209 {
210 for ( int j = 0; j < 6; j++ )
211 {
212 if ( ( dimi == dimj && i != j && ( i < 3 || j < 3 ) ) or
213 ( dimi != dimj && ( i < 3 || j < 3 ) ) )
214 {
215 boundary_mask(
216 i + num_nodes_per_wedge * dimi, j + num_nodes_per_wedge * dimj ) = 0.0;
217 }
218 }
219 }
220 }
221
222 if ( r_cell + 1 == radii_.extent( 1 ) - 1 )
223 {
224 // Outer boundary (surface).
225 for ( int i = 0; i < 6; i++ )
226 {
227 for ( int j = 0; j < 6; j++ )
228 {
229 if ( ( dimi == dimj && i != j && ( i >= 3 || j >= 3 ) ) or
230 ( dimi != dimj && ( i >= 3 || j >= 3 ) ) )
231 {
232 boundary_mask(
233 i + num_nodes_per_wedge * dimi, j + num_nodes_per_wedge * dimj ) = 0.0;
234 }
235 }
236 }
237 }
238
239 A[wedge].hadamard_product( boundary_mask );
240 }
241 }
242 }
243 }
244
245 if ( diagonal_ )
246 {
247 A[0] = A[0].diagonal();
248 A[1] = A[1].diagonal();
249 }
250
252 for ( int dimj = 0; dimj < 3; dimj++ )
253 {
255 extract_local_wedge_vector_coefficients( src_d, local_subdomain_id, x_cell, y_cell, r_cell, dimj, src_ );
256
257 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
258 {
259 for ( int i = 0; i < num_nodes_per_wedge; i++ )
260 {
261 src[wedge]( dimj * num_nodes_per_wedge + i ) = src_d[wedge]( i );
262 }
263 }
264 }
265 //extract_local_wedge_vector_coefficients( src, local_subdomain_id, x_cell, y_cell, r_cell, dimj, src_ );
266
268
269 dst[0] = A[0] * src[0];
270 dst[1] = A[1] * src[1];
271
272 //atomically_add_local_wedge_vector_coefficients( dst_, local_subdomain_id, x_cell, y_cell, r_cell, dimi, dst );
273 for ( int dimi = 0; dimi < 3; dimi++ )
274 {
276 dst_d[0] = dst[0].template slice< 6 >( dimi * num_nodes_per_wedge );
277 dst_d[1] = dst[1].template slice< 6 >( dimi * num_nodes_per_wedge );
278
280 dst_, local_subdomain_id, x_cell, y_cell, r_cell, dimi, dst_d );
281 }
282 }
283};
284
287
288} // namespace terra::fe::wedge::operators::shell
Send and receive buffers for all process-local subdomain boundaries.
Definition communication.hpp:56
void set_operator_apply_and_communication_modes(const linalg::OperatorApplyMode operator_apply_mode, const linalg::OperatorCommunicationMode operator_communication_mode)
Definition epsilon_simple.hpp:69
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition epsilon_simple.hpp:77
ScalarT ScalarType
Definition epsilon_simple.hpp:22
const grid::Grid4DDataScalar< ScalarType > & k_grid_data()
Definition epsilon_simple.hpp:67
EpsilonSimple(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_simple.hpp:44
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition epsilon_simple.hpp:116
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
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 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
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 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.
Definition mat.hpp:10
void fill(const T value)
Definition mat.hpp:201
constexpr Mat diagonal() const
Definition mat.hpp:377
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
Mat & hadamard_product(const Mat &mat)
Definition mat.hpp:213