Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
portable_tensor_product_kernels.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii__tensor_product_kernels_h
17#define dealii__tensor_product_kernels_h
18
19#include <deal.II/base/config.h>
20
22
23#include <deal.II/matrix_free/cuda_matrix_free.templates.h>
24
26
27
28namespace Portable
29{
30 namespace internal
31 {
38 // TODO: for now only the general variant is implemented
45
46
47
48#if KOKKOS_VERSION >= 40000
52 template <int n_q_points_1d,
53 typename Number,
54 int direction,
55 bool dof_to_quad,
56 bool add,
57 bool in_place,
58 typename ViewTypeIn,
59 typename ViewTypeOut>
61 apply_1d(const Kokkos::TeamPolicy<
62 MemorySpace::Default::kokkos_space::execution_space>::member_type
63 &team_member,
64 const Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
65 shape_data,
66 const ViewTypeIn in,
67 ViewTypeOut out)
68 {
69 Number t[n_q_points_1d];
70 Kokkos::parallel_for(Kokkos::TeamThreadRange(team_member, n_q_points_1d),
71 [&](const int &q) {
72 t[q] = 0;
73 // This loop simply multiplies the shape function
74 // at the quadrature point by the value finite
75 // element coefficient.
76 // FIXME check why using parallel_reduce
77 // ThreadVector is slower
78 for (int k = 0; k < n_q_points_1d; ++k)
79 {
80 const unsigned int shape_idx =
81 dof_to_quad ? (q + k * n_q_points_1d) :
82 (k + q * n_q_points_1d);
83 const unsigned int source_idx = k;
84 t[q] += shape_data[shape_idx] * in(source_idx);
85 }
86 });
87
88 if constexpr (in_place)
89 team_member.team_barrier();
90
91 Kokkos::parallel_for(Kokkos::TeamThreadRange(team_member, n_q_points_1d),
92 [&](const int &q) {
93 const unsigned int destination_idx = q;
94 if constexpr (add)
95 Kokkos::atomic_add(&out(destination_idx), t[q]);
96 else
97 out(destination_idx) = t[q];
98 });
99 }
100
101
102
106 template <int n_q_points_1d,
107 typename Number,
108 int direction,
109 bool dof_to_quad,
110 bool add,
111 bool in_place,
112 typename ViewTypeIn,
113 typename ViewTypeOut>
115 apply_2d(const Kokkos::TeamPolicy<
116 MemorySpace::Default::kokkos_space::execution_space>::member_type
117 &team_member,
118 const Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
119 shape_data,
120 const ViewTypeIn in,
121 ViewTypeOut out)
122 {
123 using TeamType = Kokkos::TeamPolicy<
124 MemorySpace::Default::kokkos_space::execution_space>::member_type;
125 constexpr unsigned int n_q_points = Utilities::pow(n_q_points_1d, 2);
126
127 Number t[n_q_points];
128 auto thread_policy =
129 Kokkos::TeamThreadMDRange<Kokkos::Rank<2>, TeamType>(team_member,
130 n_q_points_1d,
131 n_q_points_1d);
132 Kokkos::parallel_for(thread_policy, [&](const int i, const int j) {
133 int q_point = i + j * n_q_points_1d;
134
135 // This loop simply multiplies the shape function at the quadrature
136 // point by the value finite element coefficient.
137 // FIXME check why using parallel_reduce ThreadVector is slower
138 const int base_shape = dof_to_quad ? j : j * n_q_points_1d;
139 const int stride_shape = dof_to_quad ? n_q_points_1d : 1;
140 const int base_in = (direction == 0) ? (n_q_points_1d * i) : i;
141 const int stride_in = Utilities::pow(n_q_points_1d, direction);
142 Number sum = shape_data[base_shape] * in(base_in);
143 for (int k = 1; k < n_q_points_1d; ++k)
144 {
145 sum += shape_data[base_shape + k * stride_shape] *
146 in(base_in + k * stride_in);
147 }
148 t[q_point] = sum;
149 });
150
151 if (in_place)
152 team_member.team_barrier();
153
154 Kokkos::parallel_for(thread_policy, [&](const int i, const int j) {
155 const int q_point = i + j * n_q_points_1d;
156 const unsigned int destination_idx =
157 (direction == 0) ? (j + n_q_points_1d * i) : (i + n_q_points_1d * j);
158
159 if (add)
160 Kokkos::atomic_add(&out(destination_idx), t[q_point]);
161 else
162 out(destination_idx) = t[q_point];
163 });
164 }
165
166
167
171 template <int n_q_points_1d,
172 typename Number,
173 int direction,
174 bool dof_to_quad,
175 bool add,
176 bool in_place,
177 typename ViewTypeIn,
178 typename ViewTypeOut>
180 apply_3d(const Kokkos::TeamPolicy<
181 MemorySpace::Default::kokkos_space::execution_space>::member_type
182 &team_member,
183 const Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
184 shape_data,
185 const ViewTypeIn in,
186 ViewTypeOut out)
187 {
188 using TeamType = Kokkos::TeamPolicy<
189 MemorySpace::Default::kokkos_space::execution_space>::member_type;
190 constexpr unsigned int n_q_points = Utilities::pow(n_q_points_1d, 3);
191
192 Number t[n_q_points];
193 auto thread_policy = Kokkos::TeamThreadMDRange<Kokkos::Rank<3>, TeamType>(
194 team_member, n_q_points_1d, n_q_points_1d, n_q_points_1d);
195 Kokkos::parallel_for(
196 thread_policy, [&](const int i, const int j, const int q) {
197 const int q_point =
198 i + j * n_q_points_1d + q * n_q_points_1d * n_q_points_1d;
199
200 // This loop simply multiplies the shape function at the quadrature
201 // point by the value finite element coefficient.
202 // FIXME check why using parallel_reduce ThreadVector is slower
203 const int base_shape = dof_to_quad ? q : q * n_q_points_1d;
204 const int stride_shape = dof_to_quad ? n_q_points_1d : 1;
205 const int base_in =
206 (direction == 0 ?
207 (n_q_points_1d * (i + n_q_points_1d * j)) :
208 (direction == 1 ? (i + n_q_points_1d * n_q_points_1d * j) :
209 (i + n_q_points_1d * j)));
210 const int stride_in = Utilities::pow(n_q_points_1d, direction);
211 Number sum = shape_data[base_shape] * in(base_in);
212 for (int k = 1; k < n_q_points_1d; ++k)
213 {
214 sum += shape_data[base_shape + k * stride_shape] *
215 in(base_in + k * stride_in);
216 }
217 t[q_point] = sum;
218 });
219
220 if (in_place)
221 team_member.team_barrier();
222
223 Kokkos::parallel_for(
224 thread_policy, [&](const int i, const int j, const int q) {
225 const int q_point =
226 i + j * n_q_points_1d + q * n_q_points_1d * n_q_points_1d;
227 const unsigned int destination_idx =
228 (direction == 0) ? (q + n_q_points_1d * (i + n_q_points_1d * j)) :
229 (direction == 1) ? (i + n_q_points_1d * (q + n_q_points_1d * j)) :
230 (i + n_q_points_1d * (j + n_q_points_1d * q));
231
232 if (add)
233 Kokkos::atomic_add(&out(destination_idx), t[q_point]);
234 else
235 out(destination_idx) = t[q_point];
236 });
237 }
238#endif
239
240
241
245 template <int dim,
246 int n_q_points_1d,
247 typename Number,
248 int direction,
249 bool dof_to_quad,
250 bool add,
251 bool in_place,
252 typename ViewTypeIn,
253 typename ViewTypeOut>
255 apply(const Kokkos::TeamPolicy<
256 MemorySpace::Default::kokkos_space::execution_space>::member_type
257 &team_member,
258 const Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
259 shape_data,
260 const ViewTypeIn in,
261 ViewTypeOut out)
262 {
263#if KOKKOS_VERSION >= 40000
264 if constexpr (dim == 1)
265 apply_1d<n_q_points_1d, Number, direction, dof_to_quad, add, in_place>(
266 team_member, shape_data, in, out);
267 if constexpr (dim == 2)
268 apply_2d<n_q_points_1d, Number, direction, dof_to_quad, add, in_place>(
269 team_member, shape_data, in, out);
270 if constexpr (dim == 3)
271 apply_3d<n_q_points_1d, Number, direction, dof_to_quad, add, in_place>(
272 team_member, shape_data, in, out);
273#else
274 constexpr unsigned int n_q_points = Utilities::pow(n_q_points_1d, dim);
275
276 Number t[n_q_points];
277 Kokkos::parallel_for(
278 Kokkos::TeamThreadRange(team_member, n_q_points),
279 [&](const int &q_point) {
280 const unsigned int i = (dim == 1) ? 0 : q_point % n_q_points_1d;
281 const unsigned int j =
282 (dim == 3) ? (q_point / n_q_points_1d) % n_q_points_1d : 0;
283 const unsigned int q =
284 (dim == 1) ? q_point :
285 (dim == 2) ? (q_point / n_q_points_1d) % n_q_points_1d :
286 q_point / (n_q_points_1d * n_q_points_1d);
287
288 // This loop simply multiplies the shape function at the quadrature
289 // point by the value finite element coefficient.
290 const int stride_shape = dof_to_quad ? n_q_points_1d : 1;
291 const int stride = Utilities::pow(n_q_points_1d, direction);
292 const int base_shape = dof_to_quad ? q : (q * n_q_points_1d);
293 const int base =
294 (direction == 0) ? (n_q_points_1d * (i + n_q_points_1d * j)) :
295 (direction == 1) ? (i + n_q_points_1d * (n_q_points_1d * j)) :
296 (i + n_q_points_1d * j);
297 Number sum =
298 shape_data[base_shape] * (in_place ? out(base) : in(base));
299 for (int k = 1; k < n_q_points_1d; ++k)
300 {
301 sum +=
302 shape_data[base_shape + k * stride_shape] *
303 (in_place ? out(base + k * stride) : in(base + k * stride));
304 }
305 t[q_point] = sum;
306 });
307
308 if (in_place)
309 team_member.team_barrier();
310
311 Kokkos::parallel_for(
312 Kokkos::TeamThreadRange(team_member, n_q_points),
313 [&](const int &q_point) {
314 const unsigned int i = (dim == 1) ? 0 : q_point % n_q_points_1d;
315 const unsigned int j =
316 (dim == 3) ? (q_point / n_q_points_1d) % n_q_points_1d : 0;
317 const unsigned int q =
318 (dim == 1) ? q_point :
319 (dim == 2) ? (q_point / n_q_points_1d) % n_q_points_1d :
320 q_point / (n_q_points_1d * n_q_points_1d);
321
322 const unsigned int destination_idx =
323 (direction == 0) ? (q + n_q_points_1d * (i + n_q_points_1d * j)) :
324 (direction == 1) ? (i + n_q_points_1d * (q + n_q_points_1d * j)) :
325 (i + n_q_points_1d * (j + n_q_points_1d * q));
326
327 if (add)
328 Kokkos::atomic_add(&out(destination_idx), t[q_point]);
329 else
330 out(destination_idx) = t[q_point];
331 });
332#endif
333 }
334
335
341 template <EvaluatorVariant variant,
342 int dim,
343 int fe_degree,
344 int n_q_points_1d,
345 typename Number>
348
349
350
357 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
359 dim,
360 fe_degree,
361 n_q_points_1d,
362 Number>
363 {
364 public:
365 using TeamHandle = Kokkos::TeamPolicy<
366 MemorySpace::Default::kokkos_space::execution_space>::member_type;
367
370 const TeamHandle &team_member,
371 Kokkos::View<Number *, MemorySpace::Default::kokkos_space> shape_values,
372 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
373 shape_gradients,
374 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
375 co_shape_gradients);
376
380 template <typename ViewType>
382 evaluate_values(ViewType u);
383
388 template <typename ViewTypeIn, typename ViewTypeOut>
390 evaluate_gradients(const ViewTypeIn u, ViewTypeOut grad_u);
391
396 template <typename ViewType1, typename ViewType2>
398 evaluate_values_and_gradients(ViewType1 u, ViewType2 grad_u);
399
403 template <typename ViewType>
405 integrate_values(ViewType u);
406
411 template <bool add, typename ViewType1, typename ViewType2>
413 integrate_gradients(ViewType1 u, ViewType2 grad_u);
414
419 template <typename ViewType1, typename ViewType2>
421 integrate_values_and_gradients(ViewType1 u, ViewType2 grad_u);
422
427 template <int direction,
428 bool dof_to_quad,
429 bool add,
430 bool in_place,
431 typename ViewTypeIn,
432 typename ViewTypeOut>
434 values(const ViewTypeIn in, ViewTypeOut out) const;
435
440 template <int direction,
441 bool dof_to_quad,
442 bool add,
443 bool in_place,
444 typename ViewTypeIn,
445 typename ViewTypeOut>
447 gradients(const ViewTypeIn in, ViewTypeOut out) const;
448
449 public:
454 template <int direction,
455 bool dof_to_quad,
456 bool add,
457 bool in_place,
458 typename ViewTypeIn,
459 typename ViewTypeOut>
461 co_gradients(const ViewTypeIn in, ViewTypeOut out) const;
462
467
471 Kokkos::View<Number *, MemorySpace::Default::kokkos_space> shape_values;
472
476 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
478
482 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
484 };
485
486
487
488 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
491 dim,
492 fe_degree,
493 n_q_points_1d,
494 Number>::
495 EvaluatorTensorProduct(
496 const TeamHandle &team_member,
497 Kokkos::View<Number *, MemorySpace::Default::kokkos_space> shape_values,
498 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
499 shape_gradients,
500 Kokkos::View<Number *, MemorySpace::Default::kokkos_space>
501 co_shape_gradients)
502 : team_member(team_member)
503 , shape_values(shape_values)
504 , shape_gradients(shape_gradients)
505 , co_shape_gradients(co_shape_gradients)
506 {}
507
508
509
510 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
511 template <int direction,
512 bool dof_to_quad,
513 bool add,
514 bool in_place,
515 typename ViewTypeIn,
516 typename ViewTypeOut>
519 dim,
520 fe_degree,
521 n_q_points_1d,
522 Number>::values(const ViewTypeIn in,
523 ViewTypeOut out) const
524 {
526 team_member, shape_values, in, out);
527 }
528
529
530
531 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
532 template <int direction,
533 bool dof_to_quad,
534 bool add,
535 bool in_place,
536 typename ViewTypeIn,
537 typename ViewTypeOut>
540 dim,
541 fe_degree,
542 n_q_points_1d,
543 Number>::gradients(const ViewTypeIn in,
544 ViewTypeOut out) const
545 {
547 team_member, shape_gradients, in, out);
548 }
549
550
551
552 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
553 template <int direction,
554 bool dof_to_quad,
555 bool add,
556 bool in_place,
557 typename ViewTypeIn,
558 typename ViewTypeOut>
561 dim,
562 fe_degree,
563 n_q_points_1d,
564 Number>::co_gradients(const ViewTypeIn in,
565 ViewTypeOut out) const
566 {
568 team_member, co_shape_gradients, in, out);
569 }
570
571
572
573 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
574 template <typename ViewType>
575 DEAL_II_HOST_DEVICE inline void
577 dim,
578 fe_degree,
579 n_q_points_1d,
580 Number>::evaluate_values(ViewType u)
581 {
582 if constexpr (dim == 1)
583 values<0, true, false, true>(u, u);
584 else if constexpr (dim == 2)
585 {
586 values<0, true, false, true>(u, u);
587 team_member.team_barrier();
588 values<1, true, false, true>(u, u);
589 }
590 else if constexpr (dim == 3)
591 {
592 values<0, true, false, true>(u, u);
593 team_member.team_barrier();
594 values<1, true, false, true>(u, u);
595 team_member.team_barrier();
596 values<2, true, false, true>(u, u);
597 }
598 else
599 Kokkos::abort("dim must not exceed 3!");
600 }
601
602
603
604 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
605 template <typename ViewType>
606 DEAL_II_HOST_DEVICE inline void
608 dim,
609 fe_degree,
610 n_q_points_1d,
611 Number>::integrate_values(ViewType u)
612 {
613 if constexpr (dim == 1)
614 values<0, false, false, true>(u, u);
615 else if constexpr (dim == 2)
616 {
617 values<0, false, false, true>(u, u);
618 team_member.team_barrier();
619 values<1, false, false, true>(u, u);
620 }
621 else if constexpr (dim == 3)
622 {
623 values<0, false, false, true>(u, u);
624 team_member.team_barrier();
625 values<1, false, false, true>(u, u);
626 team_member.team_barrier();
627 values<2, false, false, true>(u, u);
628 }
629 else
630 Kokkos::abort("dim must not exceed 3!");
631 }
632
633
634
635 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
636 template <typename ViewTypeIn, typename ViewTypeOut>
637 DEAL_II_HOST_DEVICE inline void
639 dim,
640 fe_degree,
641 n_q_points_1d,
642 Number>::evaluate_gradients(const ViewTypeIn u,
643 ViewTypeOut grad_u)
644 {
645 if constexpr (dim == 1)
646 {
647 gradients<0, true, false, false>(
648 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
649 }
650 else if constexpr (dim == 2)
651 {
652 gradients<0, true, false, false>(
653 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
654 values<0, true, false, false>(
655 u, Kokkos::subview(grad_u, Kokkos::ALL, 1));
656
657 team_member.team_barrier();
658
659 values<1, true, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
660 Kokkos::subview(grad_u, Kokkos::ALL, 0));
661 gradients<1, true, false, true>(
662 Kokkos::subview(grad_u, Kokkos::ALL, 1),
663 Kokkos::subview(grad_u, Kokkos::ALL, 1));
664 }
665 else if constexpr (dim == 3)
666 {
667 gradients<0, true, false, false>(
668 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
669 values<0, true, false, false>(
670 u, Kokkos::subview(grad_u, Kokkos::ALL, 1));
671 values<0, true, false, false>(
672 u, Kokkos::subview(grad_u, Kokkos::ALL, 2));
673
674 team_member.team_barrier();
675
676 values<1, true, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
677 Kokkos::subview(grad_u, Kokkos::ALL, 0));
678 gradients<1, true, false, true>(
679 Kokkos::subview(grad_u, Kokkos::ALL, 1),
680 Kokkos::subview(grad_u, Kokkos::ALL, 1));
681 values<1, true, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 2),
682 Kokkos::subview(grad_u, Kokkos::ALL, 2));
683
684 team_member.team_barrier();
685
686 values<2, true, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
687 Kokkos::subview(grad_u, Kokkos::ALL, 0));
688 values<2, true, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 1),
689 Kokkos::subview(grad_u, Kokkos::ALL, 1));
690 gradients<2, true, false, true>(
691 Kokkos::subview(grad_u, Kokkos::ALL, 2),
692 Kokkos::subview(grad_u, Kokkos::ALL, 2));
693 }
694 else
695 Kokkos::abort("dim must not exceed 3!");
696 }
697
698
699
700 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
701 template <typename ViewType1, typename ViewType2>
702 DEAL_II_HOST_DEVICE inline void
704 dim,
705 fe_degree,
706 n_q_points_1d,
707 Number>::evaluate_values_and_gradients(ViewType1 u,
708 ViewType2
709 grad_u)
710 {
711 if constexpr (dim == 1)
712 {
713 values<0, true, false, true>(u, u);
714 team_member.team_barrier();
715
716 co_gradients<0, true, false, false>(
717 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
718 }
719 else if constexpr (dim == 2)
720 {
721 values<0, true, false, true>(u, u);
722 team_member.team_barrier();
723 values<1, true, false, true>(u, u);
724 team_member.team_barrier();
725
726 co_gradients<0, true, false, false>(
727 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
728 co_gradients<1, true, false, false>(
729 u, Kokkos::subview(grad_u, Kokkos::ALL, 1));
730 }
731 else if constexpr (dim == 3)
732 {
733 values<0, true, false, true>(u, u);
734 team_member.team_barrier();
735 values<1, true, false, true>(u, u);
736 team_member.team_barrier();
737 values<2, true, false, true>(u, u);
738 team_member.team_barrier();
739
740 co_gradients<0, true, false, false>(
741 u, Kokkos::subview(grad_u, Kokkos::ALL, 0));
742 co_gradients<1, true, false, false>(
743 u, Kokkos::subview(grad_u, Kokkos::ALL, 1));
744 co_gradients<2, true, false, false>(
745 u, Kokkos::subview(grad_u, Kokkos::ALL, 2));
746 }
747 else
748 Kokkos::abort("dim must not exceed 3!");
749 }
750
751
752
753 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
754 template <bool add, typename ViewType1, typename ViewType2>
755 DEAL_II_HOST_DEVICE inline void
757 dim,
758 fe_degree,
759 n_q_points_1d,
760 Number>::integrate_gradients(ViewType1 u,
761 ViewType2 grad_u)
762 {
763 if constexpr (dim == 1)
764 {
765 gradients<0, false, add, false>(
766 Kokkos::subview(grad_u, Kokkos::ALL, dim), u);
767 }
768 else if constexpr (dim == 2)
769 {
770 gradients<0, false, false, true>(
771 Kokkos::subview(grad_u, Kokkos::ALL, 0),
772 Kokkos::subview(grad_u, Kokkos::ALL, 0));
773 values<0, false, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 1),
774 Kokkos::subview(grad_u,
775 Kokkos::ALL,
776 1));
777
778 team_member.team_barrier();
779
780 values<1, false, add, false>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
781 u);
782 team_member.team_barrier();
783 gradients<1, false, true, false>(
784 Kokkos::subview(grad_u, Kokkos::ALL, 1), u);
785 }
786 else if constexpr (dim == 3)
787 {
788 gradients<0, false, false, true>(
789 Kokkos::subview(grad_u, Kokkos::ALL, 0),
790 Kokkos::subview(grad_u, Kokkos::ALL, 0));
791 values<0, false, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 1),
792 Kokkos::subview(grad_u,
793 Kokkos::ALL,
794 1));
795 values<0, false, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 2),
796 Kokkos::subview(grad_u,
797 Kokkos::ALL,
798 2));
799
800 team_member.team_barrier();
801
802 values<1, false, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
803 Kokkos::subview(grad_u,
804 Kokkos::ALL,
805 0));
806 gradients<1, false, false, true>(
807 Kokkos::subview(grad_u, Kokkos::ALL, 1),
808 Kokkos::subview(grad_u, Kokkos::ALL, 1));
809 values<1, false, false, true>(Kokkos::subview(grad_u, Kokkos::ALL, 2),
810 Kokkos::subview(grad_u,
811 Kokkos::ALL,
812 2));
813
814 team_member.team_barrier();
815
816 values<2, false, add, false>(Kokkos::subview(grad_u, Kokkos::ALL, 0),
817 u);
818 team_member.team_barrier();
819 values<2, false, true, false>(Kokkos::subview(grad_u, Kokkos::ALL, 1),
820 u);
821 team_member.team_barrier();
822 gradients<2, false, true, false>(
823 Kokkos::subview(grad_u, Kokkos::ALL, 2), u);
824 }
825 else
826 Kokkos::abort("dim must not exceed 3!");
827 }
828
829
830
831 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
832 template <typename ViewType1, typename ViewType2>
833 DEAL_II_HOST_DEVICE inline void
835 dim,
836 fe_degree,
837 n_q_points_1d,
838 Number>::integrate_values_and_gradients(ViewType1 u,
839 ViewType2
840 grad_u)
841 {
842 if constexpr (dim == 1)
843 {
844 co_gradients<0, false, true, false>(
845 Kokkos::subview(grad_u, Kokkos::ALL, 0), u);
846 team_member.team_barrier();
847
848 values<0, false, false, true>(u, u);
849 }
850 else if constexpr (dim == 2)
851 {
852 co_gradients<1, false, true, false>(
853 Kokkos::subview(grad_u, Kokkos::ALL, 1), u);
854 team_member.team_barrier();
855 co_gradients<0, false, true, false>(
856 Kokkos::subview(grad_u, Kokkos::ALL, 0), u);
857 team_member.team_barrier();
858
859 values<1, false, false, true>(u, u);
860 team_member.team_barrier();
861 values<0, false, false, true>(u, u);
862 team_member.team_barrier();
863 }
864 else if constexpr (dim == 3)
865 {
866 co_gradients<2, false, true, false>(
867 Kokkos::subview(grad_u, Kokkos::ALL, 2), u);
868 team_member.team_barrier();
869 co_gradients<1, false, true, false>(
870 Kokkos::subview(grad_u, Kokkos::ALL, 1), u);
871 team_member.team_barrier();
872 co_gradients<0, false, true, false>(
873 Kokkos::subview(grad_u, Kokkos::ALL, 0), u);
874 team_member.team_barrier();
875
876 values<2, false, false, true>(u, u);
877 team_member.team_barrier();
878 values<1, false, false, true>(u, u);
879 team_member.team_barrier();
880 values<0, false, false, true>(u, u);
881 team_member.team_barrier();
882 }
883 else
884 Kokkos::abort("dim must not exceed 3!");
885 }
886 } // namespace internal
887} // namespace Portable
888
890
891#endif
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
void apply(const Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type &team_member, const Kokkos::View< Number *, MemorySpace::Default::kokkos_space > shape_data, const ViewTypeIn in, ViewTypeOut out)
T sum(const T &t, const MPI_Comm mpi_communicator)
constexpr T pow(const T base, const int iexp)
Definition utilities.h:966
#define DEAL_II_HOST_DEVICE
Definition numbers.h:34
Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type TeamHandle