deal.II version GIT relicensing-2206-gaa53ff9447 2024-12-02 09:10:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
mapping_manifold.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
23
25
26#include <deal.II/fe/fe.h>
27#include <deal.II/fe/fe_tools.h>
30
32#include <deal.II/grid/tria.h>
34
36
37#include <algorithm>
38#include <array>
39#include <cmath>
40#include <memory>
41#include <numeric>
42
43
45
46template <int dim, int spacedim>
47std::size_t
64
65
66
67template <int dim, int spacedim>
68void
70 const UpdateFlags update_flags,
71 const Quadrature<dim> &q)
72{
73 // store the flags in the internal data object so we can access them
74 // in fill_fe_*_values()
75 this->update_each = update_flags;
76
77 const unsigned int n_q_points = q.size();
78
79 // Store the quadrature
80 this->quad.initialize(q.get_points(), q.get_weights());
81
82 // Resize the weights
83 this->vertex_weights.resize(GeometryInfo<dim>::vertices_per_cell);
84
85 // see if we need the (transformation) shape function values
86 // and/or gradients and resize the necessary arrays
87 if (this->update_each &
89 compute_manifold_quadrature_weights(q);
90
91 if (this->update_each & update_covariant_transformation)
92 covariant.resize(n_q_points);
93
94 if (this->update_each & update_contravariant_transformation)
95 contravariant.resize(n_q_points);
96
97 if (this->update_each & update_volume_elements)
98 volume_elements.resize(n_q_points);
99}
100
101
102
103template <int dim, int spacedim>
104void
106 const UpdateFlags update_flags,
107 const Quadrature<dim> &q,
108 const unsigned int n_original_q_points)
109{
110 reinit(update_flags, q);
111
112 // Set to the size of a single quadrature object for faces, as the size set
113 // in in reinit() is for all points
114 if (this->update_each & update_covariant_transformation)
115 covariant.resize(n_original_q_points);
116
117 if (this->update_each & update_contravariant_transformation)
118 contravariant.resize(n_original_q_points);
119
120 if (this->update_each & update_volume_elements)
121 volume_elements.resize(n_original_q_points);
122
123 if (dim > 1)
124 {
125 if (this->update_each & update_boundary_forms)
126 {
127 aux.resize(dim - 1,
128 std::vector<Tensor<1, spacedim>>(n_original_q_points));
129
130 // Compute tangentials to the unit cell.
131 for (const unsigned int i : GeometryInfo<dim>::face_indices())
132 {
133 unit_tangentials[i].resize(n_original_q_points);
134 std::fill(unit_tangentials[i].begin(),
135 unit_tangentials[i].end(),
137 if (dim > 2)
138 {
139 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
140 .resize(n_original_q_points);
141 std::fill(
142 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
143 .begin(),
144 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
145 .end(),
147 }
148 }
149 }
150 }
151}
152
153
154
155template <int dim, int spacedim>
156void
158 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
159{
161 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
162 vertices[i] = cell->vertex(i);
163 this->cell = cell;
164}
165
166
167
168template <int dim, int spacedim>
169void
172{
173 cell_manifold_quadrature_weights.resize(
174 quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
175 for (unsigned int q = 0; q < quad.size(); ++q)
176 {
177 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
178 {
179 cell_manifold_quadrature_weights[q][i] =
181 }
182 }
183}
184
185
186
187template <int dim, int spacedim>
191
192
193
194template <int dim, int spacedim>
195std::unique_ptr<Mapping<dim, spacedim>>
197{
198 return std::make_unique<MappingManifold<dim, spacedim>>(*this);
199}
200
201
202
203template <int dim, int spacedim>
212
213
214
215template <int dim, int spacedim>
219 const Point<dim> &p) const
220{
221 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
222 std::array<double, GeometryInfo<dim>::vertices_per_cell> weights;
223
224 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
225 {
226 vertices[v] = cell->vertex(v);
228 }
229 return cell->get_manifold().get_new_point(
230 make_array_view(vertices.begin(), vertices.end()),
231 make_array_view(weights.begin(), weights.end()));
232}
233
234
235
236// In the code below, GCC tries to instantiate MappingManifold<3,4> when
237// seeing which of the overloaded versions of
238// do_transform_real_to_unit_cell_internal() to call. This leads to bad
239// error messages and, generally, nothing very good. Avoid this by ensuring
240// that this class exists, but does not have an inner InternalData
241// type, thereby ruling out the codim-1 version of the function
242// below when doing overload resolution.
243template <>
245{};
246
247
248
249template <int dim, int spacedim>
252 const UpdateFlags in) const
253{
254 // add flags if the respective quantities are necessary to compute
255 // what we need. note that some flags appear in both the conditions
256 // and in subsequent set operations. this leads to some circular
257 // logic. the only way to treat this is to iterate. since there are
258 // 5 if-clauses in the loop, it will take at most 5 iterations to
259 // converge. do them:
260 UpdateFlags out = in;
261 for (unsigned int i = 0; i < 5; ++i)
262 {
263 // The following is a little incorrect:
264 // If not applied on a face,
265 // update_boundary_forms does not
266 // make sense. On the other hand,
267 // it is necessary on a
268 // face. Currently,
269 // update_boundary_forms is simply
270 // ignored for the interior of a
271 // cell.
274
279
280 if (out &
285
286 // The contravariant transformation used in the Piola
287 // transformation, which requires the determinant of the Jacobi
288 // matrix of the transformation. Because we have no way of
289 // knowing here whether the finite elements wants to use the
290 // contravariant of the Piola transforms, we add the JxW values
291 // to the list of flags to be updated for each cell.
293 out |= update_JxW_values;
294
295 if (out & update_normal_vectors)
296 out |= update_JxW_values;
297 }
298
299 // Now throw an exception if we stumble upon something that was not
300 // implemented yet
305
306 return out;
307}
308
309
310
311template <int dim, int spacedim>
312std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
314 const Quadrature<dim> &q) const
315{
316 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
317 std::make_unique<InternalData>();
318 data_ptr->reinit(this->requires_update_flags(update_flags), q);
319
320 return data_ptr;
321}
322
323
324
325template <int dim, int spacedim>
326std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
328 const UpdateFlags update_flags,
329 const hp::QCollection<dim - 1> &quadrature) const
330{
331 AssertDimension(quadrature.size(), 1);
332
333 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
334 std::make_unique<InternalData>();
335 auto &data = dynamic_cast<InternalData &>(*data_ptr);
336 data.initialize_face(this->requires_update_flags(update_flags),
338 ReferenceCells::get_hypercube<dim>(), quadrature[0]),
339 quadrature[0].size());
340
341 return data_ptr;
342}
343
344
345
346template <int dim, int spacedim>
347std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
349 const UpdateFlags update_flags,
350 const Quadrature<dim - 1> &quadrature) const
351{
352 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
353 std::make_unique<InternalData>();
354 auto &data = dynamic_cast<InternalData &>(*data_ptr);
355 data.initialize_face(this->requires_update_flags(update_flags),
357 ReferenceCells::get_hypercube<dim>(), quadrature),
358 quadrature.size());
359
360 return data_ptr;
361}
362
363
364
365namespace internal
366{
367 namespace MappingManifoldImplementation
368 {
369 namespace
370 {
377 template <int dim, int spacedim>
378 void
379 maybe_compute_q_points(
380 const typename QProjector<dim>::DataSetDescriptor data_set,
381 const typename ::MappingManifold<dim, spacedim>::InternalData
382 &data,
383 std::vector<Point<spacedim>> &quadrature_points)
384 {
385 const UpdateFlags update_flags = data.update_each;
386
387 AssertDimension(data.vertices.size(),
389
390 if (update_flags & update_quadrature_points)
391 {
392 for (unsigned int point = 0; point < quadrature_points.size();
393 ++point)
394 {
395 quadrature_points[point] = data.manifold->get_new_point(
396 make_array_view(data.vertices),
398 data.cell_manifold_quadrature_weights[point + data_set]));
399 }
400 }
401 }
402
403
404
411 template <int dim, int spacedim>
412 void
413 maybe_update_Jacobians(
414 const typename ::QProjector<dim>::DataSetDescriptor data_set,
415 const typename ::MappingManifold<dim, spacedim>::InternalData
416 &data)
417 {
418 const UpdateFlags update_flags = data.update_each;
419
420 if (update_flags & update_contravariant_transformation)
421 {
422 const unsigned int n_q_points = data.contravariant.size();
423
424 std::fill(data.contravariant.begin(),
425 data.contravariant.end(),
427
429 data.vertices.size());
430 for (unsigned int point = 0; point < n_q_points; ++point)
431 {
432 // Start by figuring out how to compute the direction in
433 // the reference space:
434 const Point<dim> &p = data.quad.point(point + data_set);
435
436 // And get its image on the manifold:
437 const Point<spacedim> P = data.manifold->get_new_point(
438 make_array_view(data.vertices),
440 data.cell_manifold_quadrature_weights[point + data_set]));
441
442 // To compute the Jacobian, we choose dim points aligned
443 // with the dim reference axes, which are still in the
444 // given cell, and ask for the tangent vector in these
445 // directions. Choosing the points is somewhat arbitrary,
446 // so we try to be smart and we pick points which are
447 // on the opposite quadrant w.r.t. the evaluation
448 // point.
449 for (unsigned int i = 0; i < dim; ++i)
450 {
452 const double pi = p[i];
453 Assert(pi >= 0 && pi <= 1.0,
455 "Was expecting a quadrature point "
456 "inside the unit reference element."));
457
458 // In the length L, we store also the direction sign,
459 // which is positive, if the coordinate is < .5,
460 const double L = pi > .5 ? -pi : 1 - pi;
461
462 const Point<dim> np(p + L * ei);
463
464 // Get the weights to compute the np point in real space
465 for (const unsigned int j :
467 data.vertex_weights[j] =
469
470 const Point<spacedim> NP = data.manifold->get_new_point(
471 make_array_view(data.vertices),
472 make_array_view(data.vertex_weights));
473
474 const Tensor<1, spacedim> T =
475 data.manifold->get_tangent_vector(P, NP);
476
477 for (unsigned int d = 0; d < spacedim; ++d)
478 data.contravariant[point][d][i] = T[d] / L;
479 }
480 }
481
482 if (update_flags & update_covariant_transformation)
483 {
484 const unsigned int n_q_points = data.contravariant.size();
485 for (unsigned int point = 0; point < n_q_points; ++point)
486 {
487 data.covariant[point] =
488 (data.contravariant[point]).covariant_form();
489 }
490 }
491
492 if (update_flags & update_volume_elements)
493 {
494 const unsigned int n_q_points = data.contravariant.size();
495 for (unsigned int point = 0; point < n_q_points; ++point)
496 data.volume_elements[point] =
497 data.contravariant[point].determinant();
498 }
499 }
500 }
501 } // namespace
502 } // namespace MappingManifoldImplementation
503} // namespace internal
504
505
506
507template <int dim, int spacedim>
512 const Quadrature<dim> &quadrature,
513 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
515 &output_data) const
516{
517 // ensure that the following static_cast is really correct:
518 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
520 const InternalData &data = static_cast<const InternalData &>(internal_data);
521
522 const unsigned int n_q_points = quadrature.size();
523
524 data.store_vertices(cell);
525 data.manifold = &(cell->get_manifold());
526
527 internal::MappingManifoldImplementation::maybe_compute_q_points<dim,
528 spacedim>(
530 data,
531 output_data.quadrature_points);
532
533 internal::MappingManifoldImplementation::maybe_update_Jacobians<dim,
534 spacedim>(
536
537 const UpdateFlags update_flags = data.update_each;
538 const std::vector<double> &weights = quadrature.get_weights();
539
540 // Multiply quadrature weights by absolute value of Jacobian determinants or
541 // the area element g=sqrt(DX^t DX) in case of codim > 0
542
543 if (update_flags & (update_normal_vectors | update_JxW_values))
544 {
545 AssertDimension(output_data.JxW_values.size(), n_q_points);
546
547 Assert(!(update_flags & update_normal_vectors) ||
548 (output_data.normal_vectors.size() == n_q_points),
549 ExcDimensionMismatch(output_data.normal_vectors.size(),
550 n_q_points));
551
552
553 for (unsigned int point = 0; point < n_q_points; ++point)
554 {
555 if (dim == spacedim)
556 {
557 const double det = data.contravariant[point].determinant();
558
559 // check for distorted cells.
560
561 // TODO: this allows for anisotropies of up to 1e6 in 3d and
562 // 1e12 in 2d. might want to find a finer
563 // (dimension-independent) criterion
564 Assert(det > 1e-12 * Utilities::fixed_power<dim>(
565 cell->diameter() / std::sqrt(double(dim))),
567 cell->center(), det, point)));
568
569 output_data.JxW_values[point] = weights[point] * det;
570 }
571 // if dim==spacedim, then there is no cell normal to
572 // compute. since this is for FEValues (and not FEFaceValues),
573 // there are also no face normals to compute
574 else // codim>0 case
575 {
576 Tensor<1, spacedim> DX_t[dim];
577 for (unsigned int i = 0; i < spacedim; ++i)
578 for (unsigned int j = 0; j < dim; ++j)
579 DX_t[j][i] = data.contravariant[point][i][j];
580
581 Tensor<2, dim> G; // First fundamental form
582 for (unsigned int i = 0; i < dim; ++i)
583 for (unsigned int j = 0; j < dim; ++j)
584 G[i][j] = DX_t[i] * DX_t[j];
585
586 output_data.JxW_values[point] =
587 std::sqrt(determinant(G)) * weights[point];
588
589 if (update_flags & update_normal_vectors)
590 {
591 Assert(spacedim == dim + 1,
593 "There is no (unique) cell normal for " +
595 "-dimensional cells in " +
596 Utilities::int_to_string(spacedim) +
597 "-dimensional space. This only works if the "
598 "space dimension is one greater than the "
599 "dimensionality of the mesh cells."));
600
601 if (dim == 1)
602 output_data.normal_vectors[point] =
603 cross_product_2d(-DX_t[0]);
604 else // dim == 2
605 output_data.normal_vectors[point] =
606 cross_product_3d(DX_t[0], DX_t[1]);
607
608 output_data.normal_vectors[point] /=
609 output_data.normal_vectors[point].norm();
610
611 if (cell->direction_flag() == false)
612 output_data.normal_vectors[point] *= -1.;
613 }
614 } // codim>0 case
615 }
616 }
617
618
619
620 // copy values from InternalData to vector given by reference
621 if (update_flags & update_jacobians)
622 {
623 AssertDimension(output_data.jacobians.size(), n_q_points);
624 for (unsigned int point = 0; point < n_q_points; ++point)
625 output_data.jacobians[point] = data.contravariant[point];
626 }
627
628 // copy values from InternalData to vector given by reference
629 if (update_flags & update_inverse_jacobians)
630 {
631 AssertDimension(output_data.inverse_jacobians.size(), n_q_points);
632 for (unsigned int point = 0; point < n_q_points; ++point)
633 output_data.inverse_jacobians[point] =
634 data.covariant[point].transpose();
635 }
636
638}
639
640
641
642namespace internal
643{
644 namespace MappingManifoldImplementation
645 {
646 namespace
647 {
657 template <int dim, int spacedim>
658 void
659 maybe_compute_face_data(
660 const ::MappingManifold<dim, spacedim> &mapping,
661 const typename ::Triangulation<dim, spacedim>::cell_iterator
662 &cell,
663 const unsigned int face_no,
664 const unsigned int subface_no,
665 const unsigned int n_q_points,
666 const std::vector<double> &weights,
667 const typename ::MappingManifold<dim, spacedim>::InternalData
668 &data,
670 &output_data)
671 {
672 const UpdateFlags update_flags = data.update_each;
673
674 if (update_flags & update_boundary_forms)
675 {
676 AssertDimension(output_data.boundary_forms.size(), n_q_points);
677 if (update_flags & update_normal_vectors)
678 AssertDimension(output_data.normal_vectors.size(), n_q_points);
679 if (update_flags & update_JxW_values)
680 AssertDimension(output_data.JxW_values.size(), n_q_points);
681
682 // map the unit tangentials to the real cell. checking for d!=dim-1
683 // eliminates compiler warnings regarding unsigned int expressions <
684 // 0.
685 for (unsigned int d = 0; d != dim - 1; ++d)
686 {
688 data.unit_tangentials.size(),
690 Assert(
691 data.aux[d].size() <=
692 data
693 .unit_tangentials[face_no +
695 .size(),
697
698 mapping.transform(
700 data
701 .unit_tangentials[face_no +
704 data,
705 make_array_view(data.aux[d]));
706 }
707
708 // if dim==spacedim, we can use the unit tangentials to compute the
709 // boundary form by simply taking the cross product
710 if (dim == spacedim)
711 {
712 for (unsigned int i = 0; i < n_q_points; ++i)
713 switch (dim)
714 {
715 case 1:
716 // in 1d, we don't have access to any of the data.aux
717 // fields (because it has only dim-1 components), but we
718 // can still compute the boundary form by simply
719 // looking at the number of the face
720 output_data.boundary_forms[i][0] =
721 (face_no == 0 ? -1 : +1);
722 break;
723 case 2:
724 output_data.boundary_forms[i] =
725 cross_product_2d(data.aux[0][i]);
726 break;
727 case 3:
728 output_data.boundary_forms[i] =
729 cross_product_3d(data.aux[0][i], data.aux[1][i]);
730 break;
731 default:
733 }
734 }
735 else //(dim < spacedim)
736 {
737 // in the codim-one case, the boundary form results from the
738 // cross product of all the face tangential vectors and the cell
739 // normal vector
740 //
741 // to compute the cell normal, use the same method used in
742 // fill_fe_values for cells above
743 AssertDimension(data.contravariant.size(), n_q_points);
744
745 for (unsigned int point = 0; point < n_q_points; ++point)
746 {
747 switch (dim)
748 {
749 case 1:
750 {
751 // J is a tangent vector
752 output_data.boundary_forms[point] =
753 data.contravariant[point].transpose()[0];
754 output_data.boundary_forms[point] /=
755 (face_no == 0 ? -1. : +1.) *
756 output_data.boundary_forms[point].norm();
757
758 break;
759 }
760
761 case 2:
762 {
764 data.contravariant[point].transpose();
765
766 Tensor<1, spacedim> cell_normal =
767 cross_product_3d(DX_t[0], DX_t[1]);
768 cell_normal /= cell_normal.norm();
769
770 // then compute the face normal from the face
771 // tangent and the cell normal:
772 output_data.boundary_forms[point] =
773 cross_product_3d(data.aux[0][point], cell_normal);
774
775 break;
776 }
777
778 default:
780 }
781 }
782 }
783
784 if (update_flags & (update_normal_vectors | update_JxW_values))
785 for (unsigned int i = 0; i < output_data.boundary_forms.size();
786 ++i)
787 {
788 if (update_flags & update_JxW_values)
789 {
790 output_data.JxW_values[i] =
791 output_data.boundary_forms[i].norm() * weights[i];
792
793 if (subface_no != numbers::invalid_unsigned_int)
794 {
795 const double area_ratio =
797 cell->subface_case(face_no), subface_no);
798 output_data.JxW_values[i] *= area_ratio;
799 }
800 }
801
802 if (update_flags & update_normal_vectors)
803 output_data.normal_vectors[i] =
804 Point<spacedim>(output_data.boundary_forms[i] /
805 output_data.boundary_forms[i].norm());
806 }
807
808 if (update_flags & update_jacobians)
809 for (unsigned int point = 0; point < n_q_points; ++point)
810 output_data.jacobians[point] = data.contravariant[point];
811
812 if (update_flags & update_inverse_jacobians)
813 for (unsigned int point = 0; point < n_q_points; ++point)
814 output_data.inverse_jacobians[point] =
815 data.covariant[point].transpose();
816 }
817 }
818
819
826 template <int dim, int spacedim>
827 void
829 const ::MappingManifold<dim, spacedim> &mapping,
830 const typename ::Triangulation<dim, spacedim>::cell_iterator
831 &cell,
832 const unsigned int face_no,
833 const unsigned int subface_no,
834 const typename QProjector<dim>::DataSetDescriptor data_set,
835 const Quadrature<dim - 1> &quadrature,
836 const typename ::MappingManifold<dim, spacedim>::InternalData
837 &data,
839 &output_data)
840 {
841 data.store_vertices(cell);
842
843 data.manifold = &cell->face(face_no)->get_manifold();
844
845 maybe_compute_q_points<dim, spacedim>(data_set,
846 data,
847 output_data.quadrature_points);
848 maybe_update_Jacobians<dim, spacedim>(data_set, data);
849
851 cell,
852 face_no,
853 subface_no,
854 quadrature.size(),
855 quadrature.get_weights(),
856 data,
857 output_data);
858 }
859
860 template <int dim, int spacedim, int rank>
861 void
863 const ArrayView<const Tensor<rank, dim>> &input,
864 const MappingKind mapping_kind,
865 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
866 const ArrayView<Tensor<rank, spacedim>> &output)
867 {
868 AssertDimension(input.size(), output.size());
869 Assert((dynamic_cast<const typename ::
870 MappingManifold<dim, spacedim>::InternalData *>(
871 &mapping_data) != nullptr),
873 const typename ::MappingManifold<dim, spacedim>::InternalData
874 &data =
875 static_cast<const typename ::MappingManifold<dim, spacedim>::
876 InternalData &>(mapping_data);
877
878 switch (mapping_kind)
879 {
881 {
882 Assert(
885 "update_contravariant_transformation"));
886
887 for (unsigned int i = 0; i < output.size(); ++i)
888 output[i] =
889 apply_transformation(data.contravariant[i], input[i]);
890
891 return;
892 }
893
894 case mapping_piola:
895 {
896 Assert(
899 "update_contravariant_transformation"));
900 Assert(
901 data.update_each & update_volume_elements,
903 "update_volume_elements"));
904 Assert(rank == 1, ExcMessage("Only for rank 1"));
905 if (rank != 1)
906 return;
907
908 for (unsigned int i = 0; i < output.size(); ++i)
909 {
910 output[i] =
911 apply_transformation(data.contravariant[i], input[i]);
912 output[i] /= data.volume_elements[i];
913 }
914 return;
915 }
916 // We still allow this operation as in the
917 // reference cell Derivatives are Tensor
918 // rather than DerivativeForm
920 {
921 Assert(
924 "update_covariant_transformation"));
925
926 for (unsigned int i = 0; i < output.size(); ++i)
927 output[i] = apply_transformation(data.covariant[i], input[i]);
928
929 return;
930 }
931
932 default:
934 }
935 }
936
937
938 template <int dim, int spacedim, int rank>
939 void
941 const ArrayView<const Tensor<rank, dim>> &input,
942 const MappingKind mapping_kind,
943 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
944 const ArrayView<Tensor<rank, spacedim>> &output)
945 {
946 AssertDimension(input.size(), output.size());
947 Assert((dynamic_cast<const typename ::
948 MappingManifold<dim, spacedim>::InternalData *>(
949 &mapping_data) != nullptr),
951 const typename ::MappingManifold<dim, spacedim>::InternalData
952 &data =
953 static_cast<const typename ::MappingManifold<dim, spacedim>::
954 InternalData &>(mapping_data);
955
956 switch (mapping_kind)
957 {
959 {
960 Assert(
963 "update_covariant_transformation"));
964 Assert(
967 "update_contravariant_transformation"));
968 Assert(rank == 2, ExcMessage("Only for rank 2"));
969
970 for (unsigned int i = 0; i < output.size(); ++i)
971 {
973 apply_transformation(data.contravariant[i],
974 transpose(input[i]));
975 output[i] =
976 apply_transformation(data.covariant[i], A.transpose());
977 }
978
979 return;
980 }
981
983 {
984 Assert(
987 "update_covariant_transformation"));
988 Assert(rank == 2, ExcMessage("Only for rank 2"));
989
990 for (unsigned int i = 0; i < output.size(); ++i)
991 {
993 apply_transformation(data.covariant[i],
994 transpose(input[i]));
995 output[i] =
996 apply_transformation(data.covariant[i], A.transpose());
997 }
998
999 return;
1000 }
1001
1003 {
1004 Assert(
1007 "update_covariant_transformation"));
1008 Assert(
1011 "update_contravariant_transformation"));
1012 Assert(
1013 data.update_each & update_volume_elements,
1015 "update_volume_elements"));
1016 Assert(rank == 2, ExcMessage("Only for rank 2"));
1017
1018 for (unsigned int i = 0; i < output.size(); ++i)
1019 {
1021 apply_transformation(data.covariant[i], input[i]);
1023 apply_transformation(data.contravariant[i],
1024 A.transpose());
1025
1026 output[i] = transpose(T);
1027 output[i] /= data.volume_elements[i];
1028 }
1029
1030 return;
1031 }
1032
1033 default:
1035 }
1036 }
1037
1038
1039
1040 template <int dim, int spacedim>
1041 void
1043 const ArrayView<const Tensor<3, dim>> &input,
1044 const MappingKind mapping_kind,
1045 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1046 const ArrayView<Tensor<3, spacedim>> &output)
1047 {
1048 AssertDimension(input.size(), output.size());
1049 Assert((dynamic_cast<const typename ::
1050 MappingManifold<dim, spacedim>::InternalData *>(
1051 &mapping_data) != nullptr),
1053 const typename ::MappingManifold<dim, spacedim>::InternalData
1054 &data =
1055 static_cast<const typename ::MappingManifold<dim, spacedim>::
1056 InternalData &>(mapping_data);
1057
1058 switch (mapping_kind)
1059 {
1061 {
1062 Assert(
1065 "update_covariant_transformation"));
1066 Assert(
1069 "update_contravariant_transformation"));
1070
1071 for (unsigned int q = 0; q < output.size(); ++q)
1072 for (unsigned int i = 0; i < spacedim; ++i)
1073 {
1074 double tmp1[dim][dim];
1075 for (unsigned int J = 0; J < dim; ++J)
1076 for (unsigned int K = 0; K < dim; ++K)
1077 {
1078 tmp1[J][K] =
1079 data.contravariant[q][i][0] * input[q][0][J][K];
1080 for (unsigned int I = 1; I < dim; ++I)
1081 tmp1[J][K] +=
1082 data.contravariant[q][i][I] * input[q][I][J][K];
1083 }
1084 for (unsigned int j = 0; j < spacedim; ++j)
1085 {
1086 double tmp2[dim];
1087 for (unsigned int K = 0; K < dim; ++K)
1088 {
1089 tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
1090 for (unsigned int J = 1; J < dim; ++J)
1091 tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
1092 }
1093 for (unsigned int k = 0; k < spacedim; ++k)
1094 {
1095 output[q][i][j][k] =
1096 data.covariant[q][k][0] * tmp2[0];
1097 for (unsigned int K = 1; K < dim; ++K)
1098 output[q][i][j][k] +=
1099 data.covariant[q][k][K] * tmp2[K];
1100 }
1101 }
1102 }
1103 return;
1104 }
1105
1107 {
1108 Assert(
1111 "update_covariant_transformation"));
1112
1113 for (unsigned int q = 0; q < output.size(); ++q)
1114 for (unsigned int i = 0; i < spacedim; ++i)
1115 {
1116 double tmp1[dim][dim];
1117 for (unsigned int J = 0; J < dim; ++J)
1118 for (unsigned int K = 0; K < dim; ++K)
1119 {
1120 tmp1[J][K] =
1121 data.covariant[q][i][0] * input[q][0][J][K];
1122 for (unsigned int I = 1; I < dim; ++I)
1123 tmp1[J][K] +=
1124 data.covariant[q][i][I] * input[q][I][J][K];
1125 }
1126 for (unsigned int j = 0; j < spacedim; ++j)
1127 {
1128 double tmp2[dim];
1129 for (unsigned int K = 0; K < dim; ++K)
1130 {
1131 tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
1132 for (unsigned int J = 1; J < dim; ++J)
1133 tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
1134 }
1135 for (unsigned int k = 0; k < spacedim; ++k)
1136 {
1137 output[q][i][j][k] =
1138 data.covariant[q][k][0] * tmp2[0];
1139 for (unsigned int K = 1; K < dim; ++K)
1140 output[q][i][j][k] +=
1141 data.covariant[q][k][K] * tmp2[K];
1142 }
1143 }
1144 }
1145
1146 return;
1147 }
1148
1150 {
1151 Assert(
1154 "update_covariant_transformation"));
1155 Assert(
1158 "update_contravariant_transformation"));
1159 Assert(
1160 data.update_each & update_volume_elements,
1162 "update_volume_elements"));
1163
1164 for (unsigned int q = 0; q < output.size(); ++q)
1165 for (unsigned int i = 0; i < spacedim; ++i)
1166 {
1167 double factor[dim];
1168 for (unsigned int I = 0; I < dim; ++I)
1169 factor[I] =
1170 data.contravariant[q][i][I] / data.volume_elements[q];
1171 double tmp1[dim][dim];
1172 for (unsigned int J = 0; J < dim; ++J)
1173 for (unsigned int K = 0; K < dim; ++K)
1174 {
1175 tmp1[J][K] = factor[0] * input[q][0][J][K];
1176 for (unsigned int I = 1; I < dim; ++I)
1177 tmp1[J][K] += factor[I] * input[q][I][J][K];
1178 }
1179 for (unsigned int j = 0; j < spacedim; ++j)
1180 {
1181 double tmp2[dim];
1182 for (unsigned int K = 0; K < dim; ++K)
1183 {
1184 tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
1185 for (unsigned int J = 1; J < dim; ++J)
1186 tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
1187 }
1188 for (unsigned int k = 0; k < spacedim; ++k)
1189 {
1190 output[q][i][j][k] =
1191 data.covariant[q][k][0] * tmp2[0];
1192 for (unsigned int K = 1; K < dim; ++K)
1193 output[q][i][j][k] +=
1194 data.covariant[q][k][K] * tmp2[K];
1195 }
1196 }
1197 }
1198
1199 return;
1200 }
1201
1202 default:
1204 }
1205 }
1206
1207
1208
1209 template <int dim, int spacedim, int rank>
1210 void
1213 const MappingKind mapping_kind,
1214 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1216 {
1217 AssertDimension(input.size(), output.size());
1218 Assert((dynamic_cast<const typename ::
1219 MappingManifold<dim, spacedim>::InternalData *>(
1220 &mapping_data) != nullptr),
1222 const typename ::MappingManifold<dim, spacedim>::InternalData
1223 &data =
1224 static_cast<const typename ::MappingManifold<dim, spacedim>::
1225 InternalData &>(mapping_data);
1226
1227 switch (mapping_kind)
1228 {
1229 case mapping_covariant:
1230 {
1231 Assert(
1234 "update_covariant_transformation"));
1235
1236 for (unsigned int i = 0; i < output.size(); ++i)
1237 output[i] = apply_transformation(data.covariant[i], input[i]);
1238
1239 return;
1240 }
1241 default:
1243 }
1244 }
1245 } // namespace
1246 } // namespace MappingManifoldImplementation
1247} // namespace internal
1248
1249
1250
1251template <int dim, int spacedim>
1252void
1255 const unsigned int face_no,
1256 const hp::QCollection<dim - 1> &quadrature,
1257 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1259 &output_data) const
1260{
1261 AssertDimension(quadrature.size(), 1);
1262
1263 // ensure that the following cast is really correct:
1264 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1266 const InternalData &data = static_cast<const InternalData &>(internal_data);
1267
1268 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1269 *this,
1270 cell,
1271 face_no,
1274 ReferenceCells::get_hypercube<dim>(),
1275 face_no,
1276 cell->combined_face_orientation(face_no),
1277 quadrature[0].size()),
1278 quadrature[0],
1279 data,
1280 output_data);
1281}
1282
1283
1284
1285template <int dim, int spacedim>
1286void
1289 const unsigned int face_no,
1290 const unsigned int subface_no,
1291 const Quadrature<dim - 1> &quadrature,
1292 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1294 &output_data) const
1295{
1296 // ensure that the following cast is really correct:
1297 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1299 const InternalData &data = static_cast<const InternalData &>(internal_data);
1300
1301 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1302 *this,
1303 cell,
1304 face_no,
1305 subface_no,
1307 ReferenceCells::get_hypercube<dim>(),
1308 face_no,
1309 subface_no,
1310 cell->combined_face_orientation(face_no),
1311 quadrature.size(),
1312 cell->subface_case(face_no)),
1313 quadrature,
1314 data,
1315 output_data);
1316}
1317
1318
1319
1320template <int dim, int spacedim>
1321void
1323 const ArrayView<const Tensor<1, dim>> &input,
1324 const MappingKind mapping_kind,
1325 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1326 const ArrayView<Tensor<1, spacedim>> &output) const
1327{
1328 internal::MappingManifoldImplementation::transform_fields(input,
1329 mapping_kind,
1330 mapping_data,
1331 output);
1332}
1333
1334
1335
1336template <int dim, int spacedim>
1337void
1339 const ArrayView<const DerivativeForm<1, dim, spacedim>> &input,
1340 const MappingKind mapping_kind,
1341 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1342 const ArrayView<Tensor<2, spacedim>> &output) const
1343{
1344 internal::MappingManifoldImplementation::transform_differential_forms(
1345 input, mapping_kind, mapping_data, output);
1346}
1347
1348
1349
1350template <int dim, int spacedim>
1351void
1353 const ArrayView<const Tensor<2, dim>> &input,
1354 const MappingKind mapping_kind,
1355 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1356 const ArrayView<Tensor<2, spacedim>> &output) const
1357{
1358 switch (mapping_kind)
1359 {
1361 internal::MappingManifoldImplementation::transform_fields(input,
1362 mapping_kind,
1363 mapping_data,
1364 output);
1365 return;
1366
1370 internal::MappingManifoldImplementation::transform_gradients(
1371 input, mapping_kind, mapping_data, output);
1372 return;
1373 default:
1375 }
1376}
1377
1378
1379
1380template <int dim, int spacedim>
1381void
1383 const ArrayView<const DerivativeForm<2, dim, spacedim>> &input,
1384 const MappingKind mapping_kind,
1385 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1386 const ArrayView<Tensor<3, spacedim>> &output) const
1387{
1388 AssertDimension(input.size(), output.size());
1389 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
1391 const InternalData &data = static_cast<const InternalData &>(mapping_data);
1392
1393 switch (mapping_kind)
1394 {
1396 {
1399 "update_covariant_transformation"));
1400
1401 for (unsigned int q = 0; q < output.size(); ++q)
1402 for (unsigned int i = 0; i < spacedim; ++i)
1403 for (unsigned int j = 0; j < spacedim; ++j)
1404 {
1405 double tmp[dim];
1406 for (unsigned int K = 0; K < dim; ++K)
1407 {
1408 tmp[K] = data.covariant[q][j][0] * input[q][i][0][K];
1409 for (unsigned int J = 1; J < dim; ++J)
1410 tmp[K] += data.covariant[q][j][J] * input[q][i][J][K];
1411 }
1412 for (unsigned int k = 0; k < spacedim; ++k)
1413 {
1414 output[q][i][j][k] = data.covariant[q][k][0] * tmp[0];
1415 for (unsigned int K = 1; K < dim; ++K)
1416 output[q][i][j][k] += data.covariant[q][k][K] * tmp[K];
1417 }
1418 }
1419 return;
1420 }
1421
1422 default:
1424 }
1425}
1426
1427
1428
1429template <int dim, int spacedim>
1430void
1432 const ArrayView<const Tensor<3, dim>> &input,
1433 const MappingKind mapping_kind,
1434 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1435 const ArrayView<Tensor<3, spacedim>> &output) const
1436{
1437 switch (mapping_kind)
1438 {
1442 internal::MappingManifoldImplementation::transform_hessians(
1443 input, mapping_kind, mapping_data, output);
1444 return;
1445 default:
1447 }
1448}
1449
1450//--------------------------- Explicit instantiations -----------------------
1451#include "mapping_manifold.inst"
1452
1453
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:949
DerivativeForm< 1, spacedim, dim, Number > transpose() const
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
std::vector< double > volume_elements
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< std::vector< Tensor< 1, spacedim > > > aux
virtual std::size_t memory_consumption() const override
std::vector< std::vector< double > > cell_manifold_quadrature_weights
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
std::vector< Point< spacedim > > vertices
ObserverPointer< const Manifold< dim, spacedim > > manifold
std::vector< double > vertex_weights
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Triangulation< dim, spacedim >::cell_iterator cell
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
MappingManifold()=default
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
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 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 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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) 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
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
static constexpr Point< dim, Number > unit_vector(const unsigned int i)
static DataSetDescriptor cell()
Definition qprojector.h:461
const Point< dim > & point(const unsigned int i) const
const std::vector< double > & get_weights() const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
numbers::NumberTraits< Number >::real_type norm() const
unsigned int size() const
Definition collection.h:315
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
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_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_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
MappingKind
Definition mapping.h:79
@ mapping_piola
Definition mapping.h:114
@ mapping_covariant_gradient
Definition mapping.h:100
@ mapping_covariant
Definition mapping.h:89
@ mapping_contravariant
Definition mapping.h:94
@ mapping_contravariant_hessian
Definition mapping.h:156
@ mapping_covariant_hessian
Definition mapping.h:150
@ mapping_contravariant_gradient
Definition mapping.h:106
@ mapping_piola_gradient
Definition mapping.h:120
@ mapping_piola_hessian
Definition mapping.h:162
#define DEAL_II_NOT_IMPLEMENTED()
std::vector< index_type > data
Definition mpi.cc:735
std::size_t size
Definition mpi.cc:734
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
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 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)
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 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_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_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)
static const unsigned int invalid_unsigned_int
Definition types.h:220
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
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()
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)