Loading...
Searching...
No Matches
laplace.hpp
Go to the documentation of this file.
1
2#pragma once
3
4#include "../../quadrature/quadrature.hpp"
6#include "dense/vec.hpp"
10#include "linalg/operator.hpp"
11#include "linalg/vector.hpp"
12#include "linalg/vector_q1.hpp"
13#include "util/timer.hpp"
14
16
17template < typename ScalarT >
19{
20 public:
23 using ScalarType = ScalarT;
24
25 private:
27
31
32 bool treat_boundary_;
33 bool diagonal_;
34
35 linalg::OperatorApplyMode operator_apply_mode_;
36 linalg::OperatorCommunicationMode operator_communication_mode_;
37
40
43
44 public:
50 bool treat_boundary,
51 bool diagonal,
53 linalg::OperatorCommunicationMode operator_communication_mode =
55 : domain_( domain )
56 , grid_( grid )
57 , radii_( radii )
58 , mask_( mask )
59 , treat_boundary_( treat_boundary )
60 , diagonal_( diagonal )
61 , operator_apply_mode_( operator_apply_mode )
62 , operator_communication_mode_( operator_communication_mode )
63 // TODO: we can reuse the send and recv buffers and pass in from the outside somehow
64 , send_buffers_( domain )
65 , recv_buffers_( domain )
66 {}
67
68 void apply_impl( const SrcVectorType& src, DstVectorType& dst )
69 {
70 util::Timer timer_apply( "laplace_apply" );
71
72 if ( operator_apply_mode_ == linalg::OperatorApplyMode::Replace )
73 {
74 assign( dst, 0 );
75 }
76
77 src_ = src.grid_data();
78 dst_ = dst.grid_data();
79
80 if ( src_.extent( 0 ) != dst_.extent( 0 ) || src_.extent( 1 ) != dst_.extent( 1 ) ||
81 src_.extent( 2 ) != dst_.extent( 2 ) || src_.extent( 3 ) != dst_.extent( 3 ) )
82 {
83 throw std::runtime_error( "LaplaceSimple: src/dst mismatch" );
84 }
85
86 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
87 {
88 throw std::runtime_error( "LaplaceSimple: src/dst mismatch" );
89 }
90
91 util::Timer timer_kernel( "laplace_kernel" );
92 Kokkos::parallel_for( "matvec", grid::shell::local_domain_md_range_policy_cells( domain_ ), *this );
93 Kokkos::fence();
94 timer_kernel.stop();
95
96 if ( operator_communication_mode_ == linalg::OperatorCommunicationMode::CommunicateAdditively )
97 {
98 util::Timer timer_comm( "laplace_comm" );
100 domain_, dst_, send_buffers_, recv_buffers_ );
102 }
103 }
104
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
107 {
108 // Gather surface points for each wedge.
110 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
111
112 // Gather wedge radii.
113 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
114 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
115
116 // Quadrature points.
117 constexpr auto num_quad_points = quadrature::quad_felippa_3x2_num_quad_points;
118
119 dense::Vec< ScalarT, 3 > quad_points[num_quad_points];
120 ScalarT quad_weights[num_quad_points];
121
124
125 // Compute the local element matrix.
126
127 ScalarType src_local_hex[8] = { 0 };
128 ScalarType dst_local_hex[8] = { 0 };
129
130 for ( int i = 0; i < 8; i++ )
131 {
132 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
133 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
134 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
135
136 src_local_hex[i] = src_(
137 local_subdomain_id, x_cell + hex_offset_x[i], y_cell + hex_offset_y[i], r_cell + hex_offset_r[i] );
138 }
139
140 const bool at_bot_boundary =
141 util::has_flag( mask_( local_subdomain_id, x_cell, y_cell, r_cell ), grid::shell::ShellBoundaryFlag::CMB );
142 const bool at_top_boundary = util::has_flag(
143 mask_( local_subdomain_id, x_cell, y_cell, r_cell + 1 ), grid::shell::ShellBoundaryFlag::SURFACE );
144
145 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
146 {
147 for ( int q = 0; q < num_quad_points; q++ )
148 {
149 const auto quad_point = quad_points[q];
150 const auto quad_weight = quad_weights[q];
151
152 // 1. Compute Jacobian and inverse at this quadrature point.
153
154 const auto J = jac( wedge_phy_surf[wedge], r_1, r_2, quad_points[q] );
155 const auto det = J.det();
156 const auto abs_det = Kokkos::abs( det );
157 const auto J_inv_transposed = J.inv_transposed( det );
158
159 // 2. Compute physical gradients for all nodes at this quadrature point.
161 for ( int k = 0; k < num_nodes_per_wedge; k++ )
162 {
163 grad_phy[k] = J_inv_transposed * grad_shape( k, quad_point );
164 }
165
166 if ( diagonal_ )
167 {
168 diagonal( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
169 }
170 else if ( treat_boundary_ && at_bot_boundary )
171 {
172 // Bottom boundary dirichlet
173 dirichlet_bot( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
174 }
175 else if ( treat_boundary_ && at_top_boundary )
176 {
177 // Top boundary dirichlet
178 dirichlet_top( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
179 }
180 else
181 {
182 neumann( src_local_hex, dst_local_hex, wedge, quad_weight, abs_det, grad_phy );
183 }
184 }
185 }
186
187 for ( int i = 0; i < 8; i++ )
188 {
189 constexpr int hex_offset_x[8] = { 0, 1, 0, 1, 0, 1, 0, 1 };
190 constexpr int hex_offset_y[8] = { 0, 0, 1, 1, 0, 0, 1, 1 };
191 constexpr int hex_offset_r[8] = { 0, 0, 0, 0, 1, 1, 1, 1 };
192
193 Kokkos::atomic_add(
194 &dst_(
195 local_subdomain_id, x_cell + hex_offset_x[i], y_cell + hex_offset_y[i], r_cell + hex_offset_r[i] ),
196 dst_local_hex[i] );
197 }
198 }
199
200 KOKKOS_INLINE_FUNCTION void neumann(
201 ScalarType* src_local_hex,
202 ScalarType* dst_local_hex,
203 const int wedge,
204 const ScalarType quad_weight,
205 const ScalarType abs_det,
206 const dense::Vec< ScalarType, 3 >* grad_phy ) const
207 {
208 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
209 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
210 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
211
212 // 3. Compute ∇u at this quadrature point.
214 grad_u.fill( 0.0 );
215 for ( int j = 0; j < num_nodes_per_wedge; j++ )
216 {
217 grad_u = grad_u +
218 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
219 }
220
221 // 4. Add the test function contributions.
222 for ( int i = 0; i < num_nodes_per_wedge; i++ )
223 {
224 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
225 quad_weight * grad_phy[i].dot( grad_u ) * abs_det;
226 }
227 }
228
229 KOKKOS_INLINE_FUNCTION void dirichlet_bot(
230 ScalarType* src_local_hex,
231 ScalarType* dst_local_hex,
232 const int wedge,
233 const ScalarType quad_weight,
234 const ScalarType abs_det,
235 const dense::Vec< ScalarType, 3 >* grad_phy ) const
236 {
237 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
238 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
239 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
240
241 // 3. Compute ∇u at this quadrature point.
243 grad_u.fill( 0.0 );
244 for ( int j = 3; j < num_nodes_per_wedge; j++ )
245 {
246 grad_u = grad_u +
247 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
248 }
249
250 // 4. Add the test function contributions.
251 for ( int i = 3; i < num_nodes_per_wedge; i++ )
252 {
253 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
254 quad_weight * grad_phy[i].dot( grad_u ) * abs_det;
255 }
256
257 // Diagonal for top part
258 for ( int i = 0; i < 3; i++ )
259 {
260 const auto grad_u_diag =
261 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
262
263 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
264 quad_weight * grad_phy[i].dot( grad_u_diag ) * abs_det;
265 }
266 }
267
268 KOKKOS_INLINE_FUNCTION void dirichlet_top(
269 ScalarType* src_local_hex,
270 ScalarType* dst_local_hex,
271 const int wedge,
272 const ScalarType quad_weight,
273 const ScalarType abs_det,
274 const dense::Vec< ScalarType, 3 >* grad_phy ) const
275 {
276 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
277 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
278 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
279
280 // 3. Compute ∇u at this quadrature point.
282 grad_u.fill( 0.0 );
283 for ( int j = 0; j < 3; j++ )
284 {
285 grad_u = grad_u +
286 src_local_hex[4 * offset_r[wedge][j] + 2 * offset_y[wedge][j] + offset_x[wedge][j]] * grad_phy[j];
287 }
288
289 // 4. Add the test function contributions.
290 for ( int i = 0; i < 3; i++ )
291 {
292 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
293 quad_weight * grad_phy[i].dot( grad_u ) * abs_det;
294 }
295
296 // Diagonal for top part
297 for ( int i = 3; i < num_nodes_per_wedge; i++ )
298 {
299 const auto grad_u_diag =
300 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
301
302 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
303 quad_weight * grad_phy[i].dot( grad_u_diag ) * abs_det;
304 }
305 }
306
307 KOKKOS_INLINE_FUNCTION void diagonal(
308 ScalarType* src_local_hex,
309 ScalarType* dst_local_hex,
310 const int wedge,
311 const ScalarType quad_weight,
312 const ScalarType abs_det,
313 const dense::Vec< ScalarType, 3 >* grad_phy ) const
314 {
315 constexpr int offset_x[2][6] = { { 0, 1, 0, 0, 1, 0 }, { 1, 0, 1, 1, 0, 1 } };
316 constexpr int offset_y[2][6] = { { 0, 0, 1, 0, 0, 1 }, { 1, 1, 0, 1, 1, 0 } };
317 constexpr int offset_r[2][6] = { { 0, 0, 0, 1, 1, 1 }, { 0, 0, 0, 1, 1, 1 } };
318
319 // 3. Compute ∇u at this quadrature point.
320 // 4. Add the test function contributions.
321 for ( int i = 0; i < num_nodes_per_wedge; i++ )
322 {
323 const auto grad_u =
324 src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] * grad_phy[i];
325
326 dst_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]] +=
327 quad_weight * grad_phy[i].dot( grad_u ) * abs_det;
328 }
329 }
330};
331
334
335} // namespace terra::fe::wedge::operators::shell
void neumann(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:200
void dirichlet_bot(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:229
ScalarT ScalarType
Definition laplace.hpp:23
void diagonal(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:307
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition laplace.hpp:68
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition laplace.hpp:106
Laplace(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &mask, bool treat_boundary, bool diagonal, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition laplace.hpp:45
void dirichlet_top(ScalarType *src_local_hex, ScalarType *dst_local_hex, const int wedge, const ScalarType quad_weight, const ScalarType abs_det, const dense::Vec< ScalarType, 3 > *grad_phy) const
Definition laplace.hpp:268
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
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
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 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 ****, 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 bool has_flag(E mask_value, E flag) noexcept
Checks if a bitmask value contains a specific flag.
Definition bit_masking.hpp:43
T dot(const Vec &other) const
Definition vec.hpp:39