Processing math: 0%
 deal.II version GIT relicensing-2968-g5f01c80b02 2025-03-29 13:10:00+00:00
\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\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
evaluation_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_matrix_free_evaluation_kernels_h
17#define dealii_matrix_free_evaluation_kernels_h
18
19#include <deal.II/base/config.h>
20
23
29
30
32
33
34namespace internal
35{
36 // Select evaluator type from element shape function type
37 template <MatrixFreeFunctions::ElementType element, bool is_long>
39 {};
40
41 template <bool is_long>
42 struct EvaluatorSelector<MatrixFreeFunctions::tensor_general, is_long>
43 {
44 static const EvaluatorVariant variant = evaluate_general;
45 };
46
47 template <>
48 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric, false>
49 {
50 static const EvaluatorVariant variant = evaluate_symmetric;
51 };
52
53 template <>
54 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric, true>
55 {
56 static const EvaluatorVariant variant = evaluate_evenodd;
57 };
58
59 template <bool is_long>
60 struct EvaluatorSelector<MatrixFreeFunctions::truncated_tensor, is_long>
61 {
62 static const EvaluatorVariant variant = evaluate_general;
63 };
64
65 template <>
66 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric_plus_dg0,
67 false>
68 {
69 static const EvaluatorVariant variant = evaluate_general;
70 };
71
72 template <>
73 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric_plus_dg0, true>
74 {
75 static const EvaluatorVariant variant = evaluate_evenodd;
76 };
77
78 template <bool is_long>
79 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric_collocation,
80 is_long>
81 {
82 static const EvaluatorVariant variant = evaluate_evenodd;
83 };
84
85
86
109 int dim,
110 int fe_degree,
111 int n_q_points_1d,
112 typename Number>
114 {
116 EvaluatorSelector<type, (fe_degree + n_q_points_1d > 4)>::variant;
117 using Number2 =
119
121 dim,
122 fe_degree + 1,
123 n_q_points_1d,
124 Number,
125 Number2>;
126
127 static void
128 evaluate(const unsigned int n_components,
129 const EvaluationFlags::EvaluationFlags evaluation_flag,
130 const Number *values_dofs_actual,
132
133 static void
134 integrate(const unsigned int n_components,
135 const EvaluationFlags::EvaluationFlags integration_flag,
136 Number *values_dofs_actual,
138 const bool add_into_values_array);
139
140 static Eval
143 *univariate_shape_data)
144 {
146 return Eval(univariate_shape_data->shape_values_eo,
147 univariate_shape_data->shape_gradients_eo,
148 univariate_shape_data->shape_hessians_eo,
149 univariate_shape_data->fe_degree + 1,
150 univariate_shape_data->n_q_points_1d);
151 else
152 return Eval(univariate_shape_data->shape_values,
153 univariate_shape_data->shape_gradients,
154 univariate_shape_data->shape_hessians,
155 univariate_shape_data->fe_degree + 1,
156 univariate_shape_data->n_q_points_1d);
157 }
158 };
159
160
161
166 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
167 struct FEEvaluationImpl<MatrixFreeFunctions::tensor_none,
168 dim,
169 fe_degree,
170 n_q_points_1d,
171 Number>
172 {
173 static void
174 evaluate(const unsigned int n_components,
175 const EvaluationFlags::EvaluationFlags evaluation_flag,
176 const Number *values_dofs_actual,
178
179 static void
180 integrate(const unsigned int n_components,
181 const EvaluationFlags::EvaluationFlags integration_flag,
182 Number *values_dofs_actual,
184 const bool add_into_values_array);
185 };
186
187
188
190 int dim,
191 int fe_degree,
192 int n_q_points_1d,
193 typename Number>
194 inline void
196 const unsigned int n_components,
197 const EvaluationFlags::EvaluationFlags evaluation_flag,
198 const Number *values_dofs_actual,
200 {
201 if (evaluation_flag == EvaluationFlags::nothing)
202 return;
203
204 std::array<const MatrixFreeFunctions::UnivariateShapeData<Number2> *, 3>
205 univariate_shape_data;
206
207 const auto &shape_data = fe_eval.get_shape_info().data;
208
209 univariate_shape_data.fill(&shape_data.front());
210
211 if (shape_data.size() == dim)
212 for (int i = 1; i < dim; ++i)
213 univariate_shape_data[i] = &shape_data[i];
214
215 Eval eval0 = create_evaluator_tensor_product(univariate_shape_data[0]);
216 Eval eval1 = create_evaluator_tensor_product(univariate_shape_data[1]);
217 Eval eval2 = create_evaluator_tensor_product(univariate_shape_data[2]);
218
219 const unsigned int temp_size =
220 Eval::n_rows_of_product == numbers::invalid_unsigned_int ?
221 0 :
222 (Eval::n_rows_of_product > Eval::n_columns_of_product ?
223 Eval::n_rows_of_product :
224 Eval::n_columns_of_product);
225 Number *temp1 = fe_eval.get_scratch_data().begin();
226 Number *temp2;
227 if (temp_size == 0)
228 {
229 temp2 = temp1 + std::max(Utilities::fixed_power<dim>(
230 shape_data.front().fe_degree + 1),
231 Utilities::fixed_power<dim>(
232 shape_data.front().n_q_points_1d));
233 }
234 else
235 {
236 temp2 = temp1 + temp_size;
237 }
238
239 const std::size_t n_q_points = temp_size == 0 ?
240 fe_eval.get_shape_info().n_q_points :
241 Eval::n_columns_of_product;
242 const std::size_t dofs_per_comp =
244 Utilities::pow(shape_data.front().fe_degree + 1, dim) :
246 const Number *values_dofs =
248 temp1 + 2 * (std::max<std::size_t>(
250 n_q_points)) :
251 values_dofs_actual;
252
254 embed_truncated_into_full_tensor_product<dim, fe_degree>(
255 n_components,
256 const_cast<Number *>(values_dofs),
257 values_dofs_actual,
258 fe_eval);
259
260 Number *values_quad = fe_eval.begin_values();
261 Number *gradients_quad = fe_eval.begin_gradients();
262
263 switch (dim)
264 {
265 case 1:
266 for (unsigned int c = 0; c < n_components; ++c)
267 {
268 if (evaluation_flag & EvaluationFlags::values)
269 eval0.template values<0, true, false>(values_dofs, values_quad);
270 if (evaluation_flag & EvaluationFlags::gradients)
271 eval0.template gradients<0, true, false>(values_dofs,
272 gradients_quad);
273
274 // advance the next component in 1d array
275 values_dofs += dofs_per_comp;
276 values_quad += n_q_points;
277 gradients_quad += n_q_points;
278 }
279 break;
280
281 case 2:
282 for (unsigned int c = 0; c < n_components; ++c)
283 {
284 // grad x
285 if (evaluation_flag & EvaluationFlags::gradients)
286 {
287 eval0.template gradients<0, true, false>(values_dofs, temp1);
288 eval1.template values<1, true, false, 2>(temp1,
289 gradients_quad);
290 }
291
292 // grad y
293 eval0.template values<0, true, false>(values_dofs, temp1);
294 if (evaluation_flag & EvaluationFlags::gradients)
295 eval1.template gradients<1, true, false, 2>(temp1,
296 gradients_quad + 1);
297
298 // val: can use values applied in x
299 if (evaluation_flag & EvaluationFlags::values)
300 eval1.template values<1, true, false>(temp1, values_quad);
301
302 // advance to the next component in 1d array
303 values_dofs += dofs_per_comp;
304 values_quad += n_q_points;
305 gradients_quad += 2 * n_q_points;
306 }
307 break;
308
309 case 3:
310 for (unsigned int c = 0; c < n_components; ++c)
311 {
312 if (evaluation_flag & EvaluationFlags::gradients)
313 {
314 // grad x
315 eval0.template gradients<0, true, false>(values_dofs, temp1);
316 eval1.template values<1, true, false>(temp1, temp2);
317 eval2.template values<2, true, false, 3>(temp2,
318 gradients_quad);
319 }
320
321 // grad y
322 eval0.template values<0, true, false>(values_dofs, temp1);
323 if (evaluation_flag & EvaluationFlags::gradients)
324 {
325 eval1.template gradients<1, true, false>(temp1, temp2);
326 eval2.template values<2, true, false, 3>(temp2,
327 gradients_quad + 1);
328 }
329
330 // grad z: can use the values applied in x direction stored in
331 // temp1
332 eval1.template values<1, true, false>(temp1, temp2);
333 if (evaluation_flag & EvaluationFlags::gradients)
334 eval2.template gradients<2, true, false, 3>(temp2,
335 gradients_quad + 2);
336
337 // val: can use the values applied in x & y direction stored in
338 // temp2
339 if (evaluation_flag & EvaluationFlags::values)
340 eval2.template values<2, true, false>(temp2, values_quad);
341
342 // advance to the next component in 1d array
343 values_dofs += dofs_per_comp;
344 values_quad += n_q_points;
345 gradients_quad += 3 * n_q_points;
346 }
347 break;
348
349 default:
351 }
352
353 // case additional dof for FE_Q_DG0: add values; gradients and second
354 // derivatives evaluate to zero
356 (evaluation_flag & EvaluationFlags::values))
357 {
358 values_quad -= n_components * n_q_points;
359 values_dofs -= n_components * dofs_per_comp;
360 for (std::size_t c = 0; c < n_components; ++c)
361 for (std::size_t q = 0; q < n_q_points; ++q)
362 values_quad[c * n_q_points + q] +=
363 values_dofs[(c + 1) * dofs_per_comp - 1];
364 }
365 }
366
367
368
370 int dim,
371 int fe_degree,
372 int n_q_points_1d,
373 typename Number>
374 inline void
376 const unsigned int n_components,
377 const EvaluationFlags::EvaluationFlags integration_flag,
378 Number *values_dofs_actual,
380 const bool add_into_values_array)
381 {
382 std::array<const MatrixFreeFunctions::UnivariateShapeData<Number2> *, 3>
383 univariate_shape_data;
384
385 const auto &shape_data = fe_eval.get_shape_info().data;
386 univariate_shape_data.fill(&shape_data.front());
387
388 if (shape_data.size() == dim)
389 for (int i = 1; i < dim; ++i)
390 univariate_shape_data[i] = &shape_data[i];
391
392 Eval eval0 = create_evaluator_tensor_product(univariate_shape_data[0]);
393 Eval eval1 = create_evaluator_tensor_product(univariate_shape_data[1]);
394 Eval eval2 = create_evaluator_tensor_product(univariate_shape_data[2]);
395
396 const unsigned int temp_size =
397 Eval::n_rows_of_product == numbers::invalid_unsigned_int ?
398 0 :
399 (Eval::n_rows_of_product > Eval::n_columns_of_product ?
400 Eval::n_rows_of_product :
401 Eval::n_columns_of_product);
402 Number *temp1 = fe_eval.get_scratch_data().begin();
403 Number *temp2;
404 if (temp_size == 0)
405 {
406 temp2 = temp1 + std::max(Utilities::fixed_power<dim>(
407 shape_data.front().fe_degree + 1),
408 Utilities::fixed_power<dim>(
409 shape_data.front().n_q_points_1d));
410 }
411 else
412 {
413 temp2 = temp1 + temp_size;
414 }
415
416 const std::size_t n_q_points = temp_size == 0 ?
417 fe_eval.get_shape_info().n_q_points :
418 Eval::n_columns_of_product;
419 const unsigned int dofs_per_comp =
421 Utilities::fixed_power<dim>(shape_data.front().fe_degree + 1) :
423
424 // expand dof_values to tensor product for truncated tensor products
425 Number *values_dofs =
427 temp1 + 2 * (std::max<std::size_t>(
429 n_q_points)) :
430 values_dofs_actual;
431
432 Number *values_quad = fe_eval.begin_values();
433 Number *gradients_quad = fe_eval.begin_gradients();
434
435 switch (dim)
436 {
437 case 1:
438 for (unsigned int c = 0; c < n_components; ++c)
439 {
440 if (integration_flag & EvaluationFlags::values)
441 {
442 if (add_into_values_array == false)
443 eval0.template values<0, false, false>(values_quad,
444 values_dofs);
445 else
446 eval0.template values<0, false, true>(values_quad,
447 values_dofs);
448 }
449 if (integration_flag & EvaluationFlags::gradients)
450 {
451 if (integration_flag & EvaluationFlags::values ||
452 add_into_values_array == true)
453 eval0.template gradients<0, false, true>(gradients_quad,
454 values_dofs);
455 else
456 eval0.template gradients<0, false, false>(gradients_quad,
457 values_dofs);
458 }
459
460 // advance to the next component in 1d array
461 values_dofs += dofs_per_comp;
462 values_quad += n_q_points;
463 gradients_quad += n_q_points;
464 }
465 break;
466
467 case 2:
468 for (unsigned int c = 0; c < n_components; ++c)
469 {
470 if ((integration_flag & EvaluationFlags::values) &&
471 !(integration_flag & EvaluationFlags::gradients))
472 {
473 eval1.template values<1, false, false>(values_quad, temp1);
474 if (add_into_values_array == false)
475 eval0.template values<0, false, false>(temp1, values_dofs);
476 else
477 eval0.template values<0, false, true>(temp1, values_dofs);
478 }
479 if (integration_flag & EvaluationFlags::gradients)
480 {
481 eval1.template gradients<1, false, false, 2>(gradients_quad +
482 1,
483 temp1);
484 if (integration_flag & EvaluationFlags::values)
485 eval1.template values<1, false, true>(values_quad, temp1);
486 if (add_into_values_array == false)
487 eval0.template values<0, false, false>(temp1, values_dofs);
488 else
489 eval0.template values<0, false, true>(temp1, values_dofs);
490 eval1.template values<1, false, false, 2>(gradients_quad,
491 temp1);
492 eval0.template gradients<0, false, true>(temp1, values_dofs);
493 }
494
495 // advance to the next component in 1d array
496 values_dofs += dofs_per_comp;
497 values_quad += n_q_points;
498 gradients_quad += 2 * n_q_points;
499 }
500 break;
501
502 case 3:
503 for (unsigned int c = 0; c < n_components; ++c)
504 {
505 if ((integration_flag & EvaluationFlags::values) &&
506 !(integration_flag & EvaluationFlags::gradients))
507 {
508 eval2.template values<2, false, false>(values_quad, temp1);
509 eval1.template values<1, false, false>(temp1, temp2);
510 if (add_into_values_array == false)
511 eval0.template values<0, false, false>(temp2, values_dofs);
512 else
513 eval0.template values<0, false, true>(temp2, values_dofs);
514 }
515 if (integration_flag & EvaluationFlags::gradients)
516 {
517 eval2.template gradients<2, false, false, 3>(gradients_quad +
518 2,
519 temp1);
520 if (integration_flag & EvaluationFlags::values)
521 eval2.template values<2, false, true>(values_quad, temp1);
522 eval1.template values<1, false, false>(temp1, temp2);
523 eval2.template values<2, false, false, 3>(gradients_quad + 1,
524 temp1);
525 eval1.template gradients<1, false, true>(temp1, temp2);
526 if (add_into_values_array == false)
527 eval0.template values<0, false, false>(temp2, values_dofs);
528 else
529 eval0.template values<0, false, true>(temp2, values_dofs);
530 eval2.template values<2, false, false, 3>(gradients_quad,
531 temp1);
532 eval1.template values<1, false, false>(temp1, temp2);
533 eval0.template gradients<0, false, true>(temp2, values_dofs);
534 }
535
536 // advance to the next component in 1d array
537 values_dofs += dofs_per_comp;
538 values_quad += n_q_points;
539 gradients_quad += 3 * n_q_points;
540 }
541 break;
542
543 default:
545 }
546
547 // case FE_Q_DG0: add values, gradients and second derivatives are zero
549 {
550 values_dofs -= n_components * dofs_per_comp - dofs_per_comp + 1;
551 values_quad -= n_components * n_q_points;
552 if (integration_flag & EvaluationFlags::values)
553 for (unsigned int c = 0; c < n_components; ++c)
554 {
555 values_dofs[0] = values_quad[0];
556 for (unsigned int q = 1; q < n_q_points; ++q)
557 values_dofs[0] += values_quad[q];
558 values_dofs += dofs_per_comp;
559 values_quad += n_q_points;
560 }
561 else
562 {
563 for (unsigned int c = 0; c < n_components; ++c)
564 values_dofs[c * dofs_per_comp] = Number();
565 values_dofs += n_components * dofs_per_comp;
566 }
567 }
568
570 truncate_tensor_product_to_complete_degrees<dim, fe_degree>(
571 n_components,
572 values_dofs_actual,
573 values_dofs - dofs_per_comp * n_components,
574 fe_eval);
575 }
576
577
578
579 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
580 inline void
583 dim,
584 fe_degree,
585 n_q_points_1d,
586 Number>::evaluate(const unsigned int n_components,
587 const EvaluationFlags::EvaluationFlags evaluation_flag,
588 const Number *values_dofs_actual,
590 {
591 Assert(!(evaluation_flag & EvaluationFlags::hessians), ExcNotImplemented());
592
593 const std::size_t n_dofs =
595 const std::size_t n_q_points = fe_eval.get_shape_info().n_q_points;
596
597 const auto &shape_data = fe_eval.get_shape_info().data;
598
599 using Number2 =
601
602 if (evaluation_flag & EvaluationFlags::values)
603 {
604 const auto *const shape_values = shape_data.front().shape_values.data();
605 auto *out = fe_eval.begin_values();
606 const auto *in = values_dofs_actual;
607
608 for (unsigned int c = 0; c < n_components; c += 3)
609 {
610 if (c + 1 == n_components)
613 /*transpose_matrix*/ true,
614 /*add*/ false,
615 /*consider_strides*/ false,
616 Number,
617 Number2,
618 /*n_components*/ 1>(
619 shape_values, in, out, n_dofs, n_q_points, 1, 1);
620 else if (c + 2 == n_components)
623 /*transpose_matrix*/ true,
624 /*add*/ false,
625 /*consider_strides*/ false,
626 Number,
627 Number2,
628 /*n_components*/ 2>(
629 shape_values, in, out, n_dofs, n_q_points, 1, 1);
630 else
633 /*transpose_matrix*/ true,
634 /*add*/ false,
635 /*consider_strides*/ false,
636 Number,
637 Number2,
638 /*n_components*/ 3>(
639 shape_values, in, out, n_dofs, n_q_points, 1, 1);
640
641 out += 3 * n_q_points;
642 in += 3 * n_dofs;
643 }
644 }
645
646 if (evaluation_flag & EvaluationFlags::gradients)
647 {
648 const auto *const shape_gradients =
649 shape_data.front().shape_gradients.data();
650 auto *out = fe_eval.begin_gradients();
651 const auto *in = values_dofs_actual;
652
653 for (unsigned int c = 0; c < n_components; c += 3)
654 {
655 if (c + 1 == n_components)
658 /*transpose_matrix*/ true,
659 /*add*/ false,
660 /*consider_strides*/ false,
661 Number,
662 Number2,
663 /*n_components*/ 1>(
664 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
665 else if (c + 2 == n_components)
668 /*transpose_matrix*/ true,
669 /*add*/ false,
670 /*consider_strides*/ false,
671 Number,
672 Number2,
673 /*n_components*/ 2>(
674 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
675 else
678 /*transpose_matrix*/ true,
679 /*add*/ false,
680 /*consider_strides*/ false,
681 Number,
682 Number2,
683 /*n_components*/ 3>(
684 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
685
686 out += 3 * n_q_points * dim;
687 in += 3 * n_dofs;
688 }
689 }
690 }
691
692
693
694 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
695 inline void
698 dim,
699 fe_degree,
700 n_q_points_1d,
701 Number>::integrate(const unsigned int n_components,
702 const EvaluationFlags::EvaluationFlags integration_flag,
703 Number *values_dofs_actual,
705 const bool add_into_values_array)
706 {
707 Assert(!(integration_flag & EvaluationFlags::hessians),
709
710 const std::size_t n_dofs =
712 const std::size_t n_q_points = fe_eval.get_shape_info().n_q_points;
713
714 const auto &shape_data = fe_eval.get_shape_info().data;
715
716 using Number2 =
718
719 if (integration_flag & EvaluationFlags::values)
720 {
721 const auto *const shape_values = shape_data.front().shape_values.data();
722 auto *in = fe_eval.begin_values();
723 auto *out = values_dofs_actual;
724
725 for (unsigned int c = 0; c < n_components; c += 3)
726 {
727 if (add_into_values_array == false)
728 {
729 if (c + 1 == n_components)
732 /*transpose_matrix*/ false,
733 /*add*/ false,
734 /*consider_strides*/ false,
735 Number,
736 Number2,
737 /*n_components*/ 1>(
738 shape_values, in, out, n_dofs, n_q_points, 1, 1);
739 else if (c + 2 == n_components)
742 /*transpose_matrix*/ false,
743 /*add*/ false,
744 /*consider_strides*/ false,
745 Number,
746 Number2,
747 /*n_components*/ 2>(
748 shape_values, in, out, n_dofs, n_q_points, 1, 1);
749 else
752 /*transpose_matrix*/ false,
753 /*add*/ false,
754 /*consider_strides*/ false,
755 Number,
756 Number2,
757 /*n_components*/ 3>(
758 shape_values, in, out, n_dofs, n_q_points, 1, 1);
759 }
760 else
761 {
762 if (c + 1 == n_components)
765 /*transpose_matrix*/ false,
766 /*add*/ true,
767 /*consider_strides*/ false,
768 Number,
769 Number2,
770 /*n_components*/ 1>(
771 shape_values, in, out, n_dofs, n_q_points, 1, 1);
772 else if (c + 2 == n_components)
775 /*transpose_matrix*/ false,
776 /*add*/ true,
777 /*consider_strides*/ false,
778 Number,
779 Number2,
780 /*n_components*/ 2>(
781 shape_values, in, out, n_dofs, n_q_points, 1, 1);
782 else
785 /*transpose_matrix*/ false,
786 /*add*/ true,
787 /*consider_strides*/ false,
788 Number,
789 Number2,
790 /*n_components*/ 3>(
791 shape_values, in, out, n_dofs, n_q_points, 1, 1);
792 }
793 out += 3 * n_dofs;
794 in += 3 * n_q_points;
795 }
796 }
797
798 if (integration_flag & EvaluationFlags::gradients)
799 {
800 const auto *const shape_gradients =
801 shape_data.front().shape_gradients.data();
802 auto *in = fe_eval.begin_gradients();
803 auto *out = values_dofs_actual;
804
805 for (unsigned int c = 0; c < n_components; c += 3)
806 {
807 if (add_into_values_array == false &&
808 !(integration_flag & EvaluationFlags::values))
809 {
810 if (c + 1 == n_components)
813 /*transpose_matrix*/ false,
814 /*add*/ false,
815 /*consider_strides*/ false,
816 Number,
817 Number2,
818 /*n_components*/ 1>(
819 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
820 else if (c + 2 == n_components)
823 /*transpose_matrix*/ false,
824 /*add*/ false,
825 /*consider_strides*/ false,
826 Number,
827 Number2,
828 /*n_components*/ 2>(
829 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
830 else
833 /*transpose_matrix*/ false,
834 /*add*/ false,
835 /*consider_strides*/ false,
836 Number,
837 Number2,
838 /*n_components*/ 3>(
839 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
840 }
841 else
842 {
843 if (c + 1 == n_components)
846 /*transpose_matrix*/ false,
847 /*add*/ true,
848 /*consider_strides*/ false,
849 Number,
850 Number2,
851 /*n_components*/ 1>(
852 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
853 else if (c + 2 == n_components)
856 /*transpose_matrix*/ false,
857 /*add*/ true,
858 /*consider_strides*/ false,
859 Number,
860 Number2,
861 /*n_components*/ 2>(
862 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
863 else
866 /*transpose_matrix*/ false,
867 /*add*/ true,
868 /*consider_strides*/ false,
869 Number,
870 Number2,
871 /*n_components*/ 3>(
872 shape_gradients, in, out, n_dofs, n_q_points * dim, 1, 1);
873 }
874 out += 3 * n_dofs;
875 in += 3 * n_q_points * dim;
876 }
877 }
878 }
879
880
881
891 template <EvaluatorVariant variant,
892 EvaluatorQuantity quantity,
893 int dim,
894 int basis_size_1,
895 int basis_size_2>
897 {
898 static_assert(basis_size_1 == 0 || basis_size_1 <= basis_size_2,
899 "The second dimension must not be smaller than the first");
900
923 template <typename Number, typename Number2>
924#ifndef DEBUG
926#endif
927 static void
928 do_forward(const unsigned int n_components,
929 const AlignedVector<Number2> &transformation_matrix,
930 const Number *values_in,
931 Number *values_out,
932 const unsigned int basis_size_1_variable =
934 const unsigned int basis_size_2_variable =
936 {
937 Assert(
938 basis_size_1 != 0 || basis_size_1_variable <= basis_size_2_variable,
939 ExcMessage("The second dimension must not be smaller than the first"));
940
942
943 // we do recursion until dim==1 or dim==2 and we have
944 // basis_size_1==basis_size_2. The latter optimization increases
945 // optimization possibilities for the compiler but does only work for
946 // aliased pointers if the sizes are equal.
947 constexpr int next_dim = (dim == 1 || (dim == 2 && basis_size_1 > 0 &&
948 basis_size_1 == basis_size_2)) ?
949 dim :
950 dim - 1;
951
953 dim,
954 basis_size_1,
955 (basis_size_1 == 0 ? 0 : basis_size_2),
956 Number,
957 Number2>
958 eval_val(transformation_matrix,
959 {},
960 {},
961 basis_size_1_variable,
962 basis_size_2_variable);
963 const unsigned int np_1 =
964 basis_size_1 > 0 ? basis_size_1 : basis_size_1_variable;
965 const unsigned int np_2 =
966 basis_size_1 > 0 ? basis_size_2 : basis_size_2_variable;
967 Assert(np_1 > 0 && np_1 != numbers::invalid_unsigned_int,
968 ExcMessage("Cannot transform with 0-point basis"));
969 Assert(np_2 > 0 && np_2 != numbers::invalid_unsigned_int,
970 ExcMessage("Cannot transform with 0-point basis"));
971
972 // run loop backwards to ensure correctness if values_in aliases with
973 // values_out in case with basis_size_1 < basis_size_2
974 values_in = values_in + n_components * Utilities::fixed_power<dim>(np_1);
975 values_out =
976 values_out + n_components * Utilities::fixed_power<dim>(np_2);
977 for (unsigned int c = n_components; c != 0; --c)
978 {
979 values_in -= Utilities::fixed_power<dim>(np_1);
980 values_out -= Utilities::fixed_power<dim>(np_2);
981 if (next_dim < dim)
982 for (unsigned int q = np_1; q != 0; --q)
984 quantity,
985 next_dim,
986 basis_size_1,
987 basis_size_2>::
988 do_forward(1,
989 transformation_matrix,
990 values_in +
991 (q - 1) * Utilities::fixed_power<next_dim>(np_1),
992 values_out +
993 (q - 1) * Utilities::fixed_power<next_dim>(np_2),
994 basis_size_1_variable,
995 basis_size_2_variable);
996
997 // the recursion stops if dim==1 or if dim==2 and
998 // basis_size_1==basis_size_2 (the latter is used because the
999 // compiler generates nicer code)
1000 if (basis_size_1 > 0 && basis_size_2 == basis_size_1 && dim == 2)
1001 {
1002 eval_val.template values<0, true, false>(values_in, values_out);
1003 eval_val.template values<1, true, false>(values_out, values_out);
1004 }
1005 else if (dim == 1)
1006 eval_val.template values<dim - 1, true, false>(values_in,
1007 values_out);
1008 else
1009 eval_val.template values<dim - 1, true, false>(values_out,
1010 values_out);
1011 }
1012 }
1013
1044 template <typename Number, typename Number2>
1045#ifndef DEBUG
1047#endif
1048 static void
1049 do_backward(const unsigned int n_components,
1050 const AlignedVector<Number2> &transformation_matrix,
1051 const bool add_into_result,
1052 Number *values_in,
1053 Number *values_out,
1054 const unsigned int basis_size_1_variable =
1056 const unsigned int basis_size_2_variable =
1058 {
1059 Assert(
1060 basis_size_1 != 0 || basis_size_1_variable <= basis_size_2_variable,
1061 ExcMessage("The second dimension must not be smaller than the first"));
1062 Assert(add_into_result == false || values_in != values_out,
1063 ExcMessage(
1064 "Input and output cannot alias with each other when "
1065 "adding the result of the basis change to existing data"));
1066
1067 Assert(quantity == EvaluatorQuantity::value ||
1068 quantity == EvaluatorQuantity::hessian,
1070
1071 constexpr int next_dim =
1072 (dim > 2 ||
1073 ((basis_size_1 == 0 || basis_size_2 > basis_size_1) && dim > 1)) ?
1074 dim - 1 :
1075 dim;
1076 EvaluatorTensorProduct<variant,
1077 dim,
1078 basis_size_1,
1079 (basis_size_1 == 0 ? 0 : basis_size_2),
1080 Number,
1081 Number2>
1082 eval_val(transformation_matrix,
1083 transformation_matrix,
1084 transformation_matrix,
1085 basis_size_1_variable,
1086 basis_size_2_variable);
1087 const unsigned int np_1 =
1088 basis_size_1 > 0 ? basis_size_1 : basis_size_1_variable;
1089 const unsigned int np_2 =
1090 basis_size_1 > 0 ? basis_size_2 : basis_size_2_variable;
1091 Assert(np_1 > 0 && np_1 != numbers::invalid_unsigned_int,
1092 ExcMessage("Cannot transform with 0-point basis"));
1093 Assert(np_2 > 0 && np_2 != numbers::invalid_unsigned_int,
1094 ExcMessage("Cannot transform with 0-point basis"));
1095
1096 for (unsigned int c = 0; c < n_components; ++c)
1097 {
1098 if (basis_size_1 > 0 && basis_size_2 == basis_size_1 && dim == 2)
1099 {
1100 if (quantity == EvaluatorQuantity::value)
1101 eval_val.template values<1, false, false>(values_in, values_in);
1102 else
1103 eval_val.template hessians<1, false, false>(values_in,
1104 values_in);
1105
1106 if (add_into_result)
1107 {
1108 if (quantity == EvaluatorQuantity::value)
1109 eval_val.template values<0, false, true>(values_in,
1110 values_out);
1111 else
1112 eval_val.template hessians<0, false, true>(values_in,
1113 values_out);
1114 }
1115 else
1116 {
1117 if (quantity == EvaluatorQuantity::value)
1118 eval_val.template values<0, false, false>(values_in,
1119 values_out);
1120 else
1121 eval_val.template hessians<0, false, false>(values_in,
1122 values_out);
1123 }
1124 }
1125 else
1126 {
1127 if (dim == 1 && add_into_result)
1128 {
1129 if (quantity == EvaluatorQuantity::value)
1130 eval_val.template values<0, false, true>(values_in,
1131 values_out);
1132 else
1133 eval_val.template hessians<0, false, true>(values_in,
1134 values_out);
1135 }
1136 else if (dim == 1)
1137 {
1138 if (quantity == EvaluatorQuantity::value)
1139 eval_val.template values<0, false, false>(values_in,
1140 values_out);
1141 else
1142 eval_val.template hessians<0, false, false>(values_in,
1143 values_out);
1144 }
1145 else
1146 {
1147 if (quantity == EvaluatorQuantity::value)
1148 eval_val.template values<dim - 1, false, false>(values_in,
1149 values_in);
1150 else
1151 eval_val.template hessians<dim - 1, false, false>(
1152 values_in, values_in);
1153 }
1154 }
1155 if (next_dim < dim)
1156 for (unsigned int q = 0; q < np_1; ++q)
1158 quantity,
1159 next_dim,
1160 basis_size_1,
1161 basis_size_2>::
1162 do_backward(1,
1163 transformation_matrix,
1164 add_into_result,
1165 values_in +
1166 q * Utilities::fixed_power<next_dim>(np_2),
1167 values_out +
1168 q * Utilities::fixed_power<next_dim>(np_1),
1169 basis_size_1_variable,
1170 basis_size_2_variable);
1171
1172 values_in += Utilities::fixed_power<dim>(np_2);
1173 values_out += Utilities::fixed_power<dim>(np_1);
1174 }
1175 }
1176
1197 template <typename Number, typename Number2>
1198 static void
1199 do_mass(const unsigned int n_components,
1200 const AlignedVector<Number2> &transformation_matrix,
1201 const AlignedVector<Number> &coefficients,
1202 const Number *values_in,
1203 Number *scratch_data,
1204 Number *values_out)
1205 {
1206 constexpr int next_dim = dim > 1 ? dim - 1 : dim;
1207 Number *my_scratch =
1208 basis_size_1 != basis_size_2 ? scratch_data : values_out;
1209
1210 const unsigned int size_per_component = Utilities::pow(basis_size_2, dim);
1211 Assert(coefficients.size() == size_per_component ||
1212 coefficients.size() == n_components * size_per_component,
1213 ExcDimensionMismatch(coefficients.size(), size_per_component));
1214 const unsigned int stride =
1215 coefficients.size() == size_per_component ? 0 : 1;
1216
1217 for (unsigned int q = basis_size_1; q != 0; --q)
1219 variant,
1221 next_dim,
1222 basis_size_1,
1223 basis_size_2>::do_forward(n_components,
1224 transformation_matrix,
1225 values_in +
1226 (q - 1) *
1227 Utilities::pow(basis_size_1, dim - 1),
1228 my_scratch +
1229 (q - 1) *
1230 Utilities::pow(basis_size_2, dim - 1));
1231 EvaluatorTensorProduct<variant,
1232 dim,
1233 basis_size_1,
1234 basis_size_2,
1235 Number,
1236 Number2>
1237 eval_val(transformation_matrix);
1238 const unsigned int n_inner_blocks =
1239 (dim > 1 && basis_size_2 < 10) ? basis_size_2 : 1;
1240 const unsigned int n_blocks = Utilities::pow(basis_size_2, dim - 1);
1241 for (unsigned int ii = 0; ii < n_blocks; ii += n_inner_blocks)
1242 for (unsigned int c = 0; c < n_components; ++c)
1243 {
1244 for (unsigned int i = ii; i < ii + n_inner_blocks; ++i)
1245 eval_val.template values_one_line<dim - 1, true, false>(
1246 my_scratch + i, my_scratch + i);
1247 for (unsigned int q = 0; q < basis_size_2; ++q)
1248 for (unsigned int i = ii; i < ii + n_inner_blocks; ++i)
1249 my_scratch[i + q * n_blocks + c * size_per_component] *=
1250 coefficients[i + q * n_blocks +
1251 c * stride * size_per_component];
1252 for (unsigned int i = ii; i < ii + n_inner_blocks; ++i)
1253 eval_val.template values_one_line<dim - 1, false, false>(
1254 my_scratch + i, my_scratch + i);
1255 }
1256 for (unsigned int q = 0; q < basis_size_1; ++q)
1259 next_dim,
1260 basis_size_1,
1261 basis_size_2>::
1262 do_backward(n_components,
1263 transformation_matrix,
1264 false,
1265 my_scratch + q * Utilities::pow(basis_size_2, dim - 1),
1266 values_out + q * Utilities::pow(basis_size_1, dim - 1));
1267 }
1268 };
1269
1270
1271
1279 template <int n_points_1d, int dim, typename Number, typename Number2>
1280 inline void
1283 const Number *values,
1284 Number *gradients)
1285 {
1287 (n_points_1d + 1) / 2 * n_points_1d);
1288
1290 dim,
1291 n_points_1d,
1292 n_points_1d,
1293 Number,
1294 Number2>
1295 eval({}, shape.shape_gradients_collocation_eo, {});
1297 2,
1298 n_points_1d,
1299 n_points_1d,
1300 Number,
1301 Number2>
1302 eval_2d({}, shape.shape_gradients_collocation_eo, {});
1303
1304 if (dim == 1)
1305 eval.template gradients<0, true, false>(values, gradients);
1306 else
1307 {
1308 if (dim > 2)
1309 eval.template gradients<2, true, false, dim>(values, gradients + 2);
1310 constexpr unsigned int loop_bound = (dim > 2 ? n_points_1d : 1);
1311 constexpr unsigned int n_points_2d = n_points_1d * n_points_1d;
1312 const Number *in = values + (loop_bound - 1) * n_points_2d;
1313 Number *out = gradients + (loop_bound - 1) * dim * n_points_2d;
1314 for (unsigned int l = 0; l < loop_bound; ++l)
1315 {
1316 eval_2d.template gradients<0, true, false, dim>(in, out);
1317 eval_2d.template gradients<1, true, false, dim>(in, out + 1);
1318 in -= n_points_2d;
1319 out -= dim * n_points_2d;
1320 }
1321 }
1322 }
1323
1324
1325
1333 template <int n_points_1d, int dim, typename Number, typename Number2>
1334 inline void
1337 Number *values,
1338 const Number *gradients,
1339 const bool add_into_values_array)
1340 {
1342 (n_points_1d + 1) / 2 * n_points_1d);
1343
1345 dim,
1346 n_points_1d,
1347 n_points_1d,
1348 Number,
1349 Number2>
1350 eval({}, shape.shape_gradients_collocation_eo, {});
1352 2,
1353 n_points_1d,
1354 n_points_1d,
1355 Number,
1356 Number2>
1357 eval_2d({}, shape.shape_gradients_collocation_eo, {});
1358
1359 if (dim == 1)
1360 {
1361 if (add_into_values_array)
1362 eval.template gradients<0, false, true>(gradients, values);
1363 else
1364 eval.template gradients<0, false, false>(gradients, values);
1365 }
1366 else
1367 {
1368 constexpr unsigned int loop_bound = (dim > 2 ? n_points_1d : 1);
1369 constexpr unsigned int n_points_2d = n_points_1d * n_points_1d;
1370
1371 const Number *in = gradients + (loop_bound - 1) * dim * n_points_2d;
1372 Number *out = values + (loop_bound - 1) * n_points_2d;
1373 for (unsigned int l = 0; l < loop_bound; ++l)
1374 {
1375 if (add_into_values_array)
1376 eval_2d.template gradients<0, false, true, dim>(in, out);
1377 else
1378 eval_2d.template gradients<0, false, false, dim>(in, out);
1379 eval_2d.template gradients<1, false, true, dim>(in + 1, out);
1380 in -= dim * n_points_2d;
1381 out -= n_points_2d;
1382 }
1383 }
1384 if (dim > 2)
1385 eval.template gradients<2, false, true, dim>(gradients + 2, values);
1386 }
1387
1388
1389
1396 template <int n_points_1d, int dim, typename Number>
1397 inline void
1398 evaluate_hessians_collocation(const unsigned int n_components,
1400 {
1401 using Number2 =
1403
1404 // might have non-symmetric quadrature formula, so use the more
1405 // conservative 'evaluate_general' scheme rather than 'even_odd' as the
1406 // Hessians are not used very often
1408 fe_eval.get_shape_info().data[0];
1409 AssertDimension(data.shape_gradients_collocation.size(),
1410 data.n_q_points_1d * data.n_q_points_1d);
1412 dim,
1413 n_points_1d,
1414 n_points_1d,
1415 Number,
1416 Number2>
1417 eval({},
1418 data.shape_gradients_collocation.data(),
1419 data.shape_hessians_collocation.data(),
1420 data.n_q_points_1d,
1421 data.n_q_points_1d);
1422
1423 const Number *values = fe_eval.begin_values();
1424 Number *hessians = fe_eval.begin_hessians();
1425 Number *scratch = fe_eval.get_scratch_data().begin();
1426 const std::size_t n_points = fe_eval.get_shape_info().n_q_points;
1427 for (unsigned int comp = 0; comp < n_components; ++comp)
1428 {
1429 // xx derivative
1430 eval.template hessians<0, true, false>(values, hessians);
1431 if (dim > 1)
1432 {
1433 // xy derivative: we might or might not have the gradients already
1434 // computed elsewhere, but we recompute them here since it adds
1435 // only moderate extra work (at most 25%)
1436 eval.template gradients<0, true, false>(values, scratch);
1437 eval.template gradients<1, true, false>(scratch,
1438 hessians + dim * n_points);
1439 // yy derivative
1440 eval.template hessians<1, true, false>(values, hessians + n_points);
1441 }
1442 if (dim > 2)
1443 {
1444 // xz derivative
1445 eval.template gradients<2, true, false>(scratch,
1446 hessians + 4 * n_points);
1447 // yz derivative
1448 eval.template gradients<1, true, false>(values, scratch);
1449 eval.template gradients<2, true, false>(scratch,
1450 hessians + 5 * n_points);
1451 // zz derivative
1452 eval.template hessians<2, true, false>(values,
1453 hessians + 2 * n_points);
1454 }
1455
1456 values += n_points;
1457 hessians += (dim * (dim + 1)) / 2 * n_points;
1458 }
1459 }
1460
1461
1462
1469 template <int n_q_points_1d, int dim, typename Number>
1470 inline void
1471 integrate_hessians_collocation(const unsigned int n_components,
1473 const bool add_into_values_array)
1474 {
1475 using Number2 =
1477
1479 fe_eval.get_shape_info().data[0];
1480 AssertDimension(data.shape_gradients_collocation.size(),
1481 data.n_q_points_1d * data.n_q_points_1d);
1483 dim,
1484 n_q_points_1d,
1485 n_q_points_1d,
1486 Number,
1487 Number2>
1488 eval({},
1489 data.shape_gradients_collocation.data(),
1490 data.shape_hessians_collocation.data(),
1491 data.n_q_points_1d,
1492 data.n_q_points_1d);
1493 Number *values = fe_eval.begin_values();
1494 const Number *hessians = fe_eval.begin_hessians();
1495 Number *scratch = fe_eval.get_scratch_data().begin();
1496 const std::size_t n_points = fe_eval.get_shape_info().n_q_points;
1497
1498 for (unsigned int comp = 0; comp < n_components; ++comp)
1499 {
1500 // xx derivative
1501 if (add_into_values_array == true)
1502 eval.template hessians<0, false, true>(hessians, values);
1503 else
1504 eval.template hessians<0, false, false>(hessians, values);
1505
1506 // yy derivative
1507 if (dim > 1)
1508 eval.template hessians<1, false, true>(hessians + n_points, values);
1509 if (dim > 2)
1510 {
1511 // zz derivative
1512 eval.template hessians<2, false, true>(hessians + 2 * n_points,
1513 values);
1514 // yz derivative
1515 eval.template gradients<2, false, false>(hessians + 5 * n_points,
1516 scratch);
1517 eval.template gradients<1, false, true>(scratch, values);
1518
1519 // xz derivative
1520 eval.template gradients<2, false, false>(hessians + 4 * n_points,
1521 scratch);
1522 }
1523
1524 if (dim > 1)
1525 {
1526 // xy derivative, combined with xz in 3d
1527 eval.template gradients<1, false, (dim > 2)>(hessians +
1528 dim * n_points,
1529 scratch);
1530 eval.template gradients<0, false, true>(scratch, values);
1531 }
1532
1533 values += n_points;
1534 hessians += (dim * (dim + 1)) / 2 * n_points;
1535 }
1536 }
1537
1538
1539
1546 template <int dim, typename Number>
1547 void
1548 evaluate_hessians_slow(const unsigned int n_components,
1549 const Number *values_dofs,
1551 {
1552 const auto &univariate_shape_data = fe_eval.get_shape_info().data;
1553 using Impl =
1555 using Eval = typename Impl::Eval;
1556 Eval eval0 =
1557 Impl::create_evaluator_tensor_product(&univariate_shape_data[0]);
1558 Eval eval1 = Impl::create_evaluator_tensor_product(
1559 &univariate_shape_data[std::min<int>(1,
1560 univariate_shape_data.size() - 1)]);
1561 Eval eval2 = Impl::create_evaluator_tensor_product(
1562 &univariate_shape_data[std::min<int>(2,
1563 univariate_shape_data.size() - 1)]);
1564
1565 const unsigned int n_points = fe_eval.get_shape_info().n_q_points;
1566 Number *tmp1 = fe_eval.get_scratch_data().begin();
1567 Number *tmp2 =
1568 tmp1 + std::max(Utilities::fixed_power<dim>(
1569 univariate_shape_data.front().fe_degree + 1),
1570 Utilities::fixed_power<dim>(
1571 univariate_shape_data.front().n_q_points_1d));
1572 Number *hessians = fe_eval.begin_hessians();
1573
1574 for (unsigned int comp = 0; comp < n_components;
1575 ++comp,
1576 hessians += n_points * dim * (dim + 1) / 2,
1577 values_dofs +=
1579 switch (dim)
1580 {
1581 case 1:
1582 eval0.template hessians<0, true, false>(values_dofs, hessians);
1583 break;
1584 case 2:
1585 // xx derivative
1586 eval0.template hessians<0, true, false>(values_dofs, tmp1);
1587 eval1.template values<1, true, false>(tmp1, hessians);
1588 // xy derivative
1589 eval0.template gradients<0, true, false>(values_dofs, tmp1);
1590 eval1.template gradients<1, true, false>(tmp1,
1591 hessians + 2 * n_points);
1592 // yy derivative
1593 eval0.template values<0, true, false>(values_dofs, tmp1);
1594 eval1.template hessians<1, true, false>(tmp1, hessians + n_points);
1595 break;
1596 case 3:
1597 // xx derivative
1598 eval0.template hessians<0, true, false>(values_dofs, tmp1);
1599 eval1.template values<1, true, false>(tmp1, tmp2);
1600 eval2.template values<2, true, false>(tmp2, hessians);
1601 // xy derivative
1602 eval0.template gradients<0, true, false>(values_dofs, tmp1);
1603 eval1.template gradients<1, true, false>(tmp1, tmp2);
1604 eval2.template values<2, true, false>(tmp2,
1605 hessians + 3 * n_points);
1606 // xz derivative
1607 eval1.template values<1, true, false>(tmp1, tmp2);
1608 eval2.template gradients<2, true, false>(tmp2,
1609 hessians + 4 * n_points);
1610 // yy derivative
1611 eval0.template values<0, true, false>(values_dofs, tmp1);
1612 eval1.template hessians<1, true, false>(tmp1, tmp2);
1613 eval2.template values<2, true, false>(tmp2, hessians + n_points);
1614 // yz derivative
1615 eval1.template gradients<1, true, false>(tmp1, tmp2);
1616 eval2.template gradients<2, true, false>(tmp2,
1617 hessians + 5 * n_points);
1618 // zz derivative
1619 eval1.template values<1, true, false>(tmp1, tmp2);
1620 eval2.template hessians<2, true, false>(tmp2,
1621 hessians + 2 * n_points);
1622 break;
1623
1624 default:
1625 Assert(false,
1627 "Only 1d, 2d and 3d implemented for Hessian"));
1628 }
1629 }
1630
1631
1632
1640 template <int dim, typename Number>
1641 void
1642 integrate_hessians_slow(const unsigned int n_components,
1644 Number *values_dofs,
1645 const bool add_into_values_array)
1646 {
1647 const auto &univariate_shape_data = fe_eval.get_shape_info().data;
1648 using Impl =
1650 using Eval = typename Impl::Eval;
1651 Eval eval0 =
1652 Impl::create_evaluator_tensor_product(&univariate_shape_data[0]);
1653 Eval eval1 = Impl::create_evaluator_tensor_product(
1654 &univariate_shape_data[std::min<int>(1,
1655 univariate_shape_data.size() - 1)]);
1656 Eval eval2 = Impl::create_evaluator_tensor_product(
1657 &univariate_shape_data[std::min<int>(2,
1658 univariate_shape_data.size() - 1)]);
1659
1660 const unsigned int n_points = fe_eval.get_shape_info().n_q_points;
1661 Number *tmp1 = fe_eval.get_scratch_data().begin();
1662 Number *tmp2 =
1663 tmp1 + std::max(Utilities::fixed_power<dim>(
1664 univariate_shape_data.front().fe_degree + 1),
1665 Utilities::fixed_power<dim>(
1666 univariate_shape_data.front().n_q_points_1d));
1667 const Number *hessians = fe_eval.begin_hessians();
1668
1669 for (unsigned int comp = 0; comp < n_components;
1670 ++comp,
1671 hessians += n_points * dim * (dim + 1) / 2,
1672 values_dofs +=
1674 switch (dim)
1675 {
1676 case 1:
1677 if (add_into_values_array)
1678 eval0.template hessians<0, false, true>(hessians, values_dofs);
1679 else
1680 eval0.template hessians<0, false, false>(hessians, values_dofs);
1681 break;
1682 case 2:
1683 // xx derivative
1684 eval1.template values<1, false, false>(hessians, tmp1);
1685 if (add_into_values_array)
1686 eval0.template hessians<0, false, true>(tmp1, values_dofs);
1687 else
1688 eval0.template hessians<0, false, false>(tmp1, values_dofs);
1689
1690 // xy derivative
1691 eval1.template gradients<1, false, false>(hessians + 2 * n_points,
1692 tmp1);
1693 eval0.template gradients<0, false, true>(tmp1, values_dofs);
1694 // yy derivative
1695 eval1.template hessians<1, false, false>(hessians + n_points, tmp1);
1696 eval0.template values<0, false, true>(tmp1, values_dofs);
1697 break;
1698 case 3:
1699 // xx derivative
1700 eval2.template values<2, false, false>(hessians, tmp1);
1701 eval1.template values<1, false, false>(tmp1, tmp2);
1702
1703 if (add_into_values_array)
1704 eval0.template hessians<0, false, true>(tmp2, values_dofs);
1705 else
1706 eval0.template hessians<0, false, false>(tmp2, values_dofs);
1707
1708 // xy derivative
1709 eval2.template values<2, false, false>(hessians + 3 * n_points,
1710 tmp1);
1711 eval1.template gradients<1, false, false>(tmp1, tmp2);
1712 // xz derivative
1713 eval2.template gradients<2, false, false>(hessians + 4 * n_points,
1714 tmp1);
1715 eval1.template values<1, false, true>(tmp1, tmp2);
1716 eval1.template values<0, false, true>(tmp2, values_dofs);
1717
1718 // yy derivative
1719 eval2.template values<2, false, false>(hessians + n_points, tmp1);
1720 eval1.template hessians<1, false, false>(tmp1, tmp2);
1721
1722 // yz derivative
1723 eval2.template gradients<2, false, false>(hessians + 5 * n_points,
1724 tmp1);
1725 eval1.template gradients<1, false, true>(tmp1, tmp2);
1726
1727 // zz derivative
1728 eval2.template hessians<2, false, false>(hessians + 2 * n_points,
1729 tmp1);
1730 eval1.template values<1, false, true>(tmp1, tmp2);
1731 eval0.template values<0, false, true>(tmp2, values_dofs);
1732 break;
1733
1734 default:
1735 Assert(false,
1737 "Only 1d, 2d and 3d implemented for Hessian"));
1738 }
1739 }
1740
1741
1742
1755 template <int dim, int fe_degree, typename Number>
1757 {
1758 using Number2 =
1761 dim,
1762 fe_degree + 1,
1763 fe_degree + 1,
1764 Number,
1765 Number2>;
1766
1767 static void
1768 evaluate(const unsigned int n_components,
1769 const EvaluationFlags::EvaluationFlags evaluation_flag,
1770 const Number *values_dofs,
1772 {
1773 constexpr std::size_t n_points = Utilities::pow(fe_degree + 1, dim);
1774
1775 for (unsigned int c = 0; c < n_components; ++c)
1776 {
1777 if ((evaluation_flag & EvaluationFlags::values) != 0u)
1778 for (unsigned int i = 0; i < n_points; ++i)
1779 fe_eval.begin_values()[n_points * c + i] =
1780 values_dofs[n_points * c + i];
1781
1782 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
1783 evaluate_gradients_collocation<fe_degree + 1, dim>(
1784 fe_eval.get_shape_info().data.front(),
1785 values_dofs + c * n_points,
1786 fe_eval.begin_gradients() + c * dim * n_points);
1787 }
1788 }
1789
1790 static void
1791 integrate(const unsigned int n_components,
1792 const EvaluationFlags::EvaluationFlags integration_flag,
1793 Number *values_dofs,
1795 const bool add_into_values_array)
1796 {
1797 constexpr std::size_t n_points = Utilities::pow(fe_degree + 1, dim);
1798
1799 for (unsigned int c = 0; c < n_components; ++c)
1800 {
1801 if ((integration_flag & EvaluationFlags::values) != 0u)
1802 {
1803 if (add_into_values_array)
1804 for (unsigned int i = 0; i < n_points; ++i)
1805 values_dofs[n_points * c + i] +=
1806 fe_eval.begin_values()[n_points * c + i];
1807 else
1808 for (unsigned int i = 0; i < n_points; ++i)
1809 values_dofs[n_points * c + i] =
1810 fe_eval.begin_values()[n_points * c + i];
1811 }
1812
1813 if ((integration_flag & EvaluationFlags::gradients) != 0u)
1814 integrate_gradients_collocation<fe_degree + 1, dim>(
1815 fe_eval.get_shape_info().data.front(),
1816 values_dofs + c * n_points,
1817 fe_eval.begin_gradients() + c * dim * n_points,
1818 add_into_values_array ||
1819 ((integration_flag & EvaluationFlags::values) != 0u));
1820 }
1821 }
1822 };
1823
1824
1825
1836 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
1838 {
1839 static void
1840 evaluate(const unsigned int n_components,
1841 const EvaluationFlags::EvaluationFlags evaluation_flag,
1842 const Number *values_dofs,
1844 {
1845 const auto &shape_data = fe_eval.get_shape_info().data.front();
1846
1847 Assert(n_q_points_1d > fe_degree,
1848 ExcMessage("You lose information when going to a collocation "
1849 "space of lower degree, so the evaluation results "
1850 "would be wrong. Thus, this class does not permit "
1851 "the chosen operation."));
1852 constexpr std::size_t n_dofs = Utilities::pow(fe_degree + 1, dim);
1853 constexpr std::size_t n_q_points = Utilities::pow(n_q_points_1d, dim);
1854
1855 for (unsigned int c = 0; c < n_components; ++c)
1856 {
1860 dim,
1861 (fe_degree >= n_q_points_1d ? n_q_points_1d : fe_degree + 1),
1862 n_q_points_1d>::do_forward(1,
1863 shape_data.shape_values_eo,
1864 values_dofs + c * n_dofs,
1865 fe_eval.begin_values() + c * n_q_points);
1866
1867 // apply derivatives in the collocation space
1868 if (evaluation_flag & EvaluationFlags::gradients)
1869 evaluate_gradients_collocation<n_q_points_1d, dim>(
1870 shape_data,
1871 fe_eval.begin_values() + c * n_q_points,
1872 fe_eval.begin_gradients() + c * dim * n_q_points);
1873 }
1874 }
1875
1876 static void
1877 integrate(const unsigned int n_components,
1878 const EvaluationFlags::EvaluationFlags integration_flag,
1879 Number *values_dofs,
1881 const bool add_into_values_array)
1882 {
1883 const auto &shape_data = fe_eval.get_shape_info().data.front();
1884
1885 Assert(n_q_points_1d > fe_degree,
1886 ExcMessage("You lose information when going to a collocation "
1887 "space of lower degree, so the evaluation results "
1888 "would be wrong. Thus, this class does not permit "
1889 "the chosen operation."));
1890 constexpr std::size_t n_q_points = Utilities::pow(n_q_points_1d, dim);
1891
1892 for (unsigned int c = 0; c < n_components; ++c)
1893 {
1894 // apply derivatives in collocation space
1895 if (integration_flag & EvaluationFlags::gradients)
1896 integrate_gradients_collocation<n_q_points_1d, dim>(
1897 shape_data,
1898 fe_eval.begin_values() + c * n_q_points,
1899 fe_eval.begin_gradients() + c * dim * n_q_points,
1900 /*add_into_values_array=*/
1901 integration_flag & EvaluationFlags::values);
1902
1903 // transform back to the original space
1907 dim,
1908 (fe_degree >= n_q_points_1d ? n_q_points_1d : fe_degree + 1),
1909 n_q_points_1d>::do_backward(1,
1910 shape_data.shape_values_eo,
1911 add_into_values_array,
1912 fe_eval.begin_values() + c * n_q_points,
1913 values_dofs +
1914 c *
1915 Utilities::pow(fe_degree + 1, dim));
1916 }
1917 }
1918 };
1919
1920
1921
1926 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
1927 struct FEEvaluationImpl<MatrixFreeFunctions::tensor_raviart_thomas,
1928 dim,
1929 fe_degree,
1930 n_q_points_1d,
1931 Number>
1932 {
1933 using Number2 =
1935
1936 template <bool integrate>
1937 static void
1938 evaluate_or_integrate(
1939 const EvaluationFlags::EvaluationFlags evaluation_flag,
1940 Number *values_dofs_actual,
1942 const bool add_into_values_array = false);
1943 };
1944
1945
1946
1947 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
1948 template <bool integrate>
1949 inline void
1951 dim,
1952 fe_degree,
1953 n_q_points_1d,
1954 Number>::
1955 evaluate_or_integrate(
1956 const EvaluationFlags::EvaluationFlags evaluation_flag,
1957 Number *values_dofs,
1959 const bool add)
1960 {
1961 Assert(dim == 2 || dim == 3,
1962 ExcMessage("Only dim = 2,3 implemented for Raviart-Thomas "
1963 "evaluation/integration"));
1964
1965 if (evaluation_flag == EvaluationFlags::nothing)
1966 return;
1967
1968 AssertDimension(fe_eval.get_shape_info().data.size(), 2);
1969 AssertDimension(n_q_points_1d,
1970 fe_eval.get_shape_info().data[0].n_q_points_1d);
1971 AssertDimension(n_q_points_1d,
1972 fe_eval.get_shape_info().data[1].n_q_points_1d);
1973 AssertDimension(fe_degree, fe_eval.get_shape_info().data[0].fe_degree);
1974 AssertDimension(fe_degree, fe_eval.get_shape_info().data[1].fe_degree + 1);
1975
1976 const auto &shape_data = fe_eval.get_shape_info().data;
1977 const unsigned int dofs_per_component =
1978 Utilities::pow(fe_degree, dim - 1) * (fe_degree + 1);
1979 const unsigned int n_points = Utilities::pow(n_q_points_1d, dim);
1980 Number *gradients = fe_eval.begin_gradients();
1981 Number *values = fe_eval.begin_values();
1982
1983 if (integrate)
1984 {
1986 eval;
1987
1988 const bool do_values = evaluation_flag & EvaluationFlags::values;
1989 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
1990 integrate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
1991 values,
1992 gradients,
1993 do_values);
1994 if constexpr (dim > 2)
1995 eval.template tangential<2, 0>(shape_data[1], values, values);
1996 eval.template tangential<1, 0>(shape_data[1], values, values);
1997 eval.template normal<0>(shape_data[0], values, values_dofs, add);
1998
1999 values += n_points;
2000 gradients += n_points * dim;
2001 values_dofs += dofs_per_component;
2002
2003 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
2004 integrate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
2005 values,
2006 gradients,
2007 do_values);
2008 if constexpr (dim > 2)
2009 eval.template tangential<2, 1>(shape_data[1], values, values);
2010 eval.template tangential<0, 1>(shape_data[1], values, values);
2011 eval.template normal<1>(shape_data[0], values, values_dofs, add);
2012
2013 if constexpr (dim > 2)
2014 {
2015 values += n_points;
2016 gradients += n_points * dim;
2017 values_dofs += dofs_per_component;
2018
2019 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
2020 integrate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
2021 values,
2022 gradients,
2023 do_values);
2024 eval.template tangential<1, 2>(shape_data[1], values, values);
2025 eval.template tangential<0, 2>(shape_data[1], values, values);
2026 eval.template normal<2>(shape_data[0], values, values_dofs, add);
2027 }
2028 }
2029 else
2030 {
2032 eval;
2033 eval.template normal<0>(shape_data[0], values_dofs, values);
2034 eval.template tangential<1, 0>(shape_data[1], values, values);
2035 if constexpr (dim > 2)
2036 eval.template tangential<2, 0>(shape_data[1], values, values);
2037 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
2038 evaluate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
2039 values,
2040 gradients);
2041
2042 values += n_points;
2043 gradients += n_points * dim;
2044 values_dofs += dofs_per_component;
2045
2046 eval.template normal<1>(shape_data[0], values_dofs, values);
2047 eval.template tangential<0, 1>(shape_data[1], values, values);
2048 if constexpr (dim > 2)
2049 eval.template tangential<2, 1>(shape_data[1], values, values);
2050 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
2051 evaluate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
2052 values,
2053 gradients);
2054
2055 if constexpr (dim > 2)
2056 {
2057 values += n_points;
2058 gradients += n_points * dim;
2059 values_dofs += dofs_per_component;
2060
2061 eval.template normal<2>(shape_data[0], values_dofs, values);
2062 eval.template tangential<0, 2>(shape_data[1], values, values);
2063 eval.template tangential<1, 2>(shape_data[1], values, values);
2064 if ((evaluation_flag & EvaluationFlags::gradients) != 0u)
2065 evaluate_gradients_collocation<n_q_points_1d, dim>(shape_data[0],
2066 values,
2067 gradients);
2068 }
2069 }
2070 }
2071
2072
2073
2089 template <int dim, typename Number, bool do_integrate>
2091 {
2092 template <int fe_degree, int n_q_points_1d, typename OtherNumber>
2093 static bool
2094 run(const unsigned int n_components,
2095 const EvaluationFlags::EvaluationFlags evaluation_flag,
2096 OtherNumber *values_dofs,
2098 const bool sum_into_values_array_in = false)
2099 {
2100 // `OtherNumber` is either `const Number` (evaluate()) or `Number`
2101 // (integrate())
2102 static_assert(std::is_same_v<Number, std::remove_const_t<OtherNumber>>,
2103 "Type of Number and of OtherNumber do not match.");
2104
2105 const auto element_type = fe_eval.get_shape_info().element_type;
2106 using ElementType = MatrixFreeFunctions::ElementType;
2107
2108 Assert(fe_eval.get_shape_info().data.size() == 1 ||
2109 (fe_eval.get_shape_info().data.size() == dim &&
2110 element_type == ElementType::tensor_general) ||
2111 element_type == ElementType::tensor_raviart_thomas,
2113
2114 EvaluationFlags::EvaluationFlags actual_flag = evaluation_flag;
2115 bool sum_into_values_array = sum_into_values_array_in;
2116 if (evaluation_flag & EvaluationFlags::hessians)
2117 {
2118 actual_flag |= EvaluationFlags::values;
2121 if constexpr (do_integrate)
2122 {
2123 if (fe_eval.get_shape_info().data[0].fe_degree <
2124 fe_eval.get_shape_info().data[0].n_q_points_1d)
2125 integrate_hessians_collocation<n_q_points_1d>(
2126 n_components,
2127 fe_eval,
2128 evaluation_flag & EvaluationFlags::values);
2129 else
2130 {
2131 integrate_hessians_slow(n_components,
2132 fe_eval,
2133 values_dofs,
2134 sum_into_values_array);
2135 sum_into_values_array = true;
2136 }
2137 }
2138 }
2139
2140 if (fe_degree >= 0 && fe_degree + 1 == n_q_points_1d &&
2141 element_type == ElementType::tensor_symmetric_collocation)
2142 {
2145 n_components,
2146 actual_flag,
2147 values_dofs,
2148 fe_eval,
2149 sum_into_values_array);
2150 }
2151 // '<=' on type means tensor_symmetric or tensor_symmetric_hermite, see
2152 // shape_info.h for more details
2153 else if (fe_degree >= 0 &&
2154 use_collocation_evaluation(fe_degree, n_q_points_1d) &&
2155 element_type <= ElementType::tensor_symmetric)
2156 {
2159 fe_degree,
2160 n_q_points_1d,
2161 Number>>(
2162 n_components,
2163 actual_flag,
2164 values_dofs,
2165 fe_eval,
2166 sum_into_values_array);
2167 }
2168 else if (fe_degree >= 0 &&
2169 element_type <= ElementType::tensor_symmetric_no_collocation)
2170 {
2171 evaluate_or_integrate<FEEvaluationImpl<ElementType::tensor_symmetric,
2172 dim,
2173 fe_degree,
2174 n_q_points_1d,
2175 Number>>(
2176 n_components,
2177 actual_flag,
2178 values_dofs,
2179 fe_eval,
2180 sum_into_values_array);
2181 }
2182 else if (element_type == ElementType::tensor_none)
2183 {
2185 FEEvaluationImpl<ElementType::tensor_none, dim, -1, 0, Number>>(
2186 n_components,
2187 actual_flag,
2188 values_dofs,
2189 fe_eval,
2190 sum_into_values_array);
2191 }
2192 else if (element_type == ElementType::tensor_symmetric_plus_dg0)
2193 {
2195 FEEvaluationImpl<ElementType::tensor_symmetric_plus_dg0,
2196 dim,
2197 fe_degree,
2198 n_q_points_1d,
2199 Number>>(n_components,
2200 actual_flag,
2201 values_dofs,
2202 fe_eval,
2203 sum_into_values_array);
2204 }
2205 else if (element_type == ElementType::truncated_tensor)
2206 {
2207 evaluate_or_integrate<FEEvaluationImpl<ElementType::truncated_tensor,
2208 dim,
2209 fe_degree,
2210 n_q_points_1d,
2211 Number>>(
2212 n_components,
2213 actual_flag,
2214 values_dofs,
2215 fe_eval,
2216 sum_into_values_array);
2217 }
2218 else if (element_type == ElementType::tensor_raviart_thomas)
2219 {
2220 if constexpr (fe_degree > 0 && n_q_points_1d > 0 && dim > 1)
2221 {
2222 FEEvaluationImpl<ElementType::tensor_raviart_thomas,
2223 dim,
2224 fe_degree,
2225 n_q_points_1d,
2226 Number>::
2227 template evaluate_or_integrate<do_integrate>(
2228 actual_flag,
2229 const_cast<Number *>(values_dofs),
2230 fe_eval,
2231 sum_into_values_array);
2232 }
2233 else
2234 {
2235 Assert(false,
2236 ExcNotImplemented("Raviart-Thomas currently only possible "
2237 "in 2d/3d and with templated degree"));
2238 }
2239 }
2240 else
2241 {
2242 evaluate_or_integrate<FEEvaluationImpl<ElementType::tensor_general,
2243 dim,
2244 fe_degree,
2245 n_q_points_1d,
2246 Number>>(
2247 n_components,
2248 actual_flag,
2249 values_dofs,
2250 fe_eval,
2251 sum_into_values_array);
2252 }
2253
2254 if ((evaluation_flag & EvaluationFlags::hessians) && !do_integrate)
2255 {
2258 if (fe_eval.get_shape_info().data[0].fe_degree <
2259 fe_eval.get_shape_info().data[0].n_q_points_1d)
2260 evaluate_hessians_collocation<n_q_points_1d>(n_components, fe_eval);
2261 else
2262 evaluate_hessians_slow(n_components, values_dofs, fe_eval);
2263 }
2264
2265 return false;
2266 }
2267
2268 private:
2269 template <typename T>
2270 static void
2272 const unsigned int n_components,
2273 const EvaluationFlags::EvaluationFlags evaluation_flag,
2274 const Number *values_dofs,
2276 const bool sum_into_values_array,
2277 std::bool_constant<false>)
2278 {
2279 (void)sum_into_values_array;
2280
2281 T::evaluate(n_components, evaluation_flag, values_dofs, fe_eval);
2282 }
2283
2284 template <typename T>
2285 static void
2287 const unsigned int n_components,
2288 const EvaluationFlags::EvaluationFlags evaluation_flag,
2289 Number *values_dofs,
2291 const bool sum_into_values_array,
2292 std::bool_constant<true>)
2293 {
2294 T::integrate(n_components,
2295 evaluation_flag,
2296 values_dofs,
2297 fe_eval,
2298 sum_into_values_array);
2299 }
2300
2301 template <typename T, typename OtherNumber>
2302 static void
2304 const unsigned int n_components,
2305 const EvaluationFlags::EvaluationFlags evaluation_flag,
2306 OtherNumber *values_dofs,
2308 const bool sum_into_values_array)
2309 {
2310 evaluate_or_integrate<T>(n_components,
2311 evaluation_flag,
2312 values_dofs,
2313 fe_eval,
2314 sum_into_values_array,
2315 std::bool_constant<do_integrate>());
2316 }
2317 };
2318
2319
2320
2325 template <int dim, typename Number>
2327 {
2328 using Number2 =
2330
2331 template <int fe_degree, int = 0>
2332 static bool
2333 run(const unsigned int n_components,
2335 const Number *in_array,
2336 Number *out_array)
2337 {
2338 const unsigned int given_degree =
2339 (fe_degree > -1) ? fe_degree :
2340 fe_eval.get_shape_info().data.front().fe_degree;
2341
2342 const unsigned int dofs_per_component =
2343 Utilities::pow(given_degree + 1, dim);
2344
2345 Assert(dim >= 1 || dim <= 3, ExcNotImplemented());
2349
2351 dim,
2352 fe_degree + 1,
2353 fe_degree + 1,
2354 Number,
2355 Number2>
2356 evaluator({},
2357 {},
2358 fe_eval.get_shape_info().data.front().inverse_shape_values_eo,
2359 given_degree + 1,
2360 given_degree + 1);
2361
2362 for (unsigned int d = 0; d < n_components; ++d)
2363 {
2364 const Number *in = in_array + d * dofs_per_component;
2365 Number *out = out_array + d * dofs_per_component;
2366 // Need to select 'apply' method with hessian slot because values
2367 // assume symmetries that do not exist in the inverse shapes
2368 evaluator.template hessians<0, true, false>(in, out);
2369 if (dim > 1)
2370 evaluator.template hessians<1, true, false>(out, out);
2371 if (dim > 2)
2372 evaluator.template hessians<2, true, false>(out, out);
2373 }
2374 for (unsigned int q = 0; q < dofs_per_component; ++q)
2375 {
2376 const Number inverse_JxW_q = Number(1.) / fe_eval.JxW(q);
2377 for (unsigned int d = 0; d < n_components; ++d)
2378 out_array[q + d * dofs_per_component] *= inverse_JxW_q;
2379 }
2380 for (unsigned int d = 0; d < n_components; ++d)
2381 {
2382 Number *out = out_array + d * dofs_per_component;
2383 if (dim > 2)
2384 evaluator.template hessians<2, false, false>(out, out);
2385 if (dim > 1)
2386 evaluator.template hessians<1, false, false>(out, out);
2387 evaluator.template hessians<0, false, false>(out, out);
2388 }
2389 return false;
2390 }
2391 };
2392
2393
2394
2401 template <int dim, typename Number>
2403 {
2404 using Number2 =
2406
2407 template <int fe_degree, int = 0>
2408 static bool
2409 run(const unsigned int n_desired_components,
2411 const ArrayView<const Number> &inverse_coefficients,
2412 const bool dyadic_coefficients,
2413 const Number *in_array,
2414 Number *out_array)
2415 {
2416 const unsigned int given_degree =
2417 (fe_degree > -1) ? fe_degree :
2418 fe_eval.get_shape_info().data.front().fe_degree;
2419
2420 const unsigned int dofs_per_component =
2421 Utilities::pow(given_degree + 1, dim);
2422
2423 Assert(inverse_coefficients.size() > 0 &&
2424 inverse_coefficients.size() % dofs_per_component == 0,
2425 ExcMessage(
2426 "Expected diagonal to be a multiple of scalar dof per cells"));
2427
2428 if (!dyadic_coefficients)
2429 {
2430 if (inverse_coefficients.size() != dofs_per_component)
2431 AssertDimension(n_desired_components * dofs_per_component,
2432 inverse_coefficients.size());
2433 }
2434 else
2435 {
2436 AssertDimension(n_desired_components * n_desired_components *
2437 dofs_per_component,
2438 inverse_coefficients.size());
2439 }
2440
2441 Assert(dim >= 1 || dim <= 3, ExcNotImplemented());
2445
2447 dim,
2448 fe_degree + 1,
2449 fe_degree + 1,
2450 Number,
2451 Number2>
2452 evaluator({},
2453 {},
2454 fe_eval.get_shape_info().data.front().inverse_shape_values_eo,
2455 given_degree + 1,
2456 given_degree + 1);
2457
2458 const Number *in = in_array;
2459 Number *out = out_array;
2460
2461 const Number *inv_coefficient = inverse_coefficients.data();
2462
2463 const unsigned int shift_coefficient =
2464 inverse_coefficients.size() > dofs_per_component ? dofs_per_component :
2465 0;
2466
2467 const auto n_comp_outer = dyadic_coefficients ? 1 : n_desired_components;
2468 const auto n_comp_inner = dyadic_coefficients ? n_desired_components : 1;
2469
2470 for (unsigned int d = 0; d < n_comp_outer; ++d)
2471 {
2472 for (unsigned int di = 0; di < n_comp_inner; ++di)
2473 {
2474 const Number *in_ = in + di * dofs_per_component;
2475 Number *out_ = out + di * dofs_per_component;
2476 evaluator.template hessians<0, true, false>(in_, out_);
2477 if (dim > 1)
2478 evaluator.template hessians<1, true, false>(out_, out_);
2479 if (dim > 2)
2480 evaluator.template hessians<2, true, false>(out_, out_);
2481 }
2482 if (dyadic_coefficients)
2483 {
2484 const auto n_coeff_components =
2485 n_desired_components * n_desired_components;
2486 if (n_desired_components == dim)
2487 {
2488 for (unsigned int q = 0; q < dofs_per_component; ++q)
2489 vmult<dim>(&inv_coefficient[q * n_coeff_components],
2490 &in[q],
2491 &out[q],
2492 dofs_per_component);
2493 }
2494 else
2495 {
2496 for (unsigned int q = 0; q < dofs_per_component; ++q)
2497 vmult<-1>(&inv_coefficient[q * n_coeff_components],
2498 &in[q],
2499 &out[q],
2500 dofs_per_component,
2501 n_desired_components);
2502 }
2503 }
2504 else
2505 for (unsigned int q = 0; q < dofs_per_component; ++q)
2506 out[q] *= inv_coefficient[q];
2507
2508 for (unsigned int di = 0; di < n_comp_inner; ++di)
2509 {
2510 Number *out_ = out + di * dofs_per_component;
2511 if (dim > 2)
2512 evaluator.template hessians<2, false, false>(out_, out_);
2513 if (dim > 1)
2514 evaluator.template hessians<1, false, false>(out_, out_);
2515 evaluator.template hessians<0, false, false>(out_, out_);
2516 }
2517
2518 in += dofs_per_component;
2519 out += dofs_per_component;
2520 inv_coefficient += shift_coefficient;
2521 }
2522
2523 return false;
2524 }
2525
2526 private:
2527 template <int n_components>
2528 static inline void
2529 vmult(const Number *inverse_coefficients,
2530 const Number *src,
2531 Number *dst,
2532 const unsigned int dofs_per_component,
2533 const unsigned int n_given_components = 0)
2534 {
2535 const unsigned int n_desired_components =
2536 (n_components > -1) ? n_components : n_given_components;
2537
2538 std::array<Number, dim + 2> tmp = {};
2539 Assert(n_desired_components <= dim + 2,
2540 ExcMessage(
2541 "Number of components larger than dim+2 not supported."));
2542
2543 for (unsigned int d = 0; d < n_desired_components; ++d)
2544 tmp[d] = src[d * dofs_per_component];
2545
2546 for (unsigned int d1 = 0; d1 < n_desired_components; ++d1)
2547 {
2548 const Number *inv_coeff_row =
2549 &inverse_coefficients[d1 * n_desired_components];
2550 Number sum = inv_coeff_row[0] * tmp[0];
2551 for (unsigned int d2 = 1; d2 < n_desired_components; ++d2)
2552 sum += inv_coeff_row[d2] * tmp[d2];
2553 dst[d1 * dofs_per_component] = sum;
2554 }
2555 }
2556 };
2557
2558
2559
2566 template <int dim, typename Number>
2568 {
2569 template <int fe_degree, int n_q_points_1d>
2570 static bool
2571 run(const unsigned int n_desired_components,
2573 const Number *in_array,
2574 Number *out_array)
2575 {
2576 static const bool do_inplace =
2577 fe_degree > -1 && (fe_degree + 1 == n_q_points_1d);
2578
2582
2583 const auto &inverse_shape =
2584 do_inplace ?
2585 fe_eval.get_shape_info().data.front().inverse_shape_values_eo :
2586 fe_eval.get_shape_info().data.front().inverse_shape_values;
2587
2588 const std::size_t dofs_per_component =
2589 do_inplace ? Utilities::pow(fe_degree + 1, dim) :
2591 const std::size_t n_q_points = do_inplace ?
2592 Utilities::pow(fe_degree + 1, dim) :
2593 fe_eval.get_shape_info().n_q_points;
2594
2595 using Number2 =
2598 dim,
2599 fe_degree + 1,
2600 n_q_points_1d,
2601 Number,
2602 Number2>
2603 evaluator({},
2604 {},
2605 inverse_shape,
2606 fe_eval.get_shape_info().data.front().fe_degree + 1,
2607 fe_eval.get_shape_info().data.front().n_q_points_1d);
2608
2609 for (unsigned int d = 0; d < n_desired_components; ++d)
2610 {
2611 const Number *in = in_array + d * n_q_points;
2612 Number *out = out_array + d * dofs_per_component;
2613
2614 auto *temp_1 = do_inplace ? out : fe_eval.get_scratch_data().begin();
2615 auto *temp_2 = do_inplace ?
2616 out :
2617 (temp_1 + std::max(n_q_points, dofs_per_component));
2618
2619 if (dim == 3)
2620 {
2621 evaluator.template hessians<2, false, false>(in, temp_1);
2622 evaluator.template hessians<1, false, false>(temp_1, temp_2);
2623 evaluator.template hessians<0, false, false>(temp_2, out);
2624 }
2625 if (dim == 2)
2626 {
2627 evaluator.template hessians<1, false, false>(in, temp_1);
2628 evaluator.template hessians<0, false, false>(temp_1, out);
2629 }
2630 if (dim == 1)
2631 evaluator.template hessians<0, false, false>(in, out);
2632 }
2633 return false;
2634 }
2635 };
2636} // end of namespace internal
2637
2638
2640
2641#endif
size_type size() const
iterator begin() const
Definition array_view.h:707
value_type * data() const noexcept
Definition array_view.h:666
std::size_t size() const
Definition array_view.h:689
ScalarNumber shape_info_number_type
const ShapeInfoType & get_shape_info() const
Number JxW(const unsigned int q_point) const
const Number * begin_gradients() const
ArrayView< Number > get_scratch_data() const
const Number * begin_values() const
const Number * begin_hessians() const
#define DEAL_II_ALWAYS_INLINE
Definition config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::vector< index_type > data
Definition mpi.cc:746
EvaluationFlags
The EvaluationFlags enum.
constexpr T pow(const T base, const int iexp)
Definition utilities.h:967
void evaluate_hessians_collocation(const unsigned int n_components, FEEvaluationData< dim, Number, false > &fe_eval)
constexpr bool use_collocation_evaluation(const unsigned int fe_degree, const unsigned int n_q_points_1d)
void integrate_gradients_collocation(const MatrixFreeFunctions::UnivariateShapeData< Number2 > &shape, Number *values, const Number *gradients, const bool add_into_values_array)
void evaluate_hessians_slow(const unsigned int n_components, const Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval)
std::enable_if_t<(variant==evaluate_general), void > apply_matrix_vector_product(const Number2 *matrix, const Number *in, Number *out)
void integrate_hessians_collocation(const unsigned int n_components, FEEvaluationData< dim, Number, false > &fe_eval, const bool add_into_values_array)
void evaluate_gradients_collocation(const MatrixFreeFunctions::UnivariateShapeData< Number2 > &shape, const Number *values, Number *gradients)
void integrate_hessians_slow(const unsigned int n_components, const FEEvaluationData< dim, Number, false > &fe_eval, Number *values_dofs, const bool add_into_values_array)
constexpr unsigned int invalid_unsigned_int
Definition types.h:232
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
typename FEEvaluationData< dim, Number, false >::shape_info_number_type Number2
static bool run(const unsigned int n_components, const FEEvaluationData< dim, Number, false > &fe_eval, const Number *in_array, Number *out_array)
typename FEEvaluationData< dim, Number, false >::shape_info_number_type Number2
static bool run(const unsigned int n_desired_components, const FEEvaluationData< dim, Number, false > &fe_eval, const ArrayView< const Number > &inverse_coefficients, const bool dyadic_coefficients, const Number *in_array, Number *out_array)
static void vmult(const Number *inverse_coefficients, const Number *src, Number *dst, const unsigned int dofs_per_component, const unsigned int n_given_components=0)
static bool run(const unsigned int n_desired_components, const FEEvaluationData< dim, Number, false > &fe_eval, const Number *in_array, Number *out_array)
static void do_backward(const unsigned int n_components, const AlignedVector< Number2 > &transformation_matrix, const bool add_into_result, Number *values_in, Number *values_out, const unsigned int basis_size_1_variable=numbers::invalid_unsigned_int, const unsigned int basis_size_2_variable=numbers::invalid_unsigned_int)
static void do_forward(const unsigned int n_components, const AlignedVector< Number2 > &transformation_matrix, const Number *values_in, Number *values_out, const unsigned int basis_size_1_variable=numbers::invalid_unsigned_int, const unsigned int basis_size_2_variable=numbers::invalid_unsigned_int)
static void do_mass(const unsigned int n_components, const AlignedVector< Number2 > &transformation_matrix, const AlignedVector< Number > &coefficients, const Number *values_in, Number *scratch_data, Number *values_out)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval)
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool add_into_values_array)
typename FEEvaluationData< dim, Number, false >::shape_info_number_type Number2
static void evaluate_or_integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, OtherNumber *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool sum_into_values_array)
static void evaluate_or_integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool sum_into_values_array, std::bool_constant< true >)
static bool run(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, OtherNumber *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool sum_into_values_array_in=false)
static void evaluate_or_integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool sum_into_values_array, std::bool_constant< false >)
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval, const bool add_into_values_array)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval)
static const EvaluatorVariant variant
typename FEEvaluationData< dim, Number, false >::shape_info_number_type Number2
EvaluatorTensorProduct< variant, dim, fe_degree+1, n_q_points_1d, Number, Number2 > Eval
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, Number *values_dofs_actual, FEEvaluationData< dim, Number, false > &fe_eval, const bool add_into_values_array)
static Eval create_evaluator_tensor_product(const MatrixFreeFunctions::UnivariateShapeData< Number2 > *univariate_shape_data)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs_actual, FEEvaluationData< dim, Number, false > &fe_eval)
std::vector< UnivariateShapeData< Number > > data
Definition shape_info.h:485