Loading...
Searching...
No Matches
div_k_grad.hpp
Go to the documentation of this file.
1
2#pragma once
3
4#include "../../quadrature/quadrature.hpp"
6#include "dense/vec.hpp"
9#include "grid/grid_types.hpp"
11#include "linalg/operator.hpp"
13#include "linalg/vector.hpp"
14#include "linalg/vector_q1.hpp"
15
17
18template < typename ScalarT >
20{
21 public:
24 using ScalarType = ScalarT;
25 static constexpr int LocalMatrixDim = 6;
27
28 private:
29 LocalMatrixStorage local_matrix_storage_;
30
31 bool single_quadpoint_ = false;
32
34
39
40 bool treat_boundary_;
41 bool diagonal_;
42
43 linalg::OperatorApplyMode operator_apply_mode_;
44 linalg::OperatorCommunicationMode operator_communication_mode_;
45 linalg::OperatorStoredMatrixMode operator_stored_matrix_mode_;
46
49
52
54 ScalarT quad_weights_3x2_[quadrature::quad_felippa_3x2_num_quad_points];
56 ScalarT quad_weights_1x1_[quadrature::quad_felippa_1x1_num_quad_points];
57
58 public:
65 bool treat_boundary,
66 bool diagonal,
68 linalg::OperatorCommunicationMode operator_communication_mode =
71 : domain_( domain )
72 , grid_( grid )
73 , radii_( radii )
74 , mask_( mask )
75 , k_( k )
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 )
81 // TODO: we can reuse the send and recv buffers and pass in from the outside somehow
82 , send_buffers_( domain )
83 , recv_buffers_( domain )
84 {
89 }
90
91 /// @brief Getter for mask member
92 KOKKOS_INLINE_FUNCTION
94 const int local_subdomain_id,
95 const int x_cell,
96 const int y_cell,
97 const int r_cell,
99 {
100 return util::has_flag( mask_( local_subdomain_id, x_cell, y_cell, r_cell ), flag );
101 }
102
104 const linalg::OperatorApplyMode operator_apply_mode,
105 const linalg::OperatorCommunicationMode operator_communication_mode )
106 {
107 operator_apply_mode_ = operator_apply_mode;
108 operator_communication_mode_ = operator_communication_mode;
109 }
110
111 /// @brief S/Getter for diagonal member
112 void set_diagonal( bool v ) { diagonal_ = v; }
113
114 /// @brief Getter for coefficient
116
117 /// @brief Getter for domain member
118 const grid::shell::DistributedDomain& get_domain() const { return domain_; }
119
120 /// @brief Getter for radii member
122
123 /// @brief Getter for grid member
125
126 /// @brief S/Getter for quadpoint member
127 void set_single_quadpoint( bool v ) { single_quadpoint_ = v; }
128
130 linalg::OperatorStoredMatrixMode operator_stored_matrix_mode,
131 int level_range,
133 {
134 operator_stored_matrix_mode_ = operator_stored_matrix_mode;
135
136 // allocate storage if necessary
137 if ( operator_stored_matrix_mode_ != linalg::OperatorStoredMatrixMode::Off )
138 {
140 domain_, operator_stored_matrix_mode_, level_range, GCAElements );
141 }
142 }
143
144 linalg::OperatorStoredMatrixMode get_stored_matrix_mode() { return operator_stored_matrix_mode_; }
145
146 /// @brief Set the local matrix stored in the operator
147 KOKKOS_INLINE_FUNCTION
149 const int local_subdomain_id,
150 const int x_cell,
151 const int y_cell,
152 const int r_cell,
153 const int wedge,
155 {
156 // request from storage
157 KOKKOS_ASSERT( operator_stored_matrix_mode_ != linalg::OperatorStoredMatrixMode::Off );
158 local_matrix_storage_.set_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge, mat );
159 }
160
161 /// @brief Retrives the local matrix
162 /// if there is stored local matrices, the desired local matrix is loaded and returned
163 /// if not, the local matrix is assembled on-the-fly
164 KOKKOS_INLINE_FUNCTION
166 const int local_subdomain_id,
167 const int x_cell,
168 const int y_cell,
169 const int r_cell,
170 const int wedge ) const
171 {
172 // request from storage
173 if ( operator_stored_matrix_mode_ != linalg::OperatorStoredMatrixMode::Off )
174 {
175 if ( !local_matrix_storage_.has_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge ) )
176 {
177 Kokkos::abort( "No matrix found at that spatial index." );
178 }
179 return local_matrix_storage_.get_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge );
180 }
181 else
182 {
183 return assemble_local_matrix( local_subdomain_id, x_cell, y_cell, r_cell, wedge );
184 }
185 }
186
187 void apply_impl( const SrcVectorType& src, DstVectorType& dst )
188 {
189 if ( operator_apply_mode_ == linalg::OperatorApplyMode::Replace )
190 {
191 assign( dst, 0 );
192 }
193
194 src_ = src.grid_data();
195 dst_ = dst.grid_data();
196
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 ) )
199 {
200 throw std::runtime_error( "LaplaceSimple: src/dst mismatch" );
201 }
202
203 if ( src_.extent( 1 ) != grid_.extent( 1 ) || src_.extent( 2 ) != grid_.extent( 2 ) )
204 {
205 throw std::runtime_error( "LaplaceSimple: src/dst mismatch" );
206 }
207
208 Kokkos::parallel_for( "matvec", grid::shell::local_domain_md_range_policy_cells( domain_ ), *this );
209
210 if ( operator_communication_mode_ == linalg::OperatorCommunicationMode::CommunicateAdditively )
211 {
213 domain_, dst_, send_buffers_, recv_buffers_ );
215 }
216 }
217
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
220 {
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 };
224
225 // use stored matrices (at least on some elements)
226 if ( operator_stored_matrix_mode_ != linalg::OperatorStoredMatrixMode::Off )
227 {
229
230 if ( operator_stored_matrix_mode_ == linalg::OperatorStoredMatrixMode::Full )
231 {
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 );
234 }
235 else if ( operator_stored_matrix_mode_ == linalg::OperatorStoredMatrixMode::Selective )
236 {
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 ) )
239 {
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 );
242 }
243 else
244 {
245 // Kokkos::abort("Matrix not found.");
246 A[0] = assemble_local_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 0 );
247 A[1] = assemble_local_matrix( local_subdomain_id, x_cell, y_cell, r_cell, 1 );
248 }
249 }
250
251 if ( diagonal_ )
252 {
253 A[0] = A[0].diagonal();
254 A[1] = A[1].diagonal();
255 }
256
257 if ( treat_boundary_ )
258 {
259 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
260 {
261 dense::Mat< ScalarT, 6, 6 > boundary_mask;
262 boundary_mask.fill( 1.0 );
263 if ( r_cell == 0 )
264 {
265 // Inner boundary (CMB).
266 for ( int i = 0; i < 6; i++ )
267 {
268 for ( int j = 0; j < 6; j++ )
269 {
270 if ( i != j && ( i < 3 || j < 3 ) )
271 {
272 boundary_mask( i, j ) = 0.0;
273 }
274 }
275 }
276 }
277
278 if ( r_cell + 1 == radii_.extent( 1 ) - 1 )
279 {
280 // Outer boundary (surface).
281 for ( int i = 0; i < 6; i++ )
282 {
283 for ( int j = 0; j < 6; j++ )
284 {
285 if ( i != j && ( i >= 3 || j >= 3 ) )
286 {
287 boundary_mask( i, j ) = 0.0;
288 }
289 }
290 }
291 }
292
293 A[wedge].hadamard_product( boundary_mask );
294 }
295 }
296
298 extract_local_wedge_scalar_coefficients( src, local_subdomain_id, x_cell, y_cell, r_cell, src_ );
299
301
302 dst[0] = A[0] * src[0];
303 dst[1] = A[1] * src[1];
304
305 atomically_add_local_wedge_scalar_coefficients( dst_, local_subdomain_id, x_cell, y_cell, r_cell, dst );
306 }
307 else
308 {
309 // assemble on-the-fly
310
311 // Gather surface points for each wedge.
313 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
314
315 // Gather wedge radii.
316 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
317 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
318
319 // Quadrature points.
320 int num_quad_points = single_quadpoint_ ? quadrature::quad_felippa_1x1_num_quad_points :
322
324 extract_local_wedge_scalar_coefficients( k_local_hex, local_subdomain_id, x_cell, y_cell, r_cell, k_ );
325
326 ScalarType src_local_hex[8] = { 0 };
327 ScalarType dst_local_hex[8] = { 0 };
328
329 for ( int i = 0; i < 8; i++ )
330 {
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] );
333 }
334
335 // Compute the local element matrix.
336
337 for ( int q = 0; q < num_quad_points; q++ )
338 {
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];
341
342 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
343 {
345 ScalarType jdet_keval_quadweight = 0;
346
348 wedge, qp, w, r_1, r_2, wedge_phy_surf, k_local_hex, grad, jdet_keval_quadweight );
349
350 // dot of coeff dofs and element-local shape functions to evaluate the coefficent on the current element
351 ScalarType k_eval = 0.0;
352
353 for ( int k = 0; k < num_nodes_per_wedge; k++ )
354 {
355 k_eval += shape( k, qp ) * k_local_hex[wedge]( k );
356 }
357
358 jdet_keval_quadweight *= k_eval;
359
360 fused_local_mv( src_local_hex, dst_local_hex, wedge, jdet_keval_quadweight, grad, r_cell );
361 }
362 }
363
364 for ( int i = 0; i < 8; i++ )
365 {
366 Kokkos::atomic_add(
367 &dst_(
368 local_subdomain_id,
369 x_cell + hex_offset_x[i],
370 y_cell + hex_offset_y[i],
371 r_cell + hex_offset_r[i] ),
372 dst_local_hex[i] );
373 }
374 }
375 }
376
377 KOKKOS_INLINE_FUNCTION void assemble_trial_test_vecs(
378 const int wedge,
379 const dense::Vec< ScalarType, 3 >& quad_point,
380 const ScalarType quad_weight,
381 const ScalarT r_1,
382 const ScalarT r_2,
383 dense::Vec< ScalarT, 3 > ( *wedge_phy_surf )[3],
384 const dense::Vec< ScalarT, 6 >* k_local_hex,
386 ScalarType& jdet_quadweight ) const
387 {
388 dense::Mat< ScalarType, 3, 3 > J = jac( wedge_phy_surf[wedge], r_1, r_2, quad_point );
389 const auto det = J.det();
390 const auto abs_det = Kokkos::abs( det );
391 const dense::Mat< ScalarType, 3, 3 > J_inv_transposed = J.inv_transposed( det );
392
393 for ( int k = 0; k < num_nodes_per_wedge; k++ )
394 {
395 grad[k] = J_inv_transposed * grad_shape( k, quad_point );
396 }
397 jdet_quadweight = quad_weight * abs_det;
398 }
399
400 /// @brief assemble the local matrix and return it for a given element, wedge, and vectorial component
401 /// (determined by dimi, dimj)
402 KOKKOS_INLINE_FUNCTION
404 const int local_subdomain_id,
405 const int x_cell,
406 const int y_cell,
407 const int r_cell,
408 const int wedge ) const
409 {
410 // Gather surface points for each wedge.
411 // TODO gather this for only 1 wedge
413 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
414
415 // Gather wedge radii.
416 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
417 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
418
420 extract_local_wedge_scalar_coefficients( k_local_hex, local_subdomain_id, x_cell, y_cell, r_cell, k_ );
421
422 // Compute the local element matrix.
424 int num_quad_points = single_quadpoint_ ? quadrature::quad_felippa_1x1_num_quad_points :
426
427 // spatial dimensions: quadrature points and wedge
428 for ( int q = 0; q < num_quad_points; q++ )
429 {
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];
432
434 ScalarType jdet_keval_quadweight = 0;
436 wedge, qp, w, r_1, r_2, wedge_phy_surf, k_local_hex, grad, jdet_keval_quadweight );
437
438 // dot of coeff dofs and element-local shape functions to evaluate the coefficent on the current element
439 ScalarType k_eval = 0.0;
440
441 for ( int k = 0; k < num_nodes_per_wedge; k++ )
442 {
443 k_eval += shape( k, qp ) * k_local_hex[wedge]( k );
444 }
445
446 jdet_keval_quadweight *= k_eval;
447
448 // propagate on local matrix by outer product of test and trial vecs
449 for ( int i = 0; i < num_nodes_per_wedge; i++ )
450 {
451 for ( int j = 0; j < num_nodes_per_wedge; j++ )
452 {
453 A( i, j ) += jdet_keval_quadweight * grad[j].dot( grad[i] );
454 // for the div, we just extract the component from the gradient vector
455 }
456 }
457 }
458
459 if ( treat_boundary_ )
460 {
462 boundary_mask.fill( 1.0 );
463
464 if ( r_cell == 0 )
465 {
466 // Inner boundary (CMB).
467 for ( int i = 0; i < 6; i++ )
468 {
469 for ( int j = 0; j < 6; j++ )
470 {
471 // on diagonal components of the vectorial diffusion operator, we exclude the diagonal entries from elimination
472 if ( i != j && ( i < 3 || j < 3 ) )
473 {
474 boundary_mask( i, j ) = 0.0;
475 }
476 }
477 }
478 }
479
480 if ( r_cell + 1 == radii_.extent( 1 ) - 1 )
481 {
482 // Outer boundary (surface).
483 for ( int i = 0; i < 6; i++ )
484 {
485 for ( int j = 0; j < 6; j++ )
486 {
487 // on diagonal components of the vectorial diffusion operator, we exclude the diagonal entries from elimination
488 if ( i != j && ( i >= 3 || j >= 3 ) )
489 {
490 boundary_mask( i, j ) = 0.0;
491 }
492 }
493 }
494 }
495 A.hadamard_product( boundary_mask );
496 }
497
498 return A;
499 }
500
501 // executes the fused local matvec on an element, given the assembled trial and test vectors
502 KOKKOS_INLINE_FUNCTION void fused_local_mv(
503 ScalarType src_local_hex[8],
504 ScalarType dst_local_hex[8],
505 const int wedge,
506 const ScalarType jdet_keval_quadweight,
508 int r_cell ) const
509 {
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 } };
513
515 grad_u.fill( 0.0 );
516
517 const bool at_cmb = r_cell == 0;
518 const bool at_surface = r_cell + 1 == radii_.extent( 1 ) - 1;
519 int cmb_shift = 0;
520 int surface_shift = 0;
521
522 // Compute ∇u at this quadrature point.
523 if ( !diagonal_ )
524 {
525 if ( treat_boundary_ && at_cmb )
526 {
527 // at the core-mantle boundary, we exclude dofs that are lower-indexed than the dof on the boundary
528 cmb_shift = 3;
529 }
530 else if ( treat_boundary_ && at_surface )
531 {
532 // at the surface boundary, we exclude dofs that are higher-indexed than the dof on the boundary
533 surface_shift = 3;
534 }
535
536 // accumulate the element-local gradient/divergence of the trial function (loop over columns of local matrix/local dofs)
537 // by dot of trial vec and src dofs
538 for ( int i = 0 + cmb_shift; i < num_nodes_per_wedge - surface_shift; i++ )
539 {
540 grad_u = grad_u +
541 grad[i] * src_local_hex[4 * offset_r[wedge][i] + 2 * offset_y[wedge][i] + offset_x[wedge][i]];
542 }
543
544 // Add the test function contributions.
545 // for each row of the local matrix (test-functions):
546 // multiply trial part (fully assembled for the current element from loop above) with test part corresponding to the current row/dof
547 // += due to contributions from other elements
548 for ( int j = 0 + cmb_shift; j < num_nodes_per_wedge - surface_shift; j++ )
549 {
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 );
552 // for the div, we just extract the component from the gradient vector
553 }
554 }
555
556 // Dirichlet DoFs are only to be eliminated on diagonal blocks of epsilon
557 if ( diagonal_ || ( ( treat_boundary_ && ( at_cmb || at_surface ) ) ) )
558 {
559 // for the diagonal elements at the boundary, we switch the shifts
560 for ( int i = 0 + surface_shift; i < num_nodes_per_wedge - cmb_shift; i++ )
561 {
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]];
564
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 );
567 }
568 }
569 }
570};
571
574
575} // namespace terra::fe::wedge::operators::shell
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
Definition mat.hpp:10
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