Reference documentation for deal.II version 9.3.3
\(\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\}}\)
mapping_q_generic.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
23#include <deal.II/base/table.h>
25
26#include <deal.II/fe/fe_dgq.h>
27#include <deal.II/fe/fe_tools.h>
32
34#include <deal.II/grid/tria.h>
36
38#include <boost/container/small_vector.hpp>
40
41#include <algorithm>
42#include <array>
43#include <cmath>
44#include <memory>
45#include <numeric>
46
47
49
50
51template <int dim, int spacedim>
53 const unsigned int polynomial_degree)
54 : polynomial_degree(polynomial_degree)
55 , n_shape_functions(Utilities::fixed_power<dim>(polynomial_degree + 1))
56 , line_support_points(QGaussLobatto<1>(polynomial_degree + 1))
57 , tensor_product_quadrature(false)
58{}
59
60
61
62template <int dim, int spacedim>
63std::size_t
65{
66 return (
69 MemoryConsumption::memory_consumption(shape_derivatives) +
74 MemoryConsumption::memory_consumption(mapping_support_points) +
75 MemoryConsumption::memory_consumption(cell_of_current_support_points) +
77 MemoryConsumption::memory_consumption(polynomial_degree) +
78 MemoryConsumption::memory_consumption(n_shape_functions));
79}
80
81
82
83template <int dim, int spacedim>
84void
86 const UpdateFlags update_flags,
87 const Quadrature<dim> &q,
88 const unsigned int n_original_q_points)
89{
90 // store the flags in the internal data object so we can access them
91 // in fill_fe_*_values()
92 this->update_each = update_flags;
93
94 const unsigned int n_q_points = q.size();
95
96 const bool needs_higher_order_terms =
97 this->update_each &
102
103 if (this->update_each & update_covariant_transformation)
104 covariant.resize(n_original_q_points);
105
106 if (this->update_each & update_contravariant_transformation)
107 contravariant.resize(n_original_q_points);
108
109 if (this->update_each & update_volume_elements)
110 volume_elements.resize(n_original_q_points);
111
112 tensor_product_quadrature = q.is_tensor_product();
113
114 // use of MatrixFree only for higher order elements and with more than one
115 // point where tensor products do not make sense
116 if (polynomial_degree < 2 || n_q_points == 1)
117 tensor_product_quadrature = false;
118
119 if (dim > 1)
120 {
121 // find out if the one-dimensional formula is the same
122 // in all directions
123 if (tensor_product_quadrature)
124 {
125 const std::array<Quadrature<1>, dim> quad_array =
127 for (unsigned int i = 1; i < dim && tensor_product_quadrature; ++i)
128 {
129 if (quad_array[i - 1].size() != quad_array[i].size())
130 {
131 tensor_product_quadrature = false;
132 break;
133 }
134 else
135 {
136 const std::vector<Point<1>> &points_1 =
137 quad_array[i - 1].get_points();
138 const std::vector<Point<1>> &points_2 =
139 quad_array[i].get_points();
140 const std::vector<double> &weights_1 =
141 quad_array[i - 1].get_weights();
142 const std::vector<double> &weights_2 =
143 quad_array[i].get_weights();
144 for (unsigned int j = 0; j < quad_array[i].size(); ++j)
145 {
146 if (std::abs(points_1[j][0] - points_2[j][0]) > 1.e-10 ||
147 std::abs(weights_1[j] - weights_2[j]) > 1.e-10)
148 {
149 tensor_product_quadrature = false;
150 break;
151 }
152 }
153 }
154 }
155
156 if (tensor_product_quadrature)
157 {
158 // use a 1D FE_DGQ and adjust the hierarchic -> lexicographic
159 // numbering manually (building an FE_Q<dim> is relatively
160 // expensive due to constraints)
161 const FE_DGQ<1> fe(polynomial_degree);
162 shape_info.reinit(q.get_tensor_basis()[0], fe);
163 shape_info.lexicographic_numbering =
164 FETools::lexicographic_to_hierarchic_numbering<dim>(
165 polynomial_degree);
166 shape_info.n_q_points = q.size();
167 shape_info.dofs_per_component_on_cell =
168 Utilities::pow(polynomial_degree + 1, dim);
169 }
170 }
171 }
172
173 // Only fill the big arrays on demand in case we cannot use the tensor
174 // product quadrature code path
175 if (dim == 1 || !tensor_product_quadrature || needs_higher_order_terms)
176 {
177 // see if we need the (transformation) shape function values
178 // and/or gradients and resize the necessary arrays
179 if (this->update_each & update_quadrature_points)
180 shape_values.resize(n_shape_functions * n_q_points);
181
182 if (this->update_each &
192 shape_derivatives.resize(n_shape_functions * n_q_points);
193
194 if (this->update_each &
196 shape_second_derivatives.resize(n_shape_functions * n_q_points);
197
198 if (this->update_each & (update_jacobian_2nd_derivatives |
200 shape_third_derivatives.resize(n_shape_functions * n_q_points);
201
202 if (this->update_each & (update_jacobian_3rd_derivatives |
204 shape_fourth_derivatives.resize(n_shape_functions * n_q_points);
205
206 // now also fill the various fields with their correct values
207 compute_shape_function_values(q.get_points());
208 }
209}
210
211
212
213template <int dim, int spacedim>
214void
216 const UpdateFlags update_flags,
217 const Quadrature<dim> &q,
218 const unsigned int n_original_q_points)
219{
220 initialize(update_flags, q, n_original_q_points);
221
222 if (dim > 1 && tensor_product_quadrature)
223 {
224 constexpr unsigned int facedim = dim - 1;
225 const FE_DGQ<1> fe(polynomial_degree);
226 shape_info.reinit(q.get_tensor_basis()[0], fe);
227 shape_info.lexicographic_numbering =
228 FETools::lexicographic_to_hierarchic_numbering<facedim>(
229 polynomial_degree);
230 shape_info.n_q_points = n_original_q_points;
231 shape_info.dofs_per_component_on_cell =
232 Utilities::pow(polynomial_degree + 1, dim);
233 }
234
235 if (dim > 1)
236 {
237 if (this->update_each &
240 {
241 aux.resize(dim - 1,
242 std::vector<Tensor<1, spacedim>>(n_original_q_points));
243
244 // Compute tangentials to the unit cell.
245 for (const unsigned int i : GeometryInfo<dim>::face_indices())
246 {
247 unit_tangentials[i].resize(n_original_q_points);
248 std::fill(unit_tangentials[i].begin(),
249 unit_tangentials[i].end(),
251 if (dim > 2)
252 {
253 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
254 .resize(n_original_q_points);
255 std::fill(
256 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
257 .begin(),
258 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
259 .end(),
261 }
262 }
263 }
264 }
265}
266
267
268
269template <int dim, int spacedim>
270void
272 const std::vector<Point<dim>> &unit_points)
273{
274 const unsigned int n_points = unit_points.size();
275
276 // Construct the tensor product polynomials used as shape functions for
277 // the Qp mapping of cells at the boundary.
278 const TensorProductPolynomials<dim> tensor_pols(
280 line_support_points.get_points()));
281 Assert(n_shape_functions == tensor_pols.n(), ExcInternalError());
282
283 // then also construct the mapping from lexicographic to the Qp shape
284 // function numbering
285 const std::vector<unsigned int> renumber =
286 FETools::hierarchic_to_lexicographic_numbering<dim>(polynomial_degree);
287
288 std::vector<double> values;
289 std::vector<Tensor<1, dim>> grads;
290 if (shape_values.size() != 0)
291 {
292 Assert(shape_values.size() == n_shape_functions * n_points,
294 values.resize(n_shape_functions);
295 }
296 if (shape_derivatives.size() != 0)
297 {
298 Assert(shape_derivatives.size() == n_shape_functions * n_points,
300 grads.resize(n_shape_functions);
301 }
302
303 std::vector<Tensor<2, dim>> grad2;
304 if (shape_second_derivatives.size() != 0)
305 {
306 Assert(shape_second_derivatives.size() == n_shape_functions * n_points,
308 grad2.resize(n_shape_functions);
309 }
310
311 std::vector<Tensor<3, dim>> grad3;
312 if (shape_third_derivatives.size() != 0)
313 {
314 Assert(shape_third_derivatives.size() == n_shape_functions * n_points,
316 grad3.resize(n_shape_functions);
317 }
318
319 std::vector<Tensor<4, dim>> grad4;
320 if (shape_fourth_derivatives.size() != 0)
321 {
322 Assert(shape_fourth_derivatives.size() == n_shape_functions * n_points,
324 grad4.resize(n_shape_functions);
325 }
326
327
328 if (shape_values.size() != 0 || shape_derivatives.size() != 0 ||
329 shape_second_derivatives.size() != 0 ||
330 shape_third_derivatives.size() != 0 ||
331 shape_fourth_derivatives.size() != 0)
332 for (unsigned int point = 0; point < n_points; ++point)
333 {
334 tensor_pols.evaluate(
335 unit_points[point], values, grads, grad2, grad3, grad4);
336
337 if (shape_values.size() != 0)
338 for (unsigned int i = 0; i < n_shape_functions; ++i)
339 shape(point, i) = values[renumber[i]];
340
341 if (shape_derivatives.size() != 0)
342 for (unsigned int i = 0; i < n_shape_functions; ++i)
343 derivative(point, i) = grads[renumber[i]];
344
345 if (shape_second_derivatives.size() != 0)
346 for (unsigned int i = 0; i < n_shape_functions; ++i)
347 second_derivative(point, i) = grad2[renumber[i]];
348
349 if (shape_third_derivatives.size() != 0)
350 for (unsigned int i = 0; i < n_shape_functions; ++i)
351 third_derivative(point, i) = grad3[renumber[i]];
352
353 if (shape_fourth_derivatives.size() != 0)
354 for (unsigned int i = 0; i < n_shape_functions; ++i)
355 fourth_derivative(point, i) = grad4[renumber[i]];
356 }
357}
358
359
360
361template <int dim, int spacedim>
363 : polynomial_degree(p)
364 , line_support_points(
365 QGaussLobatto<1>(this->polynomial_degree + 1).get_points())
366 , polynomials_1d(
367 Polynomials::generate_complete_Lagrange_basis(line_support_points))
368 , renumber_lexicographic_to_hierarchic(
370 , unit_cell_support_points(
371 internal::MappingQGenericImplementation::unit_support_points<dim>(
372 line_support_points,
373 renumber_lexicographic_to_hierarchic))
374 , support_point_weights_perimeter_to_interior(
375 internal::MappingQGenericImplementation::
377 this->polynomial_degree,
378 dim))
379 , support_point_weights_cell(
380 internal::MappingQGenericImplementation::
381 compute_support_point_weights_cell<dim>(this->polynomial_degree))
382{
383 Assert(p >= 1,
384 ExcMessage("It only makes sense to create polynomial mappings "
385 "with a polynomial degree greater or equal to one."));
386}
387
388
389
390template <int dim, int spacedim>
392 const MappingQGeneric<dim, spacedim> &mapping)
393 : polynomial_degree(mapping.polynomial_degree)
394 , line_support_points(mapping.line_support_points)
395 , polynomials_1d(mapping.polynomials_1d)
396 , renumber_lexicographic_to_hierarchic(
397 mapping.renumber_lexicographic_to_hierarchic)
398 , support_point_weights_perimeter_to_interior(
399 mapping.support_point_weights_perimeter_to_interior)
400 , support_point_weights_cell(mapping.support_point_weights_cell)
401{}
402
403
404
405template <int dim, int spacedim>
406std::unique_ptr<Mapping<dim, spacedim>>
408{
409 return std::make_unique<MappingQGeneric<dim, spacedim>>(*this);
410}
411
412
413
414template <int dim, int spacedim>
415unsigned int
417{
418 return polynomial_degree;
419}
420
421
422
423template <int dim, int spacedim>
427 const Point<dim> & p) const
428{
432 p,
433 polynomials_1d.size() == 2,
435 .first);
436}
437
438
439// In the code below, GCC tries to instantiate MappingQGeneric<3,4> when
440// seeing which of the overloaded versions of
441// do_transform_real_to_unit_cell_internal() to call. This leads to bad
442// error messages and, generally, nothing very good. Avoid this by ensuring
443// that this class exists, but does not have an inner InternalData
444// type, thereby ruling out the codim-1 version of the function
445// below when doing overload resolution.
446template <>
448{};
449
450
451
452// visual studio freaks out when trying to determine if
453// do_transform_real_to_unit_cell_internal with dim=3 and spacedim=4 is a good
454// candidate. So instead of letting the compiler pick the correct overload, we
455// use template specialization to make sure we pick up the right function to
456// call:
457
458template <int dim, int spacedim>
462 const Point<spacedim> &,
463 const Point<dim> &) const
464{
465 // default implementation (should never be called)
466 Assert(false, ExcInternalError());
467 return {};
468}
469
470
471
472template <>
476 const Point<1> & p,
477 const Point<1> & initial_p_unit) const
478{
479 // dispatch to the various specializations for spacedim=dim,
480 // spacedim=dim+1, etc
481 return internal::MappingQGenericImplementation::
482 do_transform_real_to_unit_cell_internal<1>(
483 p,
484 initial_p_unit,
485 this->compute_mapping_support_points(cell),
486 polynomials_1d,
487 renumber_lexicographic_to_hierarchic);
488}
489
490
491
492template <>
496 const Point<2> & p,
497 const Point<2> & initial_p_unit) const
498{
499 return internal::MappingQGenericImplementation::
500 do_transform_real_to_unit_cell_internal<2>(
501 p,
502 initial_p_unit,
503 this->compute_mapping_support_points(cell),
504 polynomials_1d,
505 renumber_lexicographic_to_hierarchic);
506}
507
508
509
510template <>
514 const Point<3> & p,
515 const Point<3> & initial_p_unit) const
516{
517 return internal::MappingQGenericImplementation::
518 do_transform_real_to_unit_cell_internal<3>(
519 p,
520 initial_p_unit,
521 this->compute_mapping_support_points(cell),
522 polynomials_1d,
523 renumber_lexicographic_to_hierarchic);
524}
525
526
527
528template <>
532 const Point<2> & p,
533 const Point<1> & initial_p_unit) const
534{
535 const int dim = 1;
536 const int spacedim = 2;
537
538 const Quadrature<dim> point_quadrature(initial_p_unit);
539
541 if (spacedim > dim)
542 update_flags |= update_jacobian_grads;
543 auto mdata = Utilities::dynamic_unique_cast<InternalData>(
544 get_data(update_flags, point_quadrature));
545
546 mdata->mapping_support_points = this->compute_mapping_support_points(cell);
547
548 // dispatch to the various specializations for spacedim=dim,
549 // spacedim=dim+1, etc
550 return internal::MappingQGenericImplementation::
551 do_transform_real_to_unit_cell_internal_codim1<1>(cell,
552 p,
553 initial_p_unit,
554 *mdata);
555}
556
557
558
559template <>
563 const Point<3> & p,
564 const Point<2> & initial_p_unit) const
565{
566 const int dim = 2;
567 const int spacedim = 3;
568
569 const Quadrature<dim> point_quadrature(initial_p_unit);
570
572 if (spacedim > dim)
573 update_flags |= update_jacobian_grads;
574 auto mdata = Utilities::dynamic_unique_cast<InternalData>(
575 get_data(update_flags, point_quadrature));
576
577 mdata->mapping_support_points = this->compute_mapping_support_points(cell);
578
579 // dispatch to the various specializations for spacedim=dim,
580 // spacedim=dim+1, etc
581 return internal::MappingQGenericImplementation::
582 do_transform_real_to_unit_cell_internal_codim1<2>(cell,
583 p,
584 initial_p_unit,
585 *mdata);
586}
587
588template <>
592 const Point<3> &,
593 const Point<1> &) const
594{
595 Assert(false, ExcNotImplemented());
596 return {};
597}
598
599
600
601template <int dim, int spacedim>
605 const Point<spacedim> & p) const
606{
607 // Use an exact formula if one is available. this is only the case
608 // for Q1 mappings in 1d, and in 2d if dim==spacedim
609 if (this->preserves_vertex_locations() && (polynomial_degree == 1) &&
610 ((dim == 1) || ((dim == 2) && (dim == spacedim))))
611 {
612 // The dimension-dependent algorithms are much faster (about 25-45x in
613 // 2D) but fail most of the time when the given point (p) is not in the
614 // cell. The dimension-independent Newton algorithm given below is
615 // slower, but more robust (though it still sometimes fails). Therefore
616 // this function implements the following strategy based on the
617 // p's dimension:
618 //
619 // * In 1D this mapping is linear, so the mapping is always invertible
620 // (and the exact formula is known) as long as the cell has non-zero
621 // length.
622 // * In 2D the exact (quadratic) formula is called first. If either the
623 // exact formula does not succeed (negative discriminant in the
624 // quadratic formula) or succeeds but finds a solution outside of the
625 // unit cell, then the Newton solver is called. The rationale for the
626 // second choice is that the exact formula may provide two different
627 // answers when mapping a point outside of the real cell, but the
628 // Newton solver (if it converges) will only return one answer.
629 // Otherwise the exact formula successfully found a point in the unit
630 // cell and that value is returned.
631 // * In 3D there is no (known to the authors) exact formula, so the Newton
632 // algorithm is used.
633 const auto vertices_ = this->get_vertices(cell);
634
635 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
636 vertices;
637 for (unsigned int i = 0; i < vertices.size(); ++i)
638 vertices[i] = vertices_[i];
639
640 try
641 {
642 switch (dim)
643 {
644 case 1:
645 {
646 // formula not subject to any issues in 1d
647 if (spacedim == 1)
649 vertices, p);
650 else
651 break;
652 }
653
654 case 2:
655 {
656 const Point<dim> point =
658 p);
659
660 // formula not guaranteed to work for points outside of
661 // the cell. only take the computed point if it lies
662 // inside the reference cell
663 const double eps = 1e-15;
664 if (-eps <= point(1) && point(1) <= 1 + eps &&
665 -eps <= point(0) && point(0) <= 1 + eps)
666 {
667 return point;
668 }
669 else
670 break;
671 }
672
673 default:
674 {
675 // we should get here, based on the if-condition at the top
676 Assert(false, ExcInternalError());
677 }
678 }
679 }
680 catch (
682 {
683 // simply fall through and continue on to the standard Newton code
684 }
685 }
686 else
687 {
688 // we can't use an explicit formula,
689 }
690
691
692 // Find the initial value for the Newton iteration by a normal
693 // projection to the least square plane determined by the vertices
694 // of the cell
695 Point<dim> initial_p_unit;
696 if (this->preserves_vertex_locations())
697 {
698 initial_p_unit = cell->real_to_unit_cell_affine_approximation(p);
699 // in 1d with spacedim > 1 the affine approximation is exact
700 if (dim == 1 && polynomial_degree == 1)
701 return initial_p_unit;
702 }
703 else
704 {
705 // else, we simply use the mid point
706 for (unsigned int d = 0; d < dim; ++d)
707 initial_p_unit[d] = 0.5;
708 }
709
710 // perform the Newton iteration and return the result. note that this
711 // statement may throw an exception, which we simply pass up to the caller
712 const Point<dim> p_unit =
713 this->transform_real_to_unit_cell_internal(cell, p, initial_p_unit);
714 if (p_unit[0] == std::numeric_limits<double>::infinity())
715 AssertThrow(false,
717 return p_unit;
718}
719
720
721
722template <int dim, int spacedim>
723void
726 const ArrayView<const Point<spacedim>> & real_points,
727 const ArrayView<Point<dim>> & unit_points) const
728{
729 // Go to base class functions for dim < spacedim because it is not yet
730 // implemented with optimized code.
731 if (dim < spacedim)
732 {
734 real_points,
735 unit_points);
736 return;
737 }
738
739 AssertDimension(real_points.size(), unit_points.size());
740 const std::vector<Point<spacedim>> support_points =
742
743 // From the given (high-order) support points, now only pick the first
744 // 2^dim points and construct an affine approximation from those.
745 internal::MappingQGenericImplementation::
746 InverseQuadraticApproximation<dim, spacedim>
747 inverse_approximation(support_points, unit_cell_support_points);
748
749 const unsigned int n_points = real_points.size();
750 const unsigned int n_lanes = VectorizedArray<double>::size();
751
752 // Use the more heavy VectorizedArray code path if there is more than
753 // one point left to compute
754 for (unsigned int i = 0; i < n_points; i += n_lanes)
755 if (n_points - i > 1)
756 {
758 for (unsigned int j = 0; j < n_lanes; ++j)
759 if (i + j < n_points)
760 for (unsigned int d = 0; d < spacedim; ++d)
761 p_vec[d][j] = real_points[i + j][d];
762 else
763 for (unsigned int d = 0; d < spacedim; ++d)
764 p_vec[d][j] = real_points[i][d];
765
767 internal::MappingQGenericImplementation::
768 do_transform_real_to_unit_cell_internal<dim, spacedim>(
769 p_vec,
770 inverse_approximation.compute(p_vec),
771 support_points,
774
775 // If the vectorized computation failed, it could be that only some of
776 // the lanes failed but others would have succeeded if we had let them
777 // compute alone without interference (like negative Jacobian
778 // determinants) from other SIMD lanes. Repeat the computation in this
779 // unlikely case with scalar arguments.
780 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
781 if (unit_point[0][j] == std::numeric_limits<double>::infinity())
782 unit_points[i + j] = internal::MappingQGenericImplementation::
783 do_transform_real_to_unit_cell_internal<dim, spacedim>(
784 real_points[i + j],
785 inverse_approximation.compute(real_points[i + j]),
786 support_points,
789 else
790 for (unsigned int d = 0; d < dim; ++d)
791 unit_points[i + j][d] = unit_point[d][j];
792 }
793 else
794 unit_points[i] = internal::MappingQGenericImplementation::
795 do_transform_real_to_unit_cell_internal<dim, spacedim>(
796 real_points[i],
797 inverse_approximation.compute(real_points[i]),
798 support_points,
801}
802
803
804
805template <int dim, int spacedim>
808 const UpdateFlags in) const
809{
810 // add flags if the respective quantities are necessary to compute
811 // what we need. note that some flags appear in both the conditions
812 // and in subsequent set operations. this leads to some circular
813 // logic. the only way to treat this is to iterate. since there are
814 // 5 if-clauses in the loop, it will take at most 5 iterations to
815 // converge. do them:
816 UpdateFlags out = in;
817 for (unsigned int i = 0; i < 5; ++i)
818 {
819 // The following is a little incorrect:
820 // If not applied on a face,
821 // update_boundary_forms does not
822 // make sense. On the other hand,
823 // it is necessary on a
824 // face. Currently,
825 // update_boundary_forms is simply
826 // ignored for the interior of a
827 // cell.
830
835
836 if (out &
841
842 // The contravariant transformation is used in the Piola
843 // transformation, which requires the determinant of the Jacobi
844 // matrix of the transformation. Because we have no way of
845 // knowing here whether the finite element wants to use the
846 // contravariant or the Piola transforms, we add the JxW values
847 // to the list of flags to be updated for each cell.
850
851 // the same is true when computing normal vectors: they require
852 // the determinant of the Jacobian
853 if (out & update_normal_vectors)
855 }
856
857 return out;
858}
859
860
861
862template <int dim, int spacedim>
863std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
865 const Quadrature<dim> &q) const
866{
867 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
868 std::make_unique<InternalData>(polynomial_degree);
869 auto &data = dynamic_cast<InternalData &>(*data_ptr);
870 data.initialize(this->requires_update_flags(update_flags), q, q.size());
871
872 return data_ptr;
873}
874
875
876
877template <int dim, int spacedim>
878std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
880 const UpdateFlags update_flags,
881 const hp::QCollection<dim - 1> &quadrature) const
882{
883 AssertDimension(quadrature.size(), 1);
884
885 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
886 std::make_unique<InternalData>(polynomial_degree);
887 auto &data = dynamic_cast<InternalData &>(*data_ptr);
888 data.initialize_face(this->requires_update_flags(update_flags),
890 ReferenceCells::get_hypercube<dim>(), quadrature[0]),
891 quadrature[0].size());
892
893 return data_ptr;
894}
895
896
897
898template <int dim, int spacedim>
899std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
901 const UpdateFlags update_flags,
902 const Quadrature<dim - 1> &quadrature) const
903{
904 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
905 std::make_unique<InternalData>(polynomial_degree);
906 auto &data = dynamic_cast<InternalData &>(*data_ptr);
907 data.initialize_face(this->requires_update_flags(update_flags),
909 ReferenceCells::get_hypercube<dim>(), quadrature),
910 quadrature.size());
911
912 return data_ptr;
913}
914
915
916
917template <int dim, int spacedim>
921 const CellSimilarity::Similarity cell_similarity,
922 const Quadrature<dim> & quadrature,
923 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
925 &output_data) const
926{
927 // ensure that the following static_cast is really correct:
928 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
930 const InternalData &data = static_cast<const InternalData &>(internal_data);
931
932 const unsigned int n_q_points = quadrature.size();
933
934 // recompute the support points of the transformation of this
935 // cell. we tried to be clever here in an earlier version of the
936 // library by checking whether the cell is the same as the one we
937 // had visited last, but it turns out to be difficult to determine
938 // that because a cell for the purposes of a mapping is
939 // characterized not just by its (triangulation, level, index)
940 // triple, but also by the locations of its vertices, the manifold
941 // object attached to the cell and all of its bounding faces/edges,
942 // etc. to reliably test that the "cell" we are on is, therefore,
943 // not easily done
944 data.mapping_support_points = this->compute_mapping_support_points(cell);
945 data.cell_of_current_support_points = cell;
946
947 // if the order of the mapping is greater than 1, then do not reuse any cell
948 // similarity information. This is necessary because the cell similarity
949 // value is computed with just cell vertices and does not take into account
950 // cell curvature.
951 const CellSimilarity::Similarity computed_cell_similarity =
952 (polynomial_degree == 1 ? cell_similarity : CellSimilarity::none);
953
954 if (dim > 1 && data.tensor_product_quadrature)
955 {
956 internal::MappingQGenericImplementation::
957 maybe_update_q_points_Jacobians_and_grads_tensor<dim, spacedim>(
958 computed_cell_similarity,
959 data,
960 output_data.quadrature_points,
961 output_data.jacobian_grads);
962 }
963 else
964 {
966 spacedim>(
968 data,
969 output_data.quadrature_points);
970
972 spacedim>(
973 computed_cell_similarity,
975 data);
976
978 dim,
979 spacedim>(computed_cell_similarity,
981 data,
982 output_data.jacobian_grads);
983 }
984
985 internal::MappingQGenericImplementation::
986 maybe_update_jacobian_pushed_forward_grads<dim, spacedim>(
987 computed_cell_similarity,
989 data,
991
992 internal::MappingQGenericImplementation::
993 maybe_update_jacobian_2nd_derivatives<dim, spacedim>(
994 computed_cell_similarity,
996 data,
997 output_data.jacobian_2nd_derivatives);
998
999 internal::MappingQGenericImplementation::
1000 maybe_update_jacobian_pushed_forward_2nd_derivatives<dim, spacedim>(
1001 computed_cell_similarity,
1003 data,
1005
1006 internal::MappingQGenericImplementation::
1007 maybe_update_jacobian_3rd_derivatives<dim, spacedim>(
1008 computed_cell_similarity,
1010 data,
1011 output_data.jacobian_3rd_derivatives);
1012
1013 internal::MappingQGenericImplementation::
1014 maybe_update_jacobian_pushed_forward_3rd_derivatives<dim, spacedim>(
1015 computed_cell_similarity,
1017 data,
1019
1020 const UpdateFlags update_flags = data.update_each;
1021 const std::vector<double> &weights = quadrature.get_weights();
1022
1023 // Multiply quadrature weights by absolute value of Jacobian determinants or
1024 // the area element g=sqrt(DX^t DX) in case of codim > 0
1025
1026 if (update_flags & (update_normal_vectors | update_JxW_values))
1027 {
1028 AssertDimension(output_data.JxW_values.size(), n_q_points);
1029
1030 Assert(!(update_flags & update_normal_vectors) ||
1031 (output_data.normal_vectors.size() == n_q_points),
1032 ExcDimensionMismatch(output_data.normal_vectors.size(),
1033 n_q_points));
1034
1035
1036 if (computed_cell_similarity != CellSimilarity::translation)
1037 for (unsigned int point = 0; point < n_q_points; ++point)
1038 {
1039 if (dim == spacedim)
1040 {
1041 const double det = data.contravariant[point].determinant();
1042
1043 // check for distorted cells.
1044
1045 // TODO: this allows for anisotropies of up to 1e6 in 3D and
1046 // 1e12 in 2D. might want to find a finer
1047 // (dimension-independent) criterion
1048 Assert(det >
1049 1e-12 * Utilities::fixed_power<dim>(
1050 cell->diameter() / std::sqrt(double(dim))),
1052 cell->center(), det, point)));
1053
1054 output_data.JxW_values[point] = weights[point] * det;
1055 }
1056 // if dim==spacedim, then there is no cell normal to
1057 // compute. since this is for FEValues (and not FEFaceValues),
1058 // there are also no face normals to compute
1059 else // codim>0 case
1060 {
1061 Tensor<1, spacedim> DX_t[dim];
1062 for (unsigned int i = 0; i < spacedim; ++i)
1063 for (unsigned int j = 0; j < dim; ++j)
1064 DX_t[j][i] = data.contravariant[point][i][j];
1065
1066 Tensor<2, dim> G; // First fundamental form
1067 for (unsigned int i = 0; i < dim; ++i)
1068 for (unsigned int j = 0; j < dim; ++j)
1069 G[i][j] = DX_t[i] * DX_t[j];
1070
1071 output_data.JxW_values[point] =
1072 std::sqrt(determinant(G)) * weights[point];
1073
1074 if (computed_cell_similarity ==
1076 {
1077 // we only need to flip the normal
1078 if (update_flags & update_normal_vectors)
1079 output_data.normal_vectors[point] *= -1.;
1080 }
1081 else
1082 {
1083 if (update_flags & update_normal_vectors)
1084 {
1085 Assert(spacedim == dim + 1,
1086 ExcMessage(
1087 "There is no (unique) cell normal for " +
1089 "-dimensional cells in " +
1090 Utilities::int_to_string(spacedim) +
1091 "-dimensional space. This only works if the "
1092 "space dimension is one greater than the "
1093 "dimensionality of the mesh cells."));
1094
1095 if (dim == 1)
1096 output_data.normal_vectors[point] =
1097 cross_product_2d(-DX_t[0]);
1098 else // dim == 2
1099 output_data.normal_vectors[point] =
1100 cross_product_3d(DX_t[0], DX_t[1]);
1101
1102 output_data.normal_vectors[point] /=
1103 output_data.normal_vectors[point].norm();
1104
1105 if (cell->direction_flag() == false)
1106 output_data.normal_vectors[point] *= -1.;
1107 }
1108 }
1109 } // codim>0 case
1110 }
1111 }
1112
1113
1114
1115 // copy values from InternalData to vector given by reference
1116 if (update_flags & update_jacobians)
1117 {
1118 AssertDimension(output_data.jacobians.size(), n_q_points);
1119 if (computed_cell_similarity != CellSimilarity::translation)
1120 for (unsigned int point = 0; point < n_q_points; ++point)
1121 output_data.jacobians[point] = data.contravariant[point];
1122 }
1123
1124 // copy values from InternalData to vector given by reference
1125 if (update_flags & update_inverse_jacobians)
1126 {
1127 AssertDimension(output_data.inverse_jacobians.size(), n_q_points);
1128 if (computed_cell_similarity != CellSimilarity::translation)
1129 for (unsigned int point = 0; point < n_q_points; ++point)
1130 output_data.inverse_jacobians[point] =
1131 data.covariant[point].transpose();
1132 }
1133
1134 return computed_cell_similarity;
1135}
1136
1137
1138
1139template <int dim, int spacedim>
1140void
1143 const unsigned int face_no,
1144 const hp::QCollection<dim - 1> & quadrature,
1145 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
1147 &output_data) const
1148{
1149 AssertDimension(quadrature.size(), 1);
1150
1151 // ensure that the following cast is really correct:
1152 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1154 const InternalData &data = static_cast<const InternalData &>(internal_data);
1155
1156 // if necessary, recompute the support points of the transformation of this
1157 // cell (note that we need to first check the triangulation pointer, since
1158 // otherwise the second test might trigger an exception if the triangulations
1159 // are not the same)
1160 if ((data.mapping_support_points.size() == 0) ||
1161 (&cell->get_triangulation() !=
1162 &data.cell_of_current_support_points->get_triangulation()) ||
1163 (cell != data.cell_of_current_support_points))
1164 {
1165 data.mapping_support_points = this->compute_mapping_support_points(cell);
1166 data.cell_of_current_support_points = cell;
1167 }
1168
1170 *this,
1171 cell,
1172 face_no,
1175 ReferenceCells::get_hypercube<dim>(),
1176 face_no,
1177 cell->face_orientation(face_no),
1178 cell->face_flip(face_no),
1179 cell->face_rotation(face_no),
1180 quadrature[0].size()),
1181 quadrature[0],
1182 data,
1183 output_data);
1184}
1185
1186
1187
1188template <int dim, int spacedim>
1189void
1192 const unsigned int face_no,
1193 const unsigned int subface_no,
1194 const Quadrature<dim - 1> & quadrature,
1195 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
1197 &output_data) const
1198{
1199 // ensure that the following cast is really correct:
1200 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1202 const InternalData &data = static_cast<const InternalData &>(internal_data);
1203
1204 // if necessary, recompute the support points of the transformation of this
1205 // cell (note that we need to first check the triangulation pointer, since
1206 // otherwise the second test might trigger an exception if the triangulations
1207 // are not the same)
1208 if ((data.mapping_support_points.size() == 0) ||
1209 (&cell->get_triangulation() !=
1210 &data.cell_of_current_support_points->get_triangulation()) ||
1211 (cell != data.cell_of_current_support_points))
1212 {
1213 data.mapping_support_points = this->compute_mapping_support_points(cell);
1214 data.cell_of_current_support_points = cell;
1215 }
1216
1218 *this,
1219 cell,
1220 face_no,
1221 subface_no,
1223 ReferenceCells::get_hypercube<dim>(),
1224 face_no,
1225 subface_no,
1226 cell->face_orientation(face_no),
1227 cell->face_flip(face_no),
1228 cell->face_rotation(face_no),
1229 quadrature.size(),
1230 cell->subface_case(face_no)),
1231 quadrature,
1232 data,
1233 output_data);
1234}
1235
1236
1237
1238template <int dim, int spacedim>
1239inline void
1242 const ArrayView<const Point<dim>> & unit_points,
1243 const UpdateFlags update_flags,
1245 &output_data) const
1246{
1247 if (update_flags == update_default)
1248 return;
1249
1250 Assert(update_flags & update_inverse_jacobians ||
1251 update_flags & update_jacobians ||
1252 update_flags & update_quadrature_points,
1254
1255 output_data.initialize(unit_points.size(), update_flags);
1256 const std::vector<Point<spacedim>> support_points =
1258
1259 const unsigned int n_points = unit_points.size();
1260 const unsigned int n_lanes = VectorizedArray<double>::size();
1261
1262 // Use the more heavy VectorizedArray code path if there is more than
1263 // one point left to compute
1264 for (unsigned int i = 0; i < n_points; i += n_lanes)
1265 if (n_points - i > 1)
1266 {
1268 for (unsigned int j = 0; j < n_lanes; ++j)
1269 if (i + j < n_points)
1270 for (unsigned int d = 0; d < dim; ++d)
1271 p_vec[d][j] = unit_points[i + j][d];
1272 else
1273 for (unsigned int d = 0; d < dim; ++d)
1274 p_vec[d][j] = unit_points[i][d];
1275
1276 const auto result =
1279 support_points,
1280 p_vec,
1281 polynomial_degree == 1,
1283
1284 if (update_flags & update_quadrature_points)
1285 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1286 for (unsigned int d = 0; d < spacedim; ++d)
1287 output_data.quadrature_points[i + j][d] = result.first[d][j];
1288
1289 if (update_flags & update_jacobians)
1290 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1291 for (unsigned int d = 0; d < spacedim; ++d)
1292 for (unsigned int e = 0; e < dim; ++e)
1293 output_data.jacobians[i + j][d][e] = result.second[e][d][j];
1294
1295 if (update_flags & update_inverse_jacobians)
1296 {
1298 result.second);
1300 inv_jac = jac.covariant_form();
1301 for (unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
1302 for (unsigned int d = 0; d < dim; ++d)
1303 for (unsigned int e = 0; e < spacedim; ++e)
1304 output_data.inverse_jacobians[i + j][d][e] = inv_jac[d][e][j];
1305 }
1306 }
1307 else
1308 {
1309 const auto result =
1312 support_points,
1313 unit_points[i],
1314 polynomial_degree == 1,
1316
1317 if (update_flags & update_quadrature_points)
1318 output_data.quadrature_points[i] = result.first;
1319
1320 if (update_flags & update_jacobians)
1321 {
1322 DerivativeForm<1, spacedim, dim> jac = result.second;
1323 output_data.jacobians[i] = jac.transpose();
1324 }
1325
1326 if (update_flags & update_inverse_jacobians)
1327 {
1328 DerivativeForm<1, spacedim, dim> jac(result.second);
1330 for (unsigned int d = 0; d < dim; ++d)
1331 for (unsigned int e = 0; e < spacedim; ++e)
1332 output_data.inverse_jacobians[i][d][e] = inv_jac[d][e];
1333 }
1334 }
1335}
1336
1337
1338
1339template <int dim, int spacedim>
1340void
1342 const ArrayView<const Tensor<1, dim>> & input,
1343 const MappingKind mapping_kind,
1344 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1345 const ArrayView<Tensor<1, spacedim>> & output) const
1346{
1348 mapping_kind,
1349 mapping_data,
1350 output);
1351}
1352
1353
1354
1355template <int dim, int spacedim>
1356void
1358 const ArrayView<const DerivativeForm<1, dim, spacedim>> &input,
1359 const MappingKind mapping_kind,
1360 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1361 const ArrayView<Tensor<2, spacedim>> & output) const
1362{
1364 input, mapping_kind, mapping_data, output);
1365}
1366
1367
1368
1369template <int dim, int spacedim>
1370void
1372 const ArrayView<const Tensor<2, dim>> & input,
1373 const MappingKind mapping_kind,
1374 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1375 const ArrayView<Tensor<2, spacedim>> & output) const
1376{
1377 switch (mapping_kind)
1378 {
1381 mapping_kind,
1382 mapping_data,
1383 output);
1384 return;
1385
1390 input, mapping_kind, mapping_data, output);
1391 return;
1392 default:
1393 Assert(false, ExcNotImplemented());
1394 }
1395}
1396
1397
1398
1399template <int dim, int spacedim>
1400void
1402 const ArrayView<const DerivativeForm<2, dim, spacedim>> &input,
1403 const MappingKind mapping_kind,
1404 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1405 const ArrayView<Tensor<3, spacedim>> & output) const
1406{
1407 AssertDimension(input.size(), output.size());
1408 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
1410 const InternalData &data = static_cast<const InternalData &>(mapping_data);
1411
1412 switch (mapping_kind)
1413 {
1415 {
1418 "update_covariant_transformation"));
1419
1420 for (unsigned int q = 0; q < output.size(); ++q)
1421 for (unsigned int i = 0; i < spacedim; ++i)
1422 for (unsigned int j = 0; j < spacedim; ++j)
1423 {
1424 double tmp[dim];
1425 for (unsigned int K = 0; K < dim; ++K)
1426 {
1427 tmp[K] = data.covariant[q][j][0] * input[q][i][0][K];
1428 for (unsigned int J = 1; J < dim; ++J)
1429 tmp[K] += data.covariant[q][j][J] * input[q][i][J][K];
1430 }
1431 for (unsigned int k = 0; k < spacedim; ++k)
1432 {
1433 output[q][i][j][k] = data.covariant[q][k][0] * tmp[0];
1434 for (unsigned int K = 1; K < dim; ++K)
1435 output[q][i][j][k] += data.covariant[q][k][K] * tmp[K];
1436 }
1437 }
1438 return;
1439 }
1440
1441 default:
1442 Assert(false, ExcNotImplemented());
1443 }
1444}
1445
1446
1447
1448template <int dim, int spacedim>
1449void
1451 const ArrayView<const Tensor<3, dim>> & input,
1452 const MappingKind mapping_kind,
1453 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1454 const ArrayView<Tensor<3, spacedim>> & output) const
1455{
1456 switch (mapping_kind)
1457 {
1462 input, mapping_kind, mapping_data, output);
1463 return;
1464 default:
1465 Assert(false, ExcNotImplemented());
1466 }
1467}
1468
1469
1470
1471template <int dim, int spacedim>
1472void
1475 std::vector<Point<spacedim>> & a) const
1476{
1477 // if we only need the midpoint, then ask for it.
1478 if (this->polynomial_degree == 2)
1479 {
1480 for (unsigned int line_no = 0;
1481 line_no < GeometryInfo<dim>::lines_per_cell;
1482 ++line_no)
1483 {
1485 (dim == 1 ?
1486 static_cast<
1488 cell->line(line_no));
1489
1490 const Manifold<dim, spacedim> &manifold =
1491 ((line->manifold_id() == numbers::flat_manifold_id) &&
1492 (dim < spacedim) ?
1493 cell->get_manifold() :
1494 line->get_manifold());
1495 a.push_back(manifold.get_new_point_on_line(line));
1496 }
1497 }
1498 else
1499 // otherwise call the more complicated functions and ask for inner points
1500 // from the manifold description
1501 {
1502 std::vector<Point<spacedim>> tmp_points;
1503 for (unsigned int line_no = 0;
1504 line_no < GeometryInfo<dim>::lines_per_cell;
1505 ++line_no)
1506 {
1508 (dim == 1 ?
1509 static_cast<
1511 cell->line(line_no));
1512
1513 const Manifold<dim, spacedim> &manifold =
1514 ((line->manifold_id() == numbers::flat_manifold_id) &&
1515 (dim < spacedim) ?
1516 cell->get_manifold() :
1517 line->get_manifold());
1518
1519 const std::array<Point<spacedim>, 2> vertices{
1520 {cell->vertex(GeometryInfo<dim>::line_to_cell_vertices(line_no, 0)),
1521 cell->vertex(
1523
1524 const std::size_t n_rows =
1526 a.resize(a.size() + n_rows);
1527 auto a_view = make_array_view(a.end() - n_rows, a.end());
1528 manifold.get_new_points(
1529 make_array_view(vertices.begin(), vertices.end()),
1531 a_view);
1532 }
1533 }
1534}
1535
1536
1537
1538template <>
1539void
1542 std::vector<Point<3>> & a) const
1543{
1544 const unsigned int faces_per_cell = GeometryInfo<3>::faces_per_cell;
1545
1546 // used if face quad at boundary or entirely in the interior of the domain
1547 std::vector<Point<3>> tmp_points;
1548
1549 // loop over all faces and collect points on them
1550 for (unsigned int face_no = 0; face_no < faces_per_cell; ++face_no)
1551 {
1552 const Triangulation<3>::face_iterator face = cell->face(face_no);
1553
1554#ifdef DEBUG
1555 const bool face_orientation = cell->face_orientation(face_no),
1556 face_flip = cell->face_flip(face_no),
1557 face_rotation = cell->face_rotation(face_no);
1558 const unsigned int vertices_per_face = GeometryInfo<3>::vertices_per_face,
1559 lines_per_face = GeometryInfo<3>::lines_per_face;
1560
1561 // some sanity checks up front
1562 for (unsigned int i = 0; i < vertices_per_face; ++i)
1563 Assert(face->vertex_index(i) ==
1564 cell->vertex_index(GeometryInfo<3>::face_to_cell_vertices(
1565 face_no, i, face_orientation, face_flip, face_rotation)),
1567
1568 // indices of the lines that bound a face are given by GeometryInfo<3>::
1569 // face_to_cell_lines
1570 for (unsigned int i = 0; i < lines_per_face; ++i)
1571 Assert(face->line(i) ==
1573 face_no, i, face_orientation, face_flip, face_rotation)),
1575#endif
1576 // extract the points surrounding a quad from the points
1577 // already computed. First get the 4 vertices and then the points on
1578 // the four lines
1579 boost::container::small_vector<Point<3>, 200> tmp_points(
1581 GeometryInfo<2>::lines_per_cell * (polynomial_degree - 1));
1582 for (const unsigned int v : GeometryInfo<2>::vertex_indices())
1583 tmp_points[v] = a[GeometryInfo<3>::face_to_cell_vertices(face_no, v)];
1584 if (polynomial_degree > 1)
1585 for (unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1586 ++line)
1587 for (unsigned int i = 0; i < polynomial_degree - 1; ++i)
1588 tmp_points[4 + line * (polynomial_degree - 1) + i] =
1590 (polynomial_degree - 1) *
1592 i];
1593
1594 const std::size_t n_rows =
1595 support_point_weights_perimeter_to_interior[1].size(0);
1596 a.resize(a.size() + n_rows);
1597 auto a_view = make_array_view(a.end() - n_rows, a.end());
1598 face->get_manifold().get_new_points(
1599 make_array_view(tmp_points.begin(), tmp_points.end()),
1600 support_point_weights_perimeter_to_interior[1],
1601 a_view);
1602 }
1603}
1604
1605
1606
1607template <>
1608void
1611 std::vector<Point<3>> & a) const
1612{
1613 std::array<Point<3>, GeometryInfo<2>::vertices_per_cell> vertices;
1614 for (const unsigned int i : GeometryInfo<2>::vertex_indices())
1615 vertices[i] = cell->vertex(i);
1616
1617 Table<2, double> weights(Utilities::fixed_power<2>(polynomial_degree - 1),
1619 for (unsigned int q = 0, q2 = 0; q2 < polynomial_degree - 1; ++q2)
1620 for (unsigned int q1 = 0; q1 < polynomial_degree - 1; ++q1, ++q)
1621 {
1622 Point<2> point(line_support_points[q1 + 1][0],
1623 line_support_points[q2 + 1][0]);
1624 for (const unsigned int i : GeometryInfo<2>::vertex_indices())
1626 }
1627
1628 const std::size_t n_rows = weights.size(0);
1629 a.resize(a.size() + n_rows);
1630 auto a_view = make_array_view(a.end() - n_rows, a.end());
1631 cell->get_manifold().get_new_points(
1632 make_array_view(vertices.begin(), vertices.end()), weights, a_view);
1633}
1634
1635
1636
1637template <int dim, int spacedim>
1638void
1641 std::vector<Point<spacedim>> &) const
1642{
1643 Assert(false, ExcInternalError());
1644}
1645
1646
1647
1648template <int dim, int spacedim>
1649std::vector<Point<spacedim>>
1651 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
1652{
1653 // get the vertices first
1654 std::vector<Point<spacedim>> a;
1655 a.reserve(Utilities::fixed_power<dim>(polynomial_degree + 1));
1656 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
1657 a.push_back(cell->vertex(i));
1658
1659 if (this->polynomial_degree > 1)
1660 {
1661 // check if all entities have the same manifold id which is when we can
1662 // simply ask the manifold for all points. the transfinite manifold can
1663 // do the interpolation better than this class, so if we detect that we
1664 // do not have to change anything here
1665 Assert(dim <= 3, ExcImpossibleInDim(dim));
1666 bool all_manifold_ids_are_equal = (dim == spacedim);
1667 if (all_manifold_ids_are_equal &&
1669 &cell->get_manifold()) == nullptr)
1670 {
1671 for (auto f : GeometryInfo<dim>::face_indices())
1672 if (&cell->face(f)->get_manifold() != &cell->get_manifold())
1673 all_manifold_ids_are_equal = false;
1674
1675 if (dim == 3)
1676 for (unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
1677 if (&cell->line(l)->get_manifold() != &cell->get_manifold())
1678 all_manifold_ids_are_equal = false;
1679 }
1680
1681 if (all_manifold_ids_are_equal)
1682 {
1683 const std::size_t n_rows = support_point_weights_cell.size(0);
1684 a.resize(a.size() + n_rows);
1685 auto a_view = make_array_view(a.end() - n_rows, a.end());
1687 a.end() - n_rows),
1689 a_view);
1690 }
1691 else
1692 switch (dim)
1693 {
1694 case 1:
1695 add_line_support_points(cell, a);
1696 break;
1697 case 2:
1698 // in 2d, add the points on the four bounding lines to the
1699 // exterior (outer) points
1700 add_line_support_points(cell, a);
1701
1702 // then get the interior support points
1703 if (dim != spacedim)
1704 add_quad_support_points(cell, a);
1705 else
1706 {
1707 const std::size_t n_rows =
1709 a.resize(a.size() + n_rows);
1710 auto a_view = make_array_view(a.end() - n_rows, a.end());
1712 make_array_view(a.begin(), a.end() - n_rows),
1714 a_view);
1715 }
1716 break;
1717
1718 case 3:
1719 // in 3d also add the points located on the boundary faces
1720 add_line_support_points(cell, a);
1721 add_quad_support_points(cell, a);
1722
1723 // then compute the interior points
1724 {
1725 const std::size_t n_rows =
1727 a.resize(a.size() + n_rows);
1728 auto a_view = make_array_view(a.end() - n_rows, a.end());
1730 make_array_view(a.begin(), a.end() - n_rows),
1732 a_view);
1733 }
1734 break;
1735
1736 default:
1737 Assert(false, ExcNotImplemented());
1738 break;
1739 }
1740 }
1741
1742 return a;
1743}
1744
1745
1746
1747template <int dim, int spacedim>
1750 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
1751{
1753}
1754
1755
1756
1757template <int dim, int spacedim>
1758bool
1760 const ReferenceCell &reference_cell) const
1761{
1762 Assert(dim == reference_cell.get_dimension(),
1763 ExcMessage("The dimension of your mapping (" +
1765 ") and the reference cell cell_type (" +
1766 Utilities::to_string(reference_cell.get_dimension()) +
1767 " ) do not agree."));
1768
1769 return reference_cell.is_hyper_cube();
1770}
1771
1772
1773
1774//--------------------------- Explicit instantiations -----------------------
1775#include "mapping_q_generic.inst"
1776
1777
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:697
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:697
DerivativeForm< 1, dim, spacedim, Number > covariant_form() const
DerivativeForm< 1, spacedim, dim, Number > transpose() const
Definition: fe_dgq.h:111
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const
virtual std::size_t memory_consumption() const override
void compute_shape_function_values(const std::vector< Point< dim > > &unit_points)
InternalData(const unsigned int polynomial_degree)
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
unsigned int get_degree() const
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
const std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
const Table< 2, double > support_point_weights_cell
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
const std::vector< Point< dim > > unit_cell_support_points
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void fill_mapping_data_for_generic_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points, const UpdateFlags update_flags, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const override
const unsigned int polynomial_degree
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const override
virtual bool preserves_vertex_locations() const override
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
MappingQGeneric(const unsigned int polynomial_degree)
const std::vector< unsigned int > renumber_lexicographic_to_hierarchic
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
const std::vector< Polynomials::Polynomial< double > > polynomials_1d
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Abstract base class for mapping classes.
Definition: mapping.h:304
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: point.h:111
static DataSetDescriptor cell()
Definition: qprojector.h:563
bool is_tensor_product() const
const std::array< Quadrature< 1 >, dim > & get_tensor_basis() const
Definition: quadrature.cc:323
const std::vector< double > & get_weights() const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
Triangulation< dim, spacedim > & get_triangulation()
unsigned int size() const
Definition: collection.h:109
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
void initialize(const unsigned int n_quadrature_points, const UpdateFlags flags)
Definition: fe_values.cc:2871
std::vector< Tensor< 1, spacedim > > normal_vectors
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< Point< spacedim > > quadrature_points
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
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:416
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:454
Point< 3 > vertices[4]
UpdateFlags
@ 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_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.
@ update_default
No update.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
MappingKind
Definition: mapping.h:65
@ mapping_covariant_gradient
Definition: mapping.h:86
@ mapping_contravariant
Definition: mapping.h:80
@ mapping_contravariant_hessian
Definition: mapping.h:142
@ mapping_covariant_hessian
Definition: mapping.h:136
@ mapping_contravariant_gradient
Definition: mapping.h:92
@ mapping_piola_gradient
Definition: mapping.h:106
@ mapping_piola_hessian
Definition: mapping.h:148
std::vector< unsigned int > lexicographic_to_hierarchic_numbering(unsigned int degree)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
Definition: polynomial.cc:701
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
constexpr T pow(const T base, const int iexp)
Definition: utilities.h:461
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:482
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:473
T fixed_power(const T t)
Definition: utilities.h:1081
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
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 QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQGeneric< dim, spacedim >::InternalData &data, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads)
std::vector<::Table< 2, double > > compute_support_point_weights_perimeter_to_interior(const unsigned int polynomial_degree, const unsigned int dim)
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_Jacobians(const CellSimilarity::Similarity cell_similarity, const typename ::QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQGeneric< dim, spacedim >::InternalData &data)
void do_fill_fe_face_values(const ::MappingQGeneric< 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 ::MappingQGeneric< dim, spacedim >::InternalData &data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
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 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)
inline ::Table< 2, double > compute_support_point_weights_cell(const unsigned int polynomial_degree)
void maybe_compute_q_points(const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQGeneric< dim, spacedim >::InternalData &data, std::vector< Point< spacedim > > &quadrature_points)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
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 std::vector< Number > &values, const Point< dim, Number2 > &p, const bool d_linear=false, const std::vector< unsigned int > &renumber={})
static const unsigned int invalid_unsigned_int
Definition: types.h:196
const types::manifold_id flat_manifold_id
Definition: types.h:264
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition: tensor.h:2539
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:2564