deal.II version 9.7.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
mapping_q_internal.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2025 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#ifndef dealii_mapping_q_internal_h
16#define dealii_mapping_q_internal_h
17
18#include <deal.II/base/config.h>
19
23#include <deal.II/base/point.h>
25#include <deal.II/base/table.h>
26#include <deal.II/base/tensor.h>
27
28#include <deal.II/fe/fe_tools.h>
33
35
42
43#include <algorithm>
44#include <array>
45#include <limits>
46#include <numeric>
47
48
50
51namespace internal
52{
58 namespace MappingQ1
59 {
60 // These are left as templates on the spatial dimension (even though dim
61 // == spacedim must be true for them to make sense) because templates are
62 // expanded before the compiler eliminates code due to the 'if (dim ==
63 // spacedim)' statement (see the body of the general
64 // transform_real_to_unit_cell).
65 template <int spacedim>
66 inline Point<1>
69 &vertices,
70 const Point<spacedim> &p)
71 {
72 Assert(spacedim == 1, ExcInternalError());
73 return Point<1>((p[0] - vertices[0][0]) /
74 (vertices[1][0] - vertices[0][0]));
75 }
76
77
78
79 template <int spacedim>
80 inline Point<2>
83 &vertices,
84 const Point<spacedim> &p)
85 {
86 Assert(spacedim == 2, ExcInternalError());
87
88 // For accuracy reasons, we do all arithmetic in extended precision
89 // (long double). This has a noticeable effect on the hit rate for
90 // borderline cases and thus makes the algorithm more robust.
91 const long double x = p[0];
92 const long double y = p[1];
93
94 const long double x0 = vertices[0][0];
95 const long double x1 = vertices[1][0];
96 const long double x2 = vertices[2][0];
97 const long double x3 = vertices[3][0];
98
99 const long double y0 = vertices[0][1];
100 const long double y1 = vertices[1][1];
101 const long double y2 = vertices[2][1];
102 const long double y3 = vertices[3][1];
103
104 const long double a = (x1 - x3) * (y0 - y2) - (x0 - x2) * (y1 - y3);
105 const long double b = -(x0 - x1 - x2 + x3) * y + (x - 2 * x1 + x3) * y0 -
106 (x - 2 * x0 + x2) * y1 - (x - x1) * y2 +
107 (x - x0) * y3;
108 const long double c = (x0 - x1) * y - (x - x1) * y0 + (x - x0) * y1;
109
110 const long double discriminant = b * b - 4 * a * c;
111 // exit if the point is not in the cell (this is the only case where the
112 // discriminant is negative)
114 discriminant > 0.0,
116
117 long double eta1;
118 long double eta2;
119 const long double sqrt_discriminant = std::sqrt(discriminant);
120 // special case #1: if a is near-zero to make the discriminant exactly
121 // equal b, then use the linear formula
122 if (b != 0.0 && std::abs(b) == sqrt_discriminant)
123 {
124 eta1 = -c / b;
125 eta2 = -c / b;
126 }
127 // special case #2: a is zero for parallelograms and very small for
128 // near-parallelograms:
129 else if (std::abs(a) < 1e-8 * std::abs(b))
130 {
131 // if both a and c are very small then the root should be near
132 // zero: this first case will capture that
133 eta1 = 2 * c / (-b - sqrt_discriminant);
134 eta2 = 2 * c / (-b + sqrt_discriminant);
135 }
136 // finally, use the plain version:
137 else
138 {
139 eta1 = (-b - sqrt_discriminant) / (2 * a);
140 eta2 = (-b + sqrt_discriminant) / (2 * a);
141 }
142 // pick the one closer to the center of the cell.
143 const long double eta =
144 (std::abs(eta1 - 0.5) < std::abs(eta2 - 0.5)) ? eta1 : eta2;
145
146 /*
147 * There are two ways to compute xi from eta, but either one may have a
148 * zero denominator.
149 */
150 const long double subexpr0 = -eta * x2 + x0 * (eta - 1);
151 const long double xi_denominator0 = eta * x3 - x1 * (eta - 1) + subexpr0;
152 const long double max_x =
153 std::max({std::abs(x0), std::abs(x1), std::abs(x2), std::abs(x3)});
154
155 if (std::abs(xi_denominator0) > 1e-10 * max_x)
156 {
157 const double xi = (x + subexpr0) / xi_denominator0;
158 return {xi, static_cast<double>(eta)};
159 }
160 else
161 {
162 const long double max_y =
163 std::max({std::abs(y0), std::abs(y1), std::abs(y2), std::abs(y3)});
164 const long double subexpr1 = -eta * y2 + y0 * (eta - 1);
165 const long double xi_denominator1 =
166 eta * y3 - y1 * (eta - 1) + subexpr1;
167 if (std::abs(xi_denominator1) > 1e-10 * max_y)
168 {
169 const double xi = (subexpr1 + y) / xi_denominator1;
170 return {xi, static_cast<double>(eta)};
171 }
172 else // give up and try Newton iteration
173 {
175 false,
176 (typename Mapping<spacedim,
177 spacedim>::ExcTransformationFailed()));
178 }
179 }
180 // bogus return to placate compiler. It should not be possible to get
181 // here.
183 return {std::numeric_limits<double>::quiet_NaN(),
184 std::numeric_limits<double>::quiet_NaN()};
185 }
186
187
188
189 template <int spacedim>
190 inline Point<3>
193 & /*vertices*/,
194 const Point<spacedim> & /*p*/)
195 {
196 // It should not be possible to get here
198 return {std::numeric_limits<double>::quiet_NaN(),
199 std::numeric_limits<double>::quiet_NaN(),
200 std::numeric_limits<double>::quiet_NaN()};
201 }
202 } // namespace MappingQ1
203
204
205
212 {
217 template <int dim>
218 std::vector<Point<dim>>
219 unit_support_points(const std::vector<Point<1>> &line_support_points,
220 const std::vector<unsigned int> &renumbering)
221 {
222 AssertDimension(Utilities::pow(line_support_points.size(), dim),
223 renumbering.size());
224 std::vector<Point<dim>> points(renumbering.size());
225 const unsigned int n1 = line_support_points.size();
226 for (unsigned int q2 = 0, q = 0; q2 < (dim > 2 ? n1 : 1); ++q2)
227 for (unsigned int q1 = 0; q1 < (dim > 1 ? n1 : 1); ++q1)
228 for (unsigned int q0 = 0; q0 < n1; ++q0, ++q)
229 {
230 points[renumbering[q]][0] = line_support_points[q0][0];
231 if (dim > 1)
232 points[renumbering[q]][1] = line_support_points[q1][0];
233 if (dim > 2)
234 points[renumbering[q]][2] = line_support_points[q2][0];
235 }
236 return points;
237 }
238
239
240
248 inline ::Table<2, double>
249 compute_support_point_weights_on_quad(const unsigned int polynomial_degree)
250 {
251 ::Table<2, double> loqvs;
252
253 // we are asked to compute weights for interior support points, but
254 // there are no interior points if degree==1
255 if (polynomial_degree == 1)
256 return loqvs;
257
258 const unsigned int M = polynomial_degree - 1;
259 const unsigned int n_inner_2d = M * M;
260 const unsigned int n_outer_2d = 4 + 4 * M;
261
262 // set the weights of transfinite interpolation
263 loqvs.reinit(n_inner_2d, n_outer_2d);
264 QGaussLobatto<2> gl(polynomial_degree + 1);
265 for (unsigned int i = 0; i < M; ++i)
266 for (unsigned int j = 0; j < M; ++j)
267 {
268 const Point<2> &p =
269 gl.point((i + 1) * (polynomial_degree + 1) + (j + 1));
270 const unsigned int index_table = i * M + j;
271 for (unsigned int v = 0; v < 4; ++v)
272 loqvs(index_table, v) =
274 loqvs(index_table, 4 + i) = 1. - p[0];
275 loqvs(index_table, 4 + i + M) = p[0];
276 loqvs(index_table, 4 + j + 2 * M) = 1. - p[1];
277 loqvs(index_table, 4 + j + 3 * M) = p[1];
278 }
279
280 // the sum of weights of the points at the outer rim should be one.
281 // check this
282 for (unsigned int unit_point = 0; unit_point < n_inner_2d; ++unit_point)
283 Assert(std::fabs(std::accumulate(loqvs[unit_point].begin(),
284 loqvs[unit_point].end(),
285 0.) -
286 1) < 1e-13 * polynomial_degree,
288
289 return loqvs;
290 }
291
292
293
300 inline ::Table<2, double>
301 compute_support_point_weights_on_hex(const unsigned int polynomial_degree)
302 {
303 ::Table<2, double> lohvs;
304
305 // we are asked to compute weights for interior support points, but
306 // there are no interior points if degree==1
307 if (polynomial_degree == 1)
308 return lohvs;
309
310 const unsigned int M = polynomial_degree - 1;
311
312 const unsigned int n_inner = Utilities::fixed_power<3>(M);
313 const unsigned int n_outer = 8 + 12 * M + 6 * M * M;
314
315 // set the weights of transfinite interpolation
316 lohvs.reinit(n_inner, n_outer);
317 QGaussLobatto<3> gl(polynomial_degree + 1);
318 for (unsigned int i = 0; i < M; ++i)
319 for (unsigned int j = 0; j < M; ++j)
320 for (unsigned int k = 0; k < M; ++k)
321 {
322 const Point<3> &p = gl.point((i + 1) * (M + 2) * (M + 2) +
323 (j + 1) * (M + 2) + (k + 1));
324 const unsigned int index_table = i * M * M + j * M + k;
325
326 // vertices
327 for (unsigned int v = 0; v < 8; ++v)
328 lohvs(index_table, v) =
330
331 // lines
332 {
333 constexpr std::array<unsigned int, 4> line_coordinates_y(
334 {{0, 1, 4, 5}});
335 const Point<2> py(p[0], p[2]);
336 for (unsigned int l = 0; l < 4; ++l)
337 lohvs(index_table, 8 + line_coordinates_y[l] * M + j) =
339 }
340
341 {
342 constexpr std::array<unsigned int, 4> line_coordinates_x(
343 {{2, 3, 6, 7}});
344 const Point<2> px(p[1], p[2]);
345 for (unsigned int l = 0; l < 4; ++l)
346 lohvs(index_table, 8 + line_coordinates_x[l] * M + k) =
348 }
349
350 {
351 constexpr std::array<unsigned int, 4> line_coordinates_z(
352 {{8, 9, 10, 11}});
353 const Point<2> pz(p[0], p[1]);
354 for (unsigned int l = 0; l < 4; ++l)
355 lohvs(index_table, 8 + line_coordinates_z[l] * M + i) =
357 }
358
359 // quads
360 lohvs(index_table, 8 + 12 * M + 0 * M * M + i * M + j) =
361 1. - p[0];
362 lohvs(index_table, 8 + 12 * M + 1 * M * M + i * M + j) = p[0];
363 lohvs(index_table, 8 + 12 * M + 2 * M * M + k * M + i) =
364 1. - p[1];
365 lohvs(index_table, 8 + 12 * M + 3 * M * M + k * M + i) = p[1];
366 lohvs(index_table, 8 + 12 * M + 4 * M * M + j * M + k) =
367 1. - p[2];
368 lohvs(index_table, 8 + 12 * M + 5 * M * M + j * M + k) = p[2];
369 }
370
371 // the sum of weights of the points at the outer rim should be one.
372 // check this
373 for (unsigned int unit_point = 0; unit_point < n_inner; ++unit_point)
374 Assert(std::fabs(std::accumulate(lohvs[unit_point].begin(),
375 lohvs[unit_point].end(),
376 0.) -
377 1) < 1e-13 * polynomial_degree,
379
380 return lohvs;
381 }
382
383
384
389 inline std::vector<::Table<2, double>>
391 const unsigned int polynomial_degree,
392 const unsigned int dim)
393 {
394 Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
395 std::vector<::Table<2, double>> output(dim);
396 if (polynomial_degree <= 1)
397 return output;
398
399 // fill the 1d interior weights
400 QGaussLobatto<1> quadrature(polynomial_degree + 1);
401 output[0].reinit(polynomial_degree - 1,
403 for (unsigned int q = 0; q < polynomial_degree - 1; ++q)
404 for (const unsigned int i : GeometryInfo<1>::vertex_indices())
405 output[0](q, i) =
406 GeometryInfo<1>::d_linear_shape_function(quadrature.point(q + 1),
407 i);
408
409 if (dim > 1)
410 output[1] = compute_support_point_weights_on_quad(polynomial_degree);
411
412 if (dim > 2)
413 output[2] = compute_support_point_weights_on_hex(polynomial_degree);
414
415 return output;
416 }
417
418
419
423 template <int dim>
424 inline ::Table<2, double>
425 compute_support_point_weights_cell(const unsigned int polynomial_degree)
426 {
427 Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
428 if (polynomial_degree <= 1)
429 return ::Table<2, double>();
430
431 QGaussLobatto<dim> quadrature(polynomial_degree + 1);
432 const std::vector<unsigned int> h2l =
434
435 ::Table<2, double> output(quadrature.size() -
438 for (unsigned int q = 0; q < output.size(0); ++q)
439 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
441 quadrature.point(h2l[q + GeometryInfo<dim>::vertices_per_cell]), i);
442
443 return output;
444 }
445
446
447
455 template <int dim, int spacedim>
456 inline Point<spacedim>
458 const typename ::MappingQ<dim, spacedim>::InternalData &data)
459 {
460 AssertDimension(data.shape_values.size(),
461 data.mapping_support_points.size());
462
463 // use now the InternalData to compute the point in real space.
464 Point<spacedim> p_real;
465 for (unsigned int i = 0; i < data.mapping_support_points.size(); ++i)
466 p_real += data.mapping_support_points[i] * data.shape(0, i);
467
468 return p_real;
469 }
470
471
472
477 template <int dim, int spacedim, typename Number>
478 inline Point<dim, Number>
481 const Point<dim, Number> &initial_p_unit,
482 const ArrayView<const Point<spacedim>> &points,
483 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
484 const std::vector<unsigned int> &renumber,
485 const bool print_iterations_to_deallog = false)
486 {
487 if (print_iterations_to_deallog)
488 deallog << "Start MappingQ::do_transform_real_to_unit_cell for real "
489 << "point [ " << p << " ] " << std::endl;
490
491 AssertDimension(points.size(),
492 Utilities::pow(polynomials_1d.size(), dim));
493
494 // Newton iteration to solve
495 // f(x)=p(x)-p=0
496 // where we are looking for 'x' and p(x) is the forward transformation
497 // from unit to real cell. We solve this using a Newton iteration
498 // x_{n+1}=x_n-[f'(x)]^{-1}f(x)
499 // The start value is set to be the linear approximation to the cell
500
501 // The shape values and derivatives of the mapping at this point are
502 // previously computed.
503
504 Point<dim, Number> p_unit = initial_p_unit;
506 polynomials_1d, points, p_unit, polynomials_1d.size() == 2, renumber);
507
508 Tensor<1, spacedim, Number> f = p_real.first - p;
509
510 // early out if we already have our point in all SIMD lanes, i.e.,
511 // f.norm_square() < 1e-24 * p_real.second[0].norm_square(). To enable
512 // this step also for VectorizedArray where we do not have operator <,
513 // compare the result to zero which is defined for SIMD types
514 if (std::max(Number(0.),
515 f.norm_square() - 1e-24 * p_real.second[0].norm_square()) ==
516 Number(0.))
517 return p_unit;
518
519 // we need to compare the position of the computed p(x) against the
520 // given point 'p'. We will terminate the iteration and return 'x' if
521 // they are less than eps apart. The question is how to choose eps --
522 // or, put maybe more generally: in which norm we want these 'p' and
523 // 'p(x)' to be eps apart.
524 //
525 // the question is difficult since we may have to deal with very
526 // elongated cells where we may achieve 1e-12*h for the distance of
527 // these two points in the 'long' direction, but achieving this
528 // tolerance in the 'short' direction of the cell may not be possible
529 //
530 // what we do instead is then to terminate iterations if
531 // \| p(x) - p \|_A < eps
532 // where the A-norm is somehow induced by the transformation of the
533 // cell. in particular, we want to measure distances relative to the
534 // sizes of the cell in its principal directions.
535 //
536 // to define what exactly A should be, note that to first order we have
537 // the following (assuming that x* is the solution of the problem, i.e.,
538 // p(x*)=p):
539 // p(x) - p = p(x) - p(x*)
540 // = -grad p(x) * (x*-x) + higher order terms
541 // This suggest to measure with a norm that corresponds to
542 // A = {[grad p(x)]^T [grad p(x)]}^{-1}
543 // because then
544 // \| p(x) - p \|_A \approx \| x - x* \|
545 // Consequently, we will try to enforce that
546 // \| p(x) - p \|_A = \| f \| <= eps
547 //
548 // Note that using this norm is a bit dangerous since the norm changes
549 // in every iteration (A isn't fixed by depending on xk). However, if
550 // the cell is not too deformed (it may be stretched, but not twisted)
551 // then the mapping is almost linear and A is indeed constant or
552 // nearly so.
553 const double eps = 1.e-11;
554 const unsigned int newton_iteration_limit = 20;
555
556 Point<dim, Number> invalid_point;
557 invalid_point[0] = std::numeric_limits<double>::lowest();
558 bool tried_project_to_unit_cell = false;
559
560 unsigned int newton_iteration = 0;
561 Number f_weighted_norm_square = 1.;
562 Number last_f_weighted_norm_square = 1.;
563
564 do
565 {
566 if (print_iterations_to_deallog)
567 deallog << "Newton iteration " << newton_iteration
568 << " for unit point guess " << p_unit << std::endl;
569
570 // f'(x)
572 for (unsigned int d = 0; d < spacedim; ++d)
573 for (unsigned int e = 0; e < dim; ++e)
574 df[d][e] = p_real.second[e][d];
575
576 // Check determinant(df) > 0 on all SIMD lanes. The
577 // condition here is unreadable (but really corresponds to
578 // asking whether det(df) > 0 for all elements of the
579 // vector) because VectorizedArray does not have a member
580 // that can be used to check that all vector elements are
581 // positive. But VectorizedArray has a std::min() function,
582 // and operator==().
583 if (!(std::min(determinant(df),
584 Number(std::numeric_limits<double>::min())) ==
585 Number(std::numeric_limits<double>::min())))
586 {
587 // We allow to enter this function with an initial guess
588 // outside the unit cell. We might have run into invalid
589 // Jacobians due to that, so we should at least try once to go
590 // back to the unit cell and go on with the Newton iteration
591 // from there. Since the outside case is unlikely, we can
592 // afford spending the extra effort at this place.
593 if (tried_project_to_unit_cell == false)
594 {
597 polynomials_1d,
598 points,
599 p_unit,
600 polynomials_1d.size() == 2,
601 renumber);
602 f = p_real.first - p;
603 f_weighted_norm_square = 1.;
604 last_f_weighted_norm_square = 1;
605 tried_project_to_unit_cell = true;
606 continue;
607 }
608 else
609 return invalid_point;
610 }
611
612 // Solve [f'(x)]d=f(x)
613 const Tensor<2, spacedim, Number> df_inverse = invert(df);
614 const Tensor<1, spacedim, Number> delta = df_inverse * f;
615 last_f_weighted_norm_square = delta.norm_square();
616
617 if (print_iterations_to_deallog)
618 deallog << " delta=" << delta << std::endl;
619
620 // do a line search
621 double step_length = 1.0;
622 do
623 {
624 // update of p_unit. The spacedim-th component of transformed
625 // point is simply ignored in codimension one case. When this
626 // component is not zero, then we are projecting the point to
627 // the surface or curve identified by the cell.
628 Point<dim, Number> p_unit_trial = p_unit;
629 for (unsigned int i = 0; i < dim; ++i)
630 p_unit_trial[i] -= step_length * delta[i];
631
632 // shape values and derivatives at new p_unit point
633 const auto p_real_trial =
635 polynomials_1d,
636 points,
637 p_unit_trial,
638 polynomials_1d.size() == 2,
639 renumber);
640 const Tensor<1, spacedim, Number> f_trial =
641 p_real_trial.first - p;
642 f_weighted_norm_square = (df_inverse * f_trial).norm_square();
643
644 if (print_iterations_to_deallog)
645 {
646 deallog << " step_length=" << step_length << std::endl;
647 if (step_length == 1.0)
648 deallog << " ||f || =" << f.norm() << std::endl;
649 deallog << " ||f*|| =" << f_trial.norm() << std::endl
650 << " ||f*||_A ="
651 << std::sqrt(f_weighted_norm_square) << std::endl;
652 }
653
654 // See if we are making progress with the current step length
655 // and if not, reduce it by a factor of two and try again.
656 //
657 // Strictly speaking, we should probably use the same norm as we
658 // use for the outer algorithm. In practice, line search is just
659 // a crutch to find a "reasonable" step length, and so using the
660 // l2 norm is probably just fine.
661 //
662 // check f_trial.norm() < f.norm() in SIMD form. This is a bit
663 // more complicated because some SIMD lanes might not be doing
664 // any progress any more as they have already reached roundoff
665 // accuracy. We define that as the case of updates less than
666 // 1e-6. The tolerance might seem coarse but since we are
667 // dealing with a Newton iteration of a polynomial function we
668 // either converge quadratically or not any more. Thus, our
669 // selection is to terminate if either the norm of f is
670 // decreasing or that threshold of 1e-6 is reached.
671 if (std::max(f_weighted_norm_square - 1e-6 * 1e-6, Number(0.)) *
672 std::max(f_trial.norm_square() - f.norm_square(),
673 Number(0.)) ==
674 Number(0.))
675 {
676 p_real = p_real_trial;
677 p_unit = p_unit_trial;
678 f = f_trial;
679 break;
680 }
681 else if (step_length > 0.05)
682 step_length *= 0.5;
683 else
684 break;
685 }
686 while (true);
687
688 // In case we terminated the line search due to the step becoming
689 // too small, we give the iteration another try with the
690 // projection of the initial guess to the unit cell before we give
691 // up, just like for the negative determinant case.
692 if (step_length <= 0.05 && tried_project_to_unit_cell == false)
693 {
696 polynomials_1d,
697 points,
698 p_unit,
699 polynomials_1d.size() == 2,
700 renumber);
701 f = p_real.first - p;
702 f_weighted_norm_square = 1.;
703 last_f_weighted_norm_square = 1;
704 tried_project_to_unit_cell = true;
705 continue;
706 }
707 else if (step_length <= 0.05)
708 return invalid_point;
709
710 ++newton_iteration;
711 if (newton_iteration > newton_iteration_limit)
712 return invalid_point;
713 }
714 // Stop if f_weighted_norm_square <= eps^2 on all SIMD lanes or if the
715 // weighted norm is less than 1e-6 and the improvement against the
716 // previous step was less than a factor of 10 (in that regime, we
717 // either have quadratic convergence or roundoff errors due to a bad
718 // mapping)
719 while (
720 !(std::max(f_weighted_norm_square - eps * eps, Number(0.)) *
721 std::max(last_f_weighted_norm_square -
722 std::min(f_weighted_norm_square, Number(1e-6 * 1e-6)) *
723 100.,
724 Number(0.)) ==
725 Number(0.)));
726
727 if (print_iterations_to_deallog)
728 deallog << "Iteration converged for p_unit = [ " << p_unit
729 << " ] and iteration error "
730 << std::sqrt(f_weighted_norm_square) << std::endl;
731
732 return p_unit;
733 }
734
735
736
740 template <int dim>
741 inline Point<dim>
743 const Point<dim + 1> &p,
744 const Point<dim> &initial_p_unit,
745 const ArrayView<const Point<dim + 1>> &points,
746 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
747 const std::vector<unsigned int> &renumber)
748 {
749 const unsigned int spacedim = dim + 1;
750
751 AssertDimension(points.size(),
752 Utilities::pow(polynomials_1d.size(), dim));
753
754 Point<dim> p_unit = initial_p_unit;
755
756 const double eps = 1.e-12;
757 const unsigned int loop_limit = 10;
758
759 unsigned int loop = 0;
760 double f_weighted_norm_square = 1.;
761
762 while (f_weighted_norm_square > eps * eps && loop++ < loop_limit)
763 {
764 const auto p_real =
766 polynomials_1d,
767 points,
768 p_unit,
769 polynomials_1d.size() == 2,
770 renumber);
771 Tensor<1, spacedim> p_minus_F = p - p_real.first;
772 const DerivativeForm<1, spacedim, dim> DF = p_real.second;
773
775 polynomials_1d, points, p_unit, renumber);
776 Point<dim> f;
778 for (unsigned int j = 0; j < dim; ++j)
779 {
780 f[j] = DF[j] * p_minus_F;
781 for (unsigned int l = 0; l < dim; ++l)
782 df[j][l] = -DF[j] * DF[l] + hessian[j][l] * p_minus_F;
783 }
784
785 // Solve [df(x)]d=f(x)
786 const Tensor<1, dim> d =
787 invert(df) * static_cast<const Tensor<1, dim> &>(f);
788 f_weighted_norm_square = d.norm_square();
789 p_unit -= d;
790 }
791
792
793 // Here we check that in the last execution of while the first
794 // condition was already wrong, meaning the residual was below
795 // eps. Only if the first condition failed, loop will have been
796 // increased and tested, and thus have reached the limit.
797 AssertThrow(loop < loop_limit,
799
800 return p_unit;
801 }
802
803
804
822 template <int dim, int spacedim>
824 {
825 public:
829 static constexpr unsigned int n_functions =
830 (spacedim == 1 ? 3 : (spacedim == 2 ? 6 : 10));
831
844 const ArrayView<const Point<spacedim>> &real_support_points,
845 const std::vector<Point<dim>> &unit_support_points)
846 : normalization_shift(real_support_points[0])
848 1. / real_support_points[0].distance(real_support_points[1]))
849 , is_affine(true)
850 {
851 AssertDimension(real_support_points.size(), unit_support_points.size());
852
853 // For the bi-/trilinear approximation, we cannot build a quadratic
854 // polynomial due to a lack of points (interpolation matrix would
855 // get singular). Similarly, it is not entirely clear how to gather
856 // enough information for the case dim < spacedim.
857 //
858 // In both cases we require the vector real_support_points to
859 // contain the vertex positions and fall back to an affine
860 // approximation:
861 Assert(dim == spacedim || real_support_points.size() ==
864 if (real_support_points.size() == GeometryInfo<dim>::vertices_per_cell)
865 {
867 make_array_view(real_support_points));
869 affine.first.covariant_form().transpose();
870
871 // The code for evaluation assumes an additional transformation of
872 // the form (x - normalization_shift) * normalization_length --
873 // account for this in the definition of the coefficients.
874 coefficients[0] =
875 apply_transformation(A_inv, normalization_shift - affine.second);
876 for (unsigned int d = 0; d < spacedim; ++d)
877 for (unsigned int e = 0; e < dim; ++e)
878 coefficients[1 + d][e] =
879 A_inv[e][d] * (1.0 / normalization_length);
880 is_affine = true;
881 return;
882 }
883
885 std::array<double, n_functions> shape_values;
886 for (unsigned int q = 0; q < unit_support_points.size(); ++q)
887 {
888 // Evaluate quadratic shape functions in point, with the
889 // normalization applied in order to avoid roundoff issues with
890 // scaling far away from 1.
891 shape_values[0] = 1.;
892 const Tensor<1, spacedim> p_scaled =
894 (real_support_points[q] - normalization_shift);
895 for (unsigned int d = 0; d < spacedim; ++d)
896 shape_values[1 + d] = p_scaled[d];
897 for (unsigned int d = 0, c = 0; d < spacedim; ++d)
898 for (unsigned int e = 0; e <= d; ++e, ++c)
899 shape_values[1 + spacedim + c] = p_scaled[d] * p_scaled[e];
900
901 // Build lower diagonal of least squares matrix and rhs, the
902 // essential part being that we construct the matrix with the
903 // real points and the right hand side by comparing to the
904 // reference point positions which sets up an inverse
905 // interpolation.
906 for (unsigned int i = 0; i < n_functions; ++i)
907 for (unsigned int j = 0; j < n_functions; ++j)
908 matrix[i][j] += shape_values[i] * shape_values[j];
909 for (unsigned int i = 0; i < n_functions; ++i)
910 coefficients[i] += shape_values[i] * unit_support_points[q];
911 }
912
913 // Factorize the matrix A = L * L^T in-place with the
914 // Cholesky-Banachiewicz algorithm. The implementation is similar to
915 // FullMatrix::cholesky() but re-implemented to avoid memory
916 // allocations and some unnecessary divisions which we can do here as
917 // we only need to solve with dim right hand sides.
918 for (unsigned int i = 0; i < n_functions; ++i)
919 {
920 double Lij_sum = 0;
921 for (unsigned int j = 0; j < i; ++j)
922 {
923 double Lik_Ljk_sum = 0;
924 for (unsigned int k = 0; k < j; ++k)
925 Lik_Ljk_sum += matrix[i][k] * matrix[j][k];
926 matrix[i][j] = matrix[j][j] * (matrix[i][j] - Lik_Ljk_sum);
927 Lij_sum += matrix[i][j] * matrix[i][j];
928 }
929 AssertThrow(matrix[i][i] - Lij_sum >= 0,
930 ExcMessage("Matrix of normal equations not positive "
931 "definite"));
932
933 // Store the inverse in the diagonal since that is the quantity
934 // needed later in the factorization as well as the forward and
935 // backward substitution, minimizing the number of divisions.
936 matrix[i][i] = 1. / std::sqrt(matrix[i][i] - Lij_sum);
937 }
938
939 // Solve lower triangular part, L * y = rhs.
940 for (unsigned int i = 0; i < n_functions; ++i)
941 {
942 Point<dim> sum = coefficients[i];
943 for (unsigned int j = 0; j < i; ++j)
944 sum -= matrix[i][j] * coefficients[j];
945 coefficients[i] = sum * matrix[i][i];
946 }
947
948 // Solve upper triangular part, L^T * x = y (i.e., x = A^{-1} * rhs)
949 for (unsigned int i = n_functions; i > 0;)
950 {
951 --i;
952 Point<dim> sum = coefficients[i];
953 for (unsigned int j = i + 1; j < n_functions; ++j)
954 sum -= matrix[j][i] * coefficients[j];
955 coefficients[i] = sum * matrix[i][i];
956 }
957
958 // Check whether the approximation is indeed affine, allowing to
959 // skip the quadratic terms.
960 is_affine = true;
961 for (unsigned int i = dim + 1; i < n_functions; ++i)
962 if (coefficients[i].norm_square() > 1e-20)
963 {
964 is_affine = false;
965 break;
966 }
967 }
968
973 default;
974
978 template <typename Number>
981 {
982 Point<dim, Number> result;
983 for (unsigned int d = 0; d < dim; ++d)
984 result[d] = coefficients[0][d];
985
986 // Apply the normalization to ensure a good conditioning. Since Number
987 // might be a vectorized array whereas the normalization is a point of
988 // doubles, we cannot use the overload of operator- and must instead
989 // loop over the components of the point.
991 for (unsigned int d = 0; d < spacedim; ++d)
992 p_scaled[d] = (p[d] - normalization_shift[d]) * normalization_length;
993
994 for (unsigned int d = 0; d < spacedim; ++d)
995 result += coefficients[1 + d] * p_scaled[d];
996
997 if (!is_affine)
998 {
999 Point<dim, Number> result_affine = result;
1000 for (unsigned int d = 0, c = 0; d < spacedim; ++d)
1001 for (unsigned int e = 0; e <= d; ++e, ++c)
1002 result +=
1003 coefficients[1 + spacedim + c] * (p_scaled[d] * p_scaled[e]);
1004
1005 // Check if the quadratic approximation ends up considerably
1006 // farther outside the unit cell on some or all SIMD lanes than
1007 // the affine approximation - in that case, we switch those
1008 // components back to the affine approximation. Note that the
1009 // quadratic approximation will grow more quickly away from the
1010 // unit cell. We make the selection for each SIMD lane with a
1011 // ternary operation.
1012 const Number distance_to_unit_cell = result.distance_square(
1014 const Number affine_distance_to_unit_cell =
1015 result_affine.distance_square(
1017 for (unsigned int d = 0; d < dim; ++d)
1019 distance_to_unit_cell,
1020 affine_distance_to_unit_cell + 0.5,
1021 result_affine[d],
1022 result[d]);
1023 }
1024 return result;
1025 }
1026
1027 private:
1036
1041
1045 std::array<Point<dim>, n_functions> coefficients;
1046
1053 };
1054
1055
1056
1062 template <int dim, int spacedim>
1063 inline void
1065 const CellSimilarity::Similarity cell_similarity,
1066 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1067 std::vector<Point<spacedim>> &quadrature_points,
1068 std::vector<DerivativeForm<1, dim, spacedim>> &jacobians,
1069 std::vector<DerivativeForm<1, spacedim, dim>> &inverse_jacobians,
1070 std::vector<DerivativeForm<2, dim, spacedim>> &jacobian_grads)
1071 {
1072 const UpdateFlags update_flags = data.update_each;
1073
1074 using VectorizedArrayType =
1075 typename ::MappingQ<dim,
1076 spacedim>::InternalData::VectorizedArrayType;
1077 const unsigned int n_shape_values = data.n_shape_functions;
1078 const unsigned int n_q_points = data.shape_info.n_q_points;
1079 constexpr unsigned int n_lanes = VectorizedArrayType::size();
1080 constexpr unsigned int n_comp = 1 + (spacedim - 1) / n_lanes;
1081 constexpr unsigned int n_hessians = (dim * (dim + 1)) / 2;
1082
1083 // Since MappingQ::InternalData does not have separate arrays for the
1084 // covariant and contravariant transformations, but uses the arrays in
1085 // the `MappingRelatedData`, it can happen that vectors do not have the
1086 // right size
1087 if (update_flags & update_contravariant_transformation)
1088 jacobians.resize(n_q_points);
1089 if (update_flags & update_covariant_transformation)
1090 inverse_jacobians.resize(n_q_points);
1091
1092 EvaluationFlags::EvaluationFlags evaluation_flag =
1095 ((cell_similarity != CellSimilarity::translation) &&
1096 (update_flags & update_contravariant_transformation) ?
1099 ((cell_similarity != CellSimilarity::translation) &&
1100 (update_flags & update_jacobian_grads) ?
1103
1104 Assert(!(evaluation_flag & EvaluationFlags::values) || n_q_points > 0,
1106 Assert(!(evaluation_flag & EvaluationFlags::values) ||
1107 n_q_points == quadrature_points.size(),
1108 ExcDimensionMismatch(n_q_points, quadrature_points.size()));
1109 Assert(!(evaluation_flag & EvaluationFlags::gradients) ||
1110 data.n_shape_functions > 0,
1112 Assert(!(evaluation_flag & EvaluationFlags::hessians) ||
1113 n_q_points == jacobian_grads.size(),
1114 ExcDimensionMismatch(n_q_points, jacobian_grads.size()));
1115
1116 // shortcut in case we have an identity interpolation and only request
1117 // the quadrature points
1118 if (evaluation_flag == EvaluationFlags::values &&
1119 data.shape_info.element_type ==
1121 {
1122 for (unsigned int q = 0; q < n_q_points; ++q)
1123 quadrature_points[q] =
1124 data.mapping_support_points[data.shape_info
1125 .lexicographic_numbering[q]];
1126 return;
1127 }
1128
1130
1131 // prepare arrays
1132 if (evaluation_flag != EvaluationFlags::nothing)
1133 {
1134 eval.set_data_pointers(&data.scratch, n_comp);
1135
1136 // make sure to initialize on all lanes also when some are unused in
1137 // the code below
1138 for (unsigned int i = 0; i < n_shape_values * n_comp; ++i)
1139 eval.begin_dof_values()[i] = VectorizedArrayType();
1140
1141 const std::vector<unsigned int> &renumber_to_lexicographic =
1142 data.shape_info.lexicographic_numbering;
1143 for (unsigned int i = 0; i < n_shape_values; ++i)
1144 for (unsigned int d = 0; d < spacedim; ++d)
1145 {
1146 const unsigned int in_comp = d % n_lanes;
1147 const unsigned int out_comp = d / n_lanes;
1148 eval
1149 .begin_dof_values()[out_comp * n_shape_values + i][in_comp] =
1150 data.mapping_support_points[renumber_to_lexicographic[i]][d];
1151 }
1152
1153 // do the actual tensorized evaluation
1155 n_comp, evaluation_flag, eval.begin_dof_values(), eval);
1156 }
1157
1158 // do the postprocessing
1159 if (evaluation_flag & EvaluationFlags::values)
1160 {
1161 for (unsigned int out_comp = 0; out_comp < n_comp; ++out_comp)
1162 for (unsigned int i = 0; i < n_q_points; ++i)
1163 for (unsigned int in_comp = 0;
1164 in_comp < n_lanes && in_comp < spacedim - out_comp * n_lanes;
1165 ++in_comp)
1166 quadrature_points[i][out_comp * n_lanes + in_comp] =
1167 eval.begin_values()[out_comp * n_q_points + i][in_comp];
1168 }
1169
1170 if (evaluation_flag & EvaluationFlags::gradients)
1171 {
1172 // We need to reinterpret the data after evaluate has been applied.
1173 for (unsigned int out_comp = 0; out_comp < n_comp; ++out_comp)
1174 for (unsigned int point = 0; point < n_q_points; ++point)
1175 for (unsigned int in_comp = 0;
1176 in_comp < n_lanes && in_comp < spacedim - out_comp * n_lanes;
1177 ++in_comp)
1178 for (unsigned int j = 0; j < dim; ++j)
1179 {
1180 jacobians[point][out_comp * n_lanes + in_comp][j] =
1181 eval.begin_gradients()[(out_comp * n_q_points + point) *
1182 dim +
1183 j][in_comp];
1184 }
1185 }
1186
1187 if (update_flags & update_volume_elements)
1188 if (cell_similarity != CellSimilarity::translation)
1189 for (unsigned int point = 0; point < n_q_points; ++point)
1190 data.volume_elements[point] = jacobians[point].determinant();
1191
1192 // copy values from InternalData to vector given by reference
1193 if (update_flags & update_covariant_transformation)
1194 {
1195 AssertDimension(jacobians.size(), n_q_points);
1196 AssertDimension(inverse_jacobians.size(), n_q_points);
1197 if (cell_similarity != CellSimilarity::translation)
1198 for (unsigned int point = 0; point < n_q_points; ++point)
1199 inverse_jacobians[point] =
1200 jacobians[point].covariant_form().transpose();
1201 }
1202
1203 if (evaluation_flag & EvaluationFlags::hessians)
1204 {
1205 constexpr int desymmetrize_3d[6][2] = {
1206 {0, 0}, {1, 1}, {2, 2}, {0, 1}, {0, 2}, {1, 2}};
1207 constexpr int desymmetrize_2d[3][2] = {{0, 0}, {1, 1}, {0, 1}};
1208
1209 // We need to reinterpret the data after evaluate has been applied.
1210 for (unsigned int out_comp = 0; out_comp < n_comp; ++out_comp)
1211 for (unsigned int point = 0; point < n_q_points; ++point)
1212 for (unsigned int j = 0; j < n_hessians; ++j)
1213 for (unsigned int in_comp = 0;
1214 in_comp < n_lanes &&
1215 in_comp < spacedim - out_comp * n_lanes;
1216 ++in_comp)
1217 {
1218 const unsigned int total_number = point * n_hessians + j;
1219 const unsigned int new_point = total_number % n_q_points;
1220 const unsigned int new_hessian_comp =
1221 total_number / n_q_points;
1222 const unsigned int new_hessian_comp_i =
1223 dim == 2 ? desymmetrize_2d[new_hessian_comp][0] :
1224 desymmetrize_3d[new_hessian_comp][0];
1225 const unsigned int new_hessian_comp_j =
1226 dim == 2 ? desymmetrize_2d[new_hessian_comp][1] :
1227 desymmetrize_3d[new_hessian_comp][1];
1228 const double value =
1229 eval.begin_hessians()[(out_comp * n_q_points + point) *
1230 n_hessians +
1231 j][in_comp];
1232 jacobian_grads[new_point][out_comp * n_lanes + in_comp]
1233 [new_hessian_comp_i][new_hessian_comp_j] =
1234 value;
1235 jacobian_grads[new_point][out_comp * n_lanes + in_comp]
1236 [new_hessian_comp_j][new_hessian_comp_i] =
1237 value;
1238 }
1239 }
1240 }
1241
1242
1243
1244 template <int dim, int spacedim>
1245 inline void
1247 const CellSimilarity::Similarity cell_similarity,
1248 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1249 const ArrayView<const Point<dim>> &unit_points,
1250 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1251 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1252 std::vector<Point<spacedim>> &quadrature_points,
1253 std::vector<DerivativeForm<1, dim, spacedim>> &jacobians,
1254 std::vector<DerivativeForm<1, spacedim, dim>> &inverse_jacobians)
1255 {
1256 const UpdateFlags update_flags = data.update_each;
1257 const ArrayView<const Point<spacedim>> support_points(
1258 data.mapping_support_points);
1259
1260 const unsigned int n_points = unit_points.size();
1261 const unsigned int n_lanes = VectorizedArray<double>::size();
1262
1263 // Since MappingQ::InternalData does not have separate arrays for the
1264 // covariant and contravariant transformations, but uses the arrays in
1265 // the `MappingRelatedData`, it can happen that vectors do not have the
1266 // right size
1267 if (update_flags & update_contravariant_transformation)
1268 jacobians.resize(n_points);
1269 if (update_flags & update_covariant_transformation)
1270 inverse_jacobians.resize(n_points);
1271
1272 const bool is_translation =
1273 cell_similarity == CellSimilarity::translation;
1274
1275 const bool needs_gradient =
1276 !is_translation &&
1277 (update_flags &
1280
1281 // Use the more heavy VectorizedArray code path if there is more than
1282 // one point left to compute
1283 for (unsigned int i = 0; i < n_points; i += n_lanes)
1284 if (n_points - i > 1)
1285 {
1287 for (unsigned int j = 0; j < n_lanes; ++j)
1288 if (i + j < n_points)
1289 for (unsigned int d = 0; d < dim; ++d)
1290 p_vec[d][j] = unit_points[i + j][d];
1291 else
1292 for (unsigned int d = 0; d < dim; ++d)
1293 p_vec[d][j] = unit_points[i][d];
1294
1297 derivative;
1298 if (needs_gradient)
1299 {
1300 const auto result =
1302 polynomials_1d,
1303 support_points,
1304 p_vec,
1305 polynomials_1d.size() == 2,
1306 renumber_lexicographic_to_hierarchic);
1307
1308 value = result.first;
1309
1310 for (unsigned int d = 0; d < spacedim; ++d)
1311 for (unsigned int e = 0; e < dim; ++e)
1312 derivative[d][e] = result.second[e][d];
1313 }
1314 else
1316 polynomials_1d,
1317 support_points,
1318 p_vec,
1319 polynomials_1d.size() == 2,
1320 renumber_lexicographic_to_hierarchic);
1321
1322 if (update_flags & update_quadrature_points)
1323 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1324 for (unsigned int d = 0; d < spacedim; ++d)
1325 quadrature_points[i + j][d] = value[d][j];
1326
1327 if (is_translation)
1328 continue;
1329
1330 if (update_flags & update_contravariant_transformation)
1331 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1332 for (unsigned int d = 0; d < spacedim; ++d)
1333 for (unsigned int e = 0; e < dim; ++e)
1334 jacobians[i + j][d][e] = derivative[d][e][j];
1335
1336 if (update_flags & update_volume_elements)
1337 {
1338 const auto determinant = derivative.determinant();
1339 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1340 data.volume_elements[i + j] = determinant[j];
1341 }
1342
1343 if (update_flags & update_covariant_transformation)
1344 {
1345 const auto covariant = derivative.covariant_form();
1346 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1347 for (unsigned int d = 0; d < dim; ++d)
1348 for (unsigned int e = 0; e < spacedim; ++e)
1349 inverse_jacobians[i + j][d][e] = covariant[e][d][j];
1350 }
1351 }
1352 else
1353 {
1356 if (needs_gradient)
1357 {
1358 const auto result =
1360 polynomials_1d,
1361 support_points,
1362 unit_points[i],
1363 polynomials_1d.size() == 2,
1364 renumber_lexicographic_to_hierarchic);
1365
1366 value = result.first;
1367
1368 for (unsigned int d = 0; d < spacedim; ++d)
1369 for (unsigned int e = 0; e < dim; ++e)
1370 derivative[d][e] = result.second[e][d];
1371 }
1372 else
1374 polynomials_1d,
1375 support_points,
1376 unit_points[i],
1377 polynomials_1d.size() == 2,
1378 renumber_lexicographic_to_hierarchic);
1379
1380 if (update_flags & update_quadrature_points)
1381 quadrature_points[i] = value;
1382
1383 if (is_translation)
1384 continue;
1385
1386 if (dim == spacedim && update_flags & update_volume_elements)
1387 data.volume_elements[i] = derivative.determinant();
1388
1389 if (update_flags & update_contravariant_transformation)
1390 jacobians[i] = derivative;
1391
1392 if (update_flags & update_covariant_transformation)
1393 inverse_jacobians[i] = jacobians[i].covariant_form().transpose();
1394 }
1395 }
1396
1397
1398
1405 template <int dim, int spacedim>
1406 inline void
1408 const CellSimilarity::Similarity cell_similarity,
1409 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1410 const ArrayView<const Point<dim>> &unit_points,
1411 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1412 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1413 std::vector<DerivativeForm<2, dim, spacedim>> &jacobian_grads)
1414 {
1415 if (data.update_each & update_jacobian_grads)
1416 {
1417 const ArrayView<const Point<spacedim>> support_points(
1418 data.mapping_support_points);
1419 const unsigned int n_q_points = jacobian_grads.size();
1420
1421 if (cell_similarity != CellSimilarity::translation)
1422 for (unsigned int point = 0; point < n_q_points; ++point)
1423 {
1426 polynomials_1d,
1427 support_points,
1428 unit_points[point],
1429 renumber_lexicographic_to_hierarchic);
1430
1431 for (unsigned int i = 0; i < spacedim; ++i)
1432 for (unsigned int j = 0; j < dim; ++j)
1433 for (unsigned int l = 0; l < dim; ++l)
1434 jacobian_grads[point][i][j][l] = second[j][l][i];
1435 }
1436 }
1437 }
1438
1439
1440
1447 template <int dim, int spacedim>
1448 inline void
1450 const CellSimilarity::Similarity cell_similarity,
1451 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1452 const ArrayView<const Point<dim>> &unit_points,
1453 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1454 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1455 std::vector<Tensor<3, spacedim>> &jacobian_pushed_forward_grads)
1456 {
1457 if (data.update_each & update_jacobian_pushed_forward_grads)
1458 {
1459 const ArrayView<const Point<spacedim>> support_points(
1460 data.mapping_support_points);
1461 const unsigned int n_q_points = jacobian_pushed_forward_grads.size();
1462
1463 if (cell_similarity != CellSimilarity::translation)
1464 {
1465 double tmp[spacedim][spacedim][spacedim];
1466 for (unsigned int point = 0; point < n_q_points; ++point)
1467 {
1470 polynomials_1d,
1471 support_points,
1472 unit_points[point],
1473 renumber_lexicographic_to_hierarchic);
1474 const DerivativeForm<1, dim, spacedim> covariant =
1475 data.output_data->inverse_jacobians[point].transpose();
1476
1477 // first push forward the j-components
1478 for (unsigned int i = 0; i < spacedim; ++i)
1479 for (unsigned int j = 0; j < spacedim; ++j)
1480 for (unsigned int l = 0; l < dim; ++l)
1481 {
1482 tmp[i][j][l] = second[0][l][i] * covariant[j][0];
1483 for (unsigned int jr = 1; jr < dim; ++jr)
1484 {
1485 tmp[i][j][l] +=
1486 second[jr][l][i] * covariant[j][jr];
1487 }
1488 }
1489
1490 // now, pushing forward the l-components
1491 for (unsigned int i = 0; i < spacedim; ++i)
1492 for (unsigned int j = 0; j < spacedim; ++j)
1493 for (unsigned int l = 0; l < spacedim; ++l)
1494 {
1495 jacobian_pushed_forward_grads[point][i][j][l] =
1496 tmp[i][j][0] * covariant[l][0];
1497 for (unsigned int lr = 1; lr < dim; ++lr)
1498 {
1499 jacobian_pushed_forward_grads[point][i][j][l] +=
1500 tmp[i][j][lr] * covariant[l][lr];
1501 }
1502 }
1503 }
1504 }
1505 }
1506 }
1507
1508
1509
1510 template <int dim, int spacedim, int length_tensor>
1513 const Tensor<1, length_tensor, Tensor<1, spacedim>> &compressed)
1514 {
1515 Assert(dim >= 1 && dim <= 3, ExcNotImplemented());
1517 for (unsigned int i = 0; i < spacedim; ++i)
1518 {
1519 if (dim == 1)
1520 result[i][0][0][0] = compressed[0][i];
1521 else if (dim >= 2)
1522 {
1523 for (unsigned int d = 0; d < 2; ++d)
1524 for (unsigned int e = 0; e < 2; ++e)
1525 for (unsigned int f = 0; f < 2; ++f)
1526 result[i][d][e][f] = compressed[d + e + f][i];
1527 if (dim == 3)
1528 {
1529 AssertDimension(length_tensor, 10);
1530
1531 // the derivatives are ordered in rows by increasing z
1532 // derivative, and in each row we have x^{(n-j)} y^{(j)} as
1533 // j walks through the columns
1534 for (unsigned int d = 0; d < 2; ++d)
1535 for (unsigned int e = 0; e < 2; ++e)
1536 {
1537 result[i][d][e][2] = compressed[4 + d + e][i];
1538 result[i][d][2][e] = compressed[4 + d + e][i];
1539 result[i][2][d][e] = compressed[4 + d + e][i];
1540 }
1541 for (unsigned int d = 0; d < 2; ++d)
1542 {
1543 result[i][d][2][2] = compressed[7 + d][i];
1544 result[i][2][d][2] = compressed[7 + d][i];
1545 result[i][2][2][d] = compressed[7 + d][i];
1546 }
1547 result[i][2][2][2] = compressed[9][i];
1548 }
1549 }
1550 }
1551 return result;
1552 }
1553
1554
1555
1562 template <int dim, int spacedim>
1563 inline void
1565 const CellSimilarity::Similarity cell_similarity,
1566 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1567 const ArrayView<const Point<dim>> &unit_points,
1568 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1569 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1570 std::vector<DerivativeForm<3, dim, spacedim>> &jacobian_2nd_derivatives)
1571 {
1572 if (data.update_each & update_jacobian_2nd_derivatives)
1573 {
1574 const ArrayView<const Point<spacedim>> support_points(
1575 data.mapping_support_points);
1576 const unsigned int n_q_points = jacobian_2nd_derivatives.size();
1577
1578 if (cell_similarity != CellSimilarity::translation)
1579 {
1580 for (unsigned int point = 0; point < n_q_points; ++point)
1581 {
1582 jacobian_2nd_derivatives[point] = expand_3rd_derivatives<dim>(
1584 polynomials_1d,
1585 support_points,
1586 unit_points[point],
1587 renumber_lexicographic_to_hierarchic));
1588 }
1589 }
1590 }
1591 }
1592
1593
1594
1602 template <int dim, int spacedim>
1603 inline void
1605 const CellSimilarity::Similarity cell_similarity,
1606 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1607 const ArrayView<const Point<dim>> &unit_points,
1608 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1609 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1610 std::vector<Tensor<4, spacedim>> &jacobian_pushed_forward_2nd_derivatives)
1611 {
1613 {
1614 const ArrayView<const Point<spacedim>> support_points(
1615 data.mapping_support_points);
1616 const unsigned int n_q_points =
1617 jacobian_pushed_forward_2nd_derivatives.size();
1618
1619 if (cell_similarity != CellSimilarity::translation)
1620 {
1623 for (unsigned int point = 0; point < n_q_points; ++point)
1624 {
1628 polynomials_1d,
1629 support_points,
1630 unit_points[point],
1631 renumber_lexicographic_to_hierarchic));
1632 const DerivativeForm<1, dim, spacedim> covariant =
1633 data.output_data->inverse_jacobians[point].transpose();
1634
1635 // push forward the j-coordinate
1636 for (unsigned int i = 0; i < spacedim; ++i)
1637 for (unsigned int j = 0; j < spacedim; ++j)
1638 for (unsigned int l = 0; l < dim; ++l)
1639 for (unsigned int m = 0; m < dim; ++m)
1640 {
1641 tmp[i][j][l][m] =
1642 third[i][0][l][m] * covariant[j][0];
1643 for (unsigned int jr = 1; jr < dim; ++jr)
1644 tmp[i][j][l][m] +=
1645 third[i][jr][l][m] * covariant[j][jr];
1646 }
1647
1648 // push forward the l-coordinate
1649 for (unsigned int i = 0; i < spacedim; ++i)
1650 for (unsigned int j = 0; j < spacedim; ++j)
1651 for (unsigned int l = 0; l < spacedim; ++l)
1652 for (unsigned int m = 0; m < dim; ++m)
1653 {
1654 tmp2[i][j][l][m] =
1655 tmp[i][j][0][m] * covariant[l][0];
1656 for (unsigned int lr = 1; lr < dim; ++lr)
1657 tmp2[i][j][l][m] +=
1658 tmp[i][j][lr][m] * covariant[l][lr];
1659 }
1660
1661 // push forward the m-coordinate
1662 for (unsigned int i = 0; i < spacedim; ++i)
1663 for (unsigned int j = 0; j < spacedim; ++j)
1664 for (unsigned int l = 0; l < spacedim; ++l)
1665 for (unsigned int m = 0; m < spacedim; ++m)
1666 {
1667 jacobian_pushed_forward_2nd_derivatives
1668 [point][i][j][l][m] =
1669 tmp2[i][j][l][0] * covariant[m][0];
1670 for (unsigned int mr = 1; mr < dim; ++mr)
1671 jacobian_pushed_forward_2nd_derivatives[point][i]
1672 [j][l]
1673 [m] +=
1674 tmp2[i][j][l][mr] * covariant[m][mr];
1675 }
1676 }
1677 }
1678 }
1679 }
1680
1681
1682
1683 template <int dim, int spacedim, int length_tensor>
1686 const Tensor<1, length_tensor, Tensor<1, spacedim>> &compressed)
1687 {
1688 Assert(dim >= 1 && dim <= 3, ExcNotImplemented());
1690 for (unsigned int i = 0; i < spacedim; ++i)
1691 {
1692 if (dim == 1)
1693 result[i][0][0][0][0] = compressed[0][i];
1694 else if (dim >= 2)
1695 {
1696 for (unsigned int d = 0; d < 2; ++d)
1697 for (unsigned int e = 0; e < 2; ++e)
1698 for (unsigned int f = 0; f < 2; ++f)
1699 for (unsigned int g = 0; g < 2; ++g)
1700 result[i][d][e][f][g] = compressed[d + e + f + g][i];
1701 if (dim == 3)
1702 {
1703 AssertDimension(length_tensor, 15);
1704
1705 // the derivatives are ordered in rows by increasing z
1706 // derivative, and in each row we have x^{(n-j)} y^{(j)} as
1707 // j walks through the columns
1708 for (unsigned int d = 0; d < 2; ++d)
1709 for (unsigned int e = 0; e < 2; ++e)
1710 for (unsigned int f = 0; f < 2; ++f)
1711 {
1712 result[i][d][e][f][2] = compressed[5 + d + e + f][i];
1713 result[i][d][e][2][f] = compressed[5 + d + e + f][i];
1714 result[i][d][2][e][f] = compressed[5 + d + e + f][i];
1715 result[i][2][d][e][f] = compressed[5 + d + e + f][i];
1716 }
1717 for (unsigned int d = 0; d < 2; ++d)
1718 for (unsigned int e = 0; e < 2; ++e)
1719 {
1720 result[i][d][e][2][2] = compressed[9 + d + e][i];
1721 result[i][d][2][e][2] = compressed[9 + d + e][i];
1722 result[i][d][2][2][e] = compressed[9 + d + e][i];
1723 result[i][2][d][e][2] = compressed[9 + d + e][i];
1724 result[i][2][d][2][e] = compressed[9 + d + e][i];
1725 result[i][2][2][d][e] = compressed[9 + d + e][i];
1726 }
1727 for (unsigned int d = 0; d < 2; ++d)
1728 {
1729 result[i][d][2][2][2] = compressed[12 + d][i];
1730 result[i][2][d][2][2] = compressed[12 + d][i];
1731 result[i][2][2][d][2] = compressed[12 + d][i];
1732 result[i][2][2][2][d] = compressed[12 + d][i];
1733 }
1734 result[i][2][2][2][2] = compressed[14][i];
1735 }
1736 }
1737 }
1738 return result;
1739 }
1740
1741
1742
1749 template <int dim, int spacedim>
1750 inline void
1752 const CellSimilarity::Similarity cell_similarity,
1753 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1754 const ArrayView<const Point<dim>> &unit_points,
1755 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1756 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1757 std::vector<DerivativeForm<4, dim, spacedim>> &jacobian_3rd_derivatives)
1758 {
1759 if (data.update_each & update_jacobian_3rd_derivatives)
1760 {
1761 const ArrayView<const Point<spacedim>> support_points(
1762 data.mapping_support_points);
1763 const unsigned int n_q_points = jacobian_3rd_derivatives.size();
1764
1765 if (cell_similarity != CellSimilarity::translation)
1766 {
1767 for (unsigned int point = 0; point < n_q_points; ++point)
1768 {
1769 jacobian_3rd_derivatives[point] = expand_4th_derivatives<dim>(
1771 polynomials_1d,
1772 support_points,
1773 unit_points[point],
1774 renumber_lexicographic_to_hierarchic));
1775 }
1776 }
1777 }
1778 }
1779
1780
1781
1789 template <int dim, int spacedim>
1790 inline void
1792 const CellSimilarity::Similarity cell_similarity,
1793 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1794 const ArrayView<const Point<dim>> &unit_points,
1795 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
1796 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
1797 std::vector<Tensor<5, spacedim>> &jacobian_pushed_forward_3rd_derivatives)
1798 {
1800 {
1801 const ArrayView<const Point<spacedim>> support_points(
1802 data.mapping_support_points);
1803 const unsigned int n_q_points =
1804 jacobian_pushed_forward_3rd_derivatives.size();
1805
1806 if (cell_similarity != CellSimilarity::translation)
1807 {
1808 ::
1809 ndarray<double, spacedim, spacedim, spacedim, spacedim, dim>
1810 tmp;
1812 tmp2;
1813 for (unsigned int point = 0; point < n_q_points; ++point)
1814 {
1818 polynomials_1d,
1819 support_points,
1820 unit_points[point],
1821 renumber_lexicographic_to_hierarchic));
1822
1823 const DerivativeForm<1, dim, spacedim> covariant =
1824 data.output_data->inverse_jacobians[point].transpose();
1825
1826 // push-forward the j-coordinate
1827 for (unsigned int i = 0; i < spacedim; ++i)
1828 for (unsigned int j = 0; j < spacedim; ++j)
1829 for (unsigned int l = 0; l < dim; ++l)
1830 for (unsigned int m = 0; m < dim; ++m)
1831 for (unsigned int n = 0; n < dim; ++n)
1832 {
1833 tmp[i][j][l][m][n] =
1834 fourth[i][0][l][m][n] * covariant[j][0];
1835 for (unsigned int jr = 1; jr < dim; ++jr)
1836 tmp[i][j][l][m][n] +=
1837 fourth[i][jr][l][m][n] * covariant[j][jr];
1838 }
1839
1840 // push-forward the l-coordinate
1841 for (unsigned int i = 0; i < spacedim; ++i)
1842 for (unsigned int j = 0; j < spacedim; ++j)
1843 for (unsigned int l = 0; l < spacedim; ++l)
1844 for (unsigned int m = 0; m < dim; ++m)
1845 for (unsigned int n = 0; n < dim; ++n)
1846 {
1847 tmp2[i][j][l][m][n] =
1848 tmp[i][j][0][m][n] * covariant[l][0];
1849 for (unsigned int lr = 1; lr < dim; ++lr)
1850 tmp2[i][j][l][m][n] +=
1851 tmp[i][j][lr][m][n] * covariant[l][lr];
1852 }
1853
1854 // push-forward the m-coordinate
1855 for (unsigned int i = 0; i < spacedim; ++i)
1856 for (unsigned int j = 0; j < spacedim; ++j)
1857 for (unsigned int l = 0; l < spacedim; ++l)
1858 for (unsigned int m = 0; m < spacedim; ++m)
1859 for (unsigned int n = 0; n < dim; ++n)
1860 {
1861 tmp[i][j][l][m][n] =
1862 tmp2[i][j][l][0][n] * covariant[m][0];
1863 for (unsigned int mr = 1; mr < dim; ++mr)
1864 tmp[i][j][l][m][n] +=
1865 tmp2[i][j][l][mr][n] * covariant[m][mr];
1866 }
1867
1868 // push-forward the n-coordinate
1869 for (unsigned int i = 0; i < spacedim; ++i)
1870 for (unsigned int j = 0; j < spacedim; ++j)
1871 for (unsigned int l = 0; l < spacedim; ++l)
1872 for (unsigned int m = 0; m < spacedim; ++m)
1873 for (unsigned int n = 0; n < spacedim; ++n)
1874 {
1875 jacobian_pushed_forward_3rd_derivatives
1876 [point][i][j][l][m][n] =
1877 tmp[i][j][l][m][0] * covariant[n][0];
1878 for (unsigned int nr = 1; nr < dim; ++nr)
1879 jacobian_pushed_forward_3rd_derivatives[point]
1880 [i][j][l]
1881 [m][n] +=
1882 tmp[i][j][l][m][nr] * covariant[n][nr];
1883 }
1884 }
1885 }
1886 }
1887 }
1888
1889
1890
1900 template <int dim, int spacedim>
1901 inline void
1903 const ::MappingQ<dim, spacedim> &mapping,
1904 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
1905 const unsigned int face_no,
1906 const unsigned int subface_no,
1907 const unsigned int n_q_points,
1908 const std::vector<double> &weights,
1909 const typename ::MappingQ<dim, spacedim>::InternalData &data,
1911 &output_data)
1912 {
1913 const UpdateFlags update_flags = data.update_each;
1914
1915 if (update_flags &
1917 {
1918 if (update_flags & update_boundary_forms)
1919 AssertDimension(output_data.boundary_forms.size(), n_q_points);
1920 if (update_flags & update_normal_vectors)
1921 AssertDimension(output_data.normal_vectors.size(), n_q_points);
1922 if (update_flags & update_JxW_values)
1923 AssertDimension(output_data.JxW_values.size(), n_q_points);
1924
1925 Assert(data.aux.size() + 1 >= dim, ExcInternalError());
1926
1927 // first compute some common data that is used for evaluating
1928 // all of the flags below
1929
1930 // map the unit tangentials to the real cell. checking for d!=dim-1
1931 // eliminates compiler warnings regarding unsigned int expressions <
1932 // 0.
1933 for (unsigned int d = 0; d != dim - 1; ++d)
1934 {
1935 const unsigned int vector_index =
1937 Assert(vector_index < data.unit_tangentials.size(),
1939 Assert(data.aux[d].size() <=
1940 data.unit_tangentials[vector_index].size(),
1942 mapping.transform(make_array_view(
1943 data.unit_tangentials[vector_index]),
1945 data,
1946 make_array_view(data.aux[d]));
1947 }
1948
1949 if (update_flags & update_boundary_forms)
1950 {
1951 // if dim==spacedim, we can use the unit tangentials to compute
1952 // the boundary form by simply taking the cross product
1953 if (dim == spacedim)
1954 {
1955 for (unsigned int i = 0; i < n_q_points; ++i)
1956 switch (dim)
1957 {
1958 case 1:
1959 // in 1d, we don't have access to any of the
1960 // data.aux fields (because it has only dim-1
1961 // components), but we can still compute the
1962 // boundary form by simply looking at the number of
1963 // the face
1964 output_data.boundary_forms[i][0] =
1965 (face_no == 0 ? -1 : +1);
1966 break;
1967 case 2:
1968 output_data.boundary_forms[i] =
1969 cross_product_2d(data.aux[0][i]);
1970 break;
1971 case 3:
1972 output_data.boundary_forms[i] =
1973 cross_product_3d(data.aux[0][i], data.aux[1][i]);
1974 break;
1975 default:
1977 }
1978 }
1979 else //(dim < spacedim)
1980 {
1981 // in the codim-one case, the boundary form results from the
1982 // cross product of all the face tangential vectors and the
1983 // cell normal vector
1984 //
1985 // to compute the cell normal, use the same method used in
1986 // fill_fe_values for cells above
1987 AssertDimension(data.output_data->jacobians.size(),
1988 n_q_points);
1989
1990 for (unsigned int point = 0; point < n_q_points; ++point)
1991 {
1992 const DerivativeForm<1, dim, spacedim> contravariant =
1993 data.output_data->jacobians[point];
1994 if (dim == 1)
1995 {
1996 // J is a tangent vector
1997 output_data.boundary_forms[point] =
1998 contravariant.transpose()[0];
1999 output_data.boundary_forms[point] /=
2000 (face_no == 0 ? -1. : +1.) *
2001 output_data.boundary_forms[point].norm();
2002 }
2003
2004 if (dim == 2)
2005 {
2007 contravariant.transpose();
2008
2009 Tensor<1, spacedim> cell_normal =
2010 cross_product_3d(DX_t[0], DX_t[1]);
2011 cell_normal /= cell_normal.norm();
2012
2013 // then compute the face normal from the face
2014 // tangent and the cell normal:
2015 output_data.boundary_forms[point] =
2016 cross_product_3d(data.aux[0][point], cell_normal);
2017 }
2018 }
2019 }
2020 }
2021
2022 if (update_flags & update_JxW_values)
2023 for (unsigned int i = 0; i < output_data.boundary_forms.size(); ++i)
2024 {
2025 output_data.JxW_values[i] =
2026 output_data.boundary_forms[i].norm() * weights[i];
2027
2028 if (subface_no != numbers::invalid_unsigned_int)
2029 {
2030 const double area_ratio = GeometryInfo<dim>::subface_ratio(
2031 cell->subface_case(face_no), subface_no);
2032 output_data.JxW_values[i] *= area_ratio;
2033 }
2034 }
2035
2036 if (update_flags & update_normal_vectors)
2037 for (unsigned int i = 0; i < output_data.normal_vectors.size(); ++i)
2038 output_data.normal_vectors[i] =
2039 Point<spacedim>(output_data.boundary_forms[i] /
2040 output_data.boundary_forms[i].norm());
2041 }
2042 }
2043
2044
2045
2052 template <int dim, int spacedim>
2053 inline void
2055 const ::MappingQ<dim, spacedim> &mapping,
2056 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
2057 const unsigned int face_no,
2058 const unsigned int subface_no,
2059 const typename QProjector<dim>::DataSetDescriptor data_set,
2060 const Quadrature<dim - 1> &quadrature,
2061 const typename ::MappingQ<dim, spacedim>::InternalData &data,
2062 const std::vector<Polynomials::Polynomial<double>> &polynomials_1d,
2063 const std::vector<unsigned int> &renumber_lexicographic_to_hierarchic,
2065 &output_data)
2066 {
2067 const ArrayView<const Point<dim>> quadrature_points(
2068 &data.quadrature_points[data_set], quadrature.size());
2069 if (dim > 1 && data.tensor_product_quadrature)
2070 {
2073 data,
2074 output_data.quadrature_points,
2075 output_data.jacobians,
2076 output_data.inverse_jacobians,
2077 output_data.jacobian_grads);
2078 }
2079 else
2080 {
2084 data,
2085 quadrature_points,
2086 polynomials_1d,
2087 renumber_lexicographic_to_hierarchic,
2088 output_data.quadrature_points,
2089 output_data.jacobians,
2090 output_data.inverse_jacobians);
2093 data,
2094 quadrature_points,
2095 polynomials_1d,
2096 renumber_lexicographic_to_hierarchic,
2097 output_data.jacobian_grads);
2098 }
2101 data,
2102 quadrature_points,
2103 polynomials_1d,
2104 renumber_lexicographic_to_hierarchic,
2105 output_data.jacobian_pushed_forward_grads);
2108 data,
2109 quadrature_points,
2110 polynomials_1d,
2111 renumber_lexicographic_to_hierarchic,
2112 output_data.jacobian_2nd_derivatives);
2115 data,
2116 quadrature_points,
2117 polynomials_1d,
2118 renumber_lexicographic_to_hierarchic,
2122 data,
2123 quadrature_points,
2124 polynomials_1d,
2125 renumber_lexicographic_to_hierarchic,
2126 output_data.jacobian_3rd_derivatives);
2129 data,
2130 quadrature_points,
2131 polynomials_1d,
2132 renumber_lexicographic_to_hierarchic,
2134
2136 cell,
2137 face_no,
2138 subface_no,
2139 quadrature.size(),
2140 quadrature.get_weights(),
2141 data,
2142 output_data);
2143 }
2144
2145
2146
2150 template <int dim, int spacedim, int rank>
2151 inline void
2153 const ArrayView<const Tensor<rank, dim>> &input,
2154 const MappingKind mapping_kind,
2155 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
2156 const ArrayView<Tensor<rank, spacedim>> &output)
2157 {
2158 AssertDimension(input.size(), output.size());
2159 Assert((dynamic_cast<
2160 const typename ::MappingQ<dim, spacedim>::InternalData *>(
2161 &mapping_data) != nullptr),
2163 const typename ::MappingQ<dim, spacedim>::InternalData &data =
2164 static_cast<
2165 const typename ::MappingQ<dim, spacedim>::InternalData &>(
2166 mapping_data);
2167
2168 switch (mapping_kind)
2169 {
2171 {
2174 "update_contravariant_transformation"));
2175
2176 for (unsigned int i = 0; i < output.size(); ++i)
2177 output[i] = apply_transformation(data.output_data->jacobians[i],
2178 input[i]);
2179
2180 return;
2181 }
2182
2183 case mapping_piola:
2184 {
2187 "update_contravariant_transformation"));
2188 Assert(data.update_each & update_volume_elements,
2190 "update_volume_elements"));
2191 Assert(rank == 1, ExcMessage("Only for rank 1"));
2192 if (rank != 1)
2193 return;
2194
2195 for (unsigned int i = 0; i < output.size(); ++i)
2196 {
2197 output[i] =
2198 apply_transformation(data.output_data->jacobians[i],
2199 input[i]);
2200 output[i] /= data.volume_elements[i];
2201 }
2202 return;
2203 }
2204 // We still allow this operation as in the
2205 // reference cell Derivatives are Tensor
2206 // rather than DerivativeForm
2207 case mapping_covariant:
2208 {
2211 "update_covariant_transformation"));
2212
2213 for (unsigned int i = 0; i < output.size(); ++i)
2214 output[i] = apply_transformation(
2215 data.output_data->inverse_jacobians[i].transpose(), input[i]);
2216
2217 return;
2218 }
2219
2220 default:
2222 }
2223 }
2224
2225
2226
2230 template <int dim, int spacedim, int rank>
2231 inline void
2233 const ArrayView<const Tensor<rank, dim>> &input,
2234 const MappingKind mapping_kind,
2235 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
2236 const ArrayView<Tensor<rank, spacedim>> &output)
2237 {
2238 AssertDimension(input.size(), output.size());
2239 Assert((dynamic_cast<
2240 const typename ::MappingQ<dim, spacedim>::InternalData *>(
2241 &mapping_data) != nullptr),
2243 const typename ::MappingQ<dim, spacedim>::InternalData &data =
2244 static_cast<
2245 const typename ::MappingQ<dim, spacedim>::InternalData &>(
2246 mapping_data);
2247
2248 switch (mapping_kind)
2249 {
2251 {
2252 Assert(data.update_each & update_covariant_transformation,
2254 "update_covariant_transformation"));
2257 "update_contravariant_transformation"));
2258 Assert(rank == 2, ExcMessage("Only for rank 2"));
2259
2260 for (unsigned int i = 0; i < output.size(); ++i)
2261 {
2263 apply_transformation(data.output_data->jacobians[i],
2264 transpose(input[i]));
2265 output[i] = apply_transformation(
2266 data.output_data->inverse_jacobians[i].transpose(),
2267 A.transpose());
2268 }
2269
2270 return;
2271 }
2272
2274 {
2275 Assert(data.update_each & update_covariant_transformation,
2277 "update_covariant_transformation"));
2278 Assert(rank == 2, ExcMessage("Only for rank 2"));
2279
2280 for (unsigned int i = 0; i < output.size(); ++i)
2281 {
2282 const DerivativeForm<1, dim, spacedim> covariant =
2283 data.output_data->inverse_jacobians[i].transpose();
2285 apply_transformation(covariant, transpose(input[i]));
2286 output[i] = apply_transformation(covariant, A.transpose());
2287 }
2288
2289 return;
2290 }
2291
2293 {
2294 Assert(data.update_each & update_covariant_transformation,
2296 "update_covariant_transformation"));
2299 "update_contravariant_transformation"));
2300 Assert(data.update_each & update_volume_elements,
2302 "update_volume_elements"));
2303 Assert(rank == 2, ExcMessage("Only for rank 2"));
2304
2305 for (unsigned int i = 0; i < output.size(); ++i)
2306 {
2307 const DerivativeForm<1, dim, spacedim> covariant =
2308 data.output_data->inverse_jacobians[i].transpose();
2310 apply_transformation(covariant, input[i]);
2311 const Tensor<2, spacedim> T =
2312 apply_transformation(data.output_data->jacobians[i],
2313 A.transpose());
2314
2315 output[i] = transpose(T);
2316 output[i] /= data.volume_elements[i];
2317 }
2318
2319 return;
2320 }
2321
2322 default:
2324 }
2325 }
2326
2327
2328
2332 template <int dim, int spacedim>
2333 inline void
2335 const ArrayView<const Tensor<3, dim>> &input,
2336 const MappingKind mapping_kind,
2337 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
2338 const ArrayView<Tensor<3, spacedim>> &output)
2339 {
2340 AssertDimension(input.size(), output.size());
2341 Assert((dynamic_cast<
2342 const typename ::MappingQ<dim, spacedim>::InternalData *>(
2343 &mapping_data) != nullptr),
2345 const typename ::MappingQ<dim, spacedim>::InternalData &data =
2346 static_cast<
2347 const typename ::MappingQ<dim, spacedim>::InternalData &>(
2348 mapping_data);
2349
2350 switch (mapping_kind)
2351 {
2353 {
2354 Assert(data.update_each & update_covariant_transformation,
2356 "update_covariant_transformation"));
2359 "update_contravariant_transformation"));
2360
2361 for (unsigned int q = 0; q < output.size(); ++q)
2362 {
2363 const DerivativeForm<1, dim, spacedim> covariant =
2364 data.output_data->inverse_jacobians[q].transpose();
2365 const DerivativeForm<1, dim, spacedim> contravariant =
2366 data.output_data->jacobians[q];
2367 output[q] =
2369 contravariant,
2370 input[q]);
2371 }
2372
2373 return;
2374 }
2375
2377 {
2378 Assert(data.update_each & update_covariant_transformation,
2380 "update_covariant_transformation"));
2381
2382 for (unsigned int q = 0; q < output.size(); ++q)
2383 {
2384 const DerivativeForm<1, dim, spacedim> covariant =
2385 data.output_data->inverse_jacobians[q].transpose();
2386 output[q] =
2387 internal::apply_covariant_hessian(covariant, input[q]);
2388 }
2389
2390 return;
2391 }
2392
2394 {
2395 Assert(data.update_each & update_covariant_transformation,
2397 "update_covariant_transformation"));
2400 "update_contravariant_transformation"));
2401 Assert(data.update_each & update_volume_elements,
2403 "update_volume_elements"));
2404
2405 for (unsigned int q = 0; q < output.size(); ++q)
2406 {
2407 const DerivativeForm<1, dim, spacedim> covariant =
2408 data.output_data->inverse_jacobians[q].transpose();
2409 const DerivativeForm<1, dim, spacedim> contravariant =
2410 data.output_data->jacobians[q];
2411 output[q] =
2413 contravariant,
2414 data.volume_elements[q],
2415 input[q]);
2416 }
2417
2418 return;
2419 }
2420
2421 default:
2423 }
2424 }
2425
2426
2427
2432 template <int dim, int spacedim, int rank>
2433 inline void
2436 const MappingKind mapping_kind,
2437 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
2439 {
2440 AssertDimension(input.size(), output.size());
2441 Assert((dynamic_cast<
2442 const typename ::MappingQ<dim, spacedim>::InternalData *>(
2443 &mapping_data) != nullptr),
2445 const typename ::MappingQ<dim, spacedim>::InternalData &data =
2446 static_cast<
2447 const typename ::MappingQ<dim, spacedim>::InternalData &>(
2448 mapping_data);
2449
2450 switch (mapping_kind)
2451 {
2452 case mapping_covariant:
2453 {
2456 "update_covariant_transformation"));
2457
2458 for (unsigned int i = 0; i < output.size(); ++i)
2459 output[i] = apply_transformation(
2460 data.output_data->inverse_jacobians[i].transpose(), input[i]);
2461
2462 return;
2463 }
2464 default:
2466 }
2467 }
2468 } // namespace MappingQImplementation
2469} // namespace internal
2470
2472
2473#endif
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
std::size_t size() const
Definition array_view.h:689
Number determinant() const
DerivativeForm< 1, dim, spacedim, Number > covariant_form() const
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
DerivativeForm< 1, spacedim, dim, Number > transpose() const
void set_data_pointers(AlignedVector< Number > *scratch_data, const unsigned int n_components)
const Number * begin_gradients() const
const Number * begin_values() const
const Number * begin_dof_values() const
const Number * begin_hessians() const
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
constexpr numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to...
Definition qprojector.h:334
const std::vector< double > & get_weights() const
unsigned int size() const
numbers::NumberTraits< Number >::real_type norm() const
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
InverseQuadraticApproximation(const InverseQuadraticApproximation &)=default
Point< dim, Number > compute(const Point< spacedim, Number > &p) const
InverseQuadraticApproximation(const ArrayView< const Point< spacedim > > &real_support_points, const std::vector< Point< dim > > &unit_support_points)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Tensor< 1, spacedim, typename ProductType< Number1, Number2 >::type > apply_transformation(const DerivativeForm< 1, dim, spacedim, Number1 > &grad_F, const Tensor< 1, dim, Number2 > &d_x)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcTransformationFailed()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#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)
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_jacobian_3rd_derivatives
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_quadrature_points
Transformed quadrature points.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives
@ mapping_piola
Definition mapping.h:116
@ mapping_covariant_gradient
Definition mapping.h:102
@ mapping_covariant
Definition mapping.h:91
@ mapping_contravariant
Definition mapping.h:96
@ mapping_contravariant_hessian
Definition mapping.h:158
@ mapping_covariant_hessian
Definition mapping.h:152
@ mapping_contravariant_gradient
Definition mapping.h:108
@ mapping_piola_gradient
Definition mapping.h:122
@ mapping_piola_hessian
Definition mapping.h:164
LogStream deallog
Definition logstream.cc:38
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
EvaluationFlags
The EvaluationFlags enum.
std::vector< unsigned int > hierarchic_to_lexicographic_numbering(unsigned int degree)
std::pair< DerivativeForm< 1, dim, spacedim >, Tensor< 1, spacedim > > affine_cell_approximation(const ArrayView< const Point< spacedim > > &vertices)
constexpr T fixed_power(const T t)
Definition utilities.h:943
constexpr T pow(const T base, const int iexp)
Definition utilities.h:967
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
void maybe_update_jacobian_pushed_forward_2nd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 4, spacedim > > &jacobian_pushed_forward_2nd_derivatives)
Point< dim, Number > do_transform_real_to_unit_cell_internal(const Point< spacedim, Number > &p, const Point< dim, Number > &initial_p_unit, const ArrayView< const Point< spacedim > > &points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber, const bool print_iterations_to_deallog=false)
DerivativeForm< 3, dim, spacedim > expand_3rd_derivatives(const Tensor< 1, length_tensor, Tensor< 1, spacedim > > &compressed)
inline ::Table< 2, double > compute_support_point_weights_cell(const unsigned int polynomial_degree)
void transform_differential_forms(const ArrayView< const DerivativeForm< rank, dim, spacedim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank+1, spacedim > > &output)
void maybe_update_jacobian_grads(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads)
void do_fill_fe_face_values(const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const typename QProjector< dim >::DataSetDescriptor data_set, const Quadrature< dim - 1 > &quadrature, const typename ::MappingQ< dim, spacedim >::InternalData &data, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
inline ::Table< 2, double > compute_support_point_weights_on_hex(const unsigned int polynomial_degree)
void transform_fields(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void maybe_update_jacobian_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 4, dim, spacedim > > &jacobian_3rd_derivatives)
void maybe_update_jacobian_pushed_forward_grads(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 3, spacedim > > &jacobian_pushed_forward_grads)
inline ::Table< 2, double > compute_support_point_weights_on_quad(const unsigned int polynomial_degree)
void maybe_update_q_points_Jacobians_generic(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Point< spacedim > > &quadrature_points, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
void transform_gradients(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void maybe_update_q_points_Jacobians_and_grads_tensor(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Point< spacedim > > &quadrature_points, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads)
void transform_hessians(const ArrayView< const Tensor< 3, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< 3, spacedim > > &output)
void maybe_update_jacobian_pushed_forward_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 5, spacedim > > &jacobian_pushed_forward_3rd_derivatives)
std::vector<::Table< 2, double > > compute_support_point_weights_perimeter_to_interior(const unsigned int polynomial_degree, const unsigned int dim)
DerivativeForm< 4, dim, spacedim > expand_4th_derivatives(const Tensor< 1, length_tensor, Tensor< 1, spacedim > > &compressed)
Point< spacedim > compute_mapped_location_of_point(const typename ::MappingQ< dim, spacedim >::InternalData &data)
Point< dim > do_transform_real_to_unit_cell_internal_codim1(const Point< dim+1 > &p, const Point< dim > &initial_p_unit, const ArrayView< const Point< dim+1 > > &points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber)
void maybe_compute_face_data(const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int n_q_points, const std::vector< double > &weights, const typename ::MappingQ< dim, spacedim >::InternalData &data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
void maybe_update_jacobian_2nd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 3, dim, spacedim > > &jacobian_2nd_derivatives)
Tensor< 3, spacedim, Number > apply_contravariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_piola_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Number &volume_element, const Tensor< 3, dim, Number > &input)
SymmetricTensor< 2, dim, typename ProductTypeNoPoint< Number, Number2 >::type > evaluate_tensor_product_hessian(const std::vector< Polynomials::Polynomial< double > > &poly, const ArrayView< const Number > &values, const Point< dim, Number2 > &p, const std::vector< unsigned int > &renumber={})
Tensor< 3, spacedim, Number > apply_covariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const Tensor< 3, dim, Number > &input)
std::pair< typename ProductTypeNoPoint< Number, Number2 >::type, Tensor< 1, dim, typename ProductTypeNoPoint< Number, Number2 >::type > > evaluate_tensor_product_value_and_gradient(const std::vector< Polynomials::Polynomial< double > > &poly, const ArrayView< const Number > &values, const Point< dim, Number2 > &p, const bool d_linear=false, const std::vector< unsigned int > &renumber={})
ProductTypeNoPoint< Number, Number2 >::type evaluate_tensor_product_value(const std::vector< Polynomials::Polynomial< double > > &poly, const ArrayView< const Number > &values, const Point< dim, Number2 > &p, const bool d_linear=false, const std::vector< unsigned int > &renumber={})
Tensor< 1, 1, typename ProductTypeNoPoint< Number, Number2 >::type > evaluate_tensor_product_higher_derivatives(const std::vector< Polynomials::Polynomial< double > > &poly, const ArrayView< const Number > &values, const Point< 1, Number2 > &p, const std::vector< unsigned int > &renumber={})
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
Definition ndarray.h:107
static constexpr unsigned int vertices_per_cell
static constexpr unsigned int faces_per_cell
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
static Point< dim, Number > project_to_unit_cell(const Point< dim, Number > &p)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs, FEEvaluationData< dim, Number, false > &fe_eval)
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition tensor.h:2647
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition tensor.h:2672
Number compare_and_apply_mask(const Number &left, const Number &right, const Number &true_value, const Number &false_value)