deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
mapping_manifold.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
23
25
26#include <deal.II/fe/fe.h>
27#include <deal.II/fe/fe_tools.h>
31
33#include <deal.II/grid/tria.h>
35
37
38#include <algorithm>
39#include <array>
40#include <cmath>
41#include <memory>
42#include <numeric>
43
44
46
47template <int dim, int spacedim>
48std::size_t
65
66
67
68template <int dim, int spacedim>
69void
71 const UpdateFlags update_flags,
72 const Quadrature<dim> &q)
73{
74 // store the flags in the internal data object so we can access them
75 // in fill_fe_*_values()
76 this->update_each = update_flags;
77
78 const unsigned int n_q_points = q.size();
79
80 // Store the quadrature
81 this->quad.initialize(q.get_points(), q.get_weights());
82
83 // see if we need the (transformation) shape function values
84 // and/or gradients and resize the necessary arrays
85 if (this->update_each &
89 if (this->update_each & update_covariant_transformation)
90 covariant.resize(n_q_points);
91
92 if (this->update_each & update_contravariant_transformation)
93 contravariant.resize(n_q_points);
95 if (this->update_each & update_volume_elements)
96 volume_elements.resize(n_q_points);
97}
98
99
100
101template <int dim, int spacedim>
102void
104 const UpdateFlags update_flags,
105 const Quadrature<dim> &q,
106 const unsigned int n_original_q_points)
107{
108 reinit(update_flags, q);
110 // Set to the size of a single quadrature object for faces, as the size set
111 // in in reinit() is for all points
112 if (this->update_each & update_covariant_transformation)
113 covariant.resize(n_original_q_points);
114
115 if (this->update_each & update_contravariant_transformation)
116 contravariant.resize(n_original_q_points);
117
118 if (this->update_each & update_volume_elements)
119 volume_elements.resize(n_original_q_points);
120
121 if (dim > 1)
122 {
123 if (this->update_each & update_boundary_forms)
124 {
125 aux.resize(dim - 1,
126 std::vector<Tensor<1, spacedim>>(n_original_q_points));
127
128 // Compute tangentials to the unit cell.
129 for (const unsigned int i : GeometryInfo<dim>::face_indices())
130 {
131 unit_tangentials[i].resize(n_original_q_points);
132 std::fill(unit_tangentials[i].begin(),
133 unit_tangentials[i].end(),
135 if (dim > 2)
136 {
137 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
138 .resize(n_original_q_points);
139 std::fill(
140 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
141 .begin(),
142 unit_tangentials[GeometryInfo<dim>::faces_per_cell + i]
143 .end(),
145 }
146 }
147 }
148 }
149}
150
151
152
153template <int dim, int spacedim>
154void
157{
158 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
159 vertices[i] = cell->vertex(i);
160 this->cell = cell;
161}
162
163
164
165template <int dim, int spacedim>
166void
169{
171 for (unsigned int q = 0; q < quad.size(); ++q)
172 {
173 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
174 {
177 }
178 }
179}
180
181
182
183template <int dim, int spacedim>
187
188
189
190template <int dim, int spacedim>
191std::unique_ptr<Mapping<dim, spacedim>>
193{
194 return std::make_unique<MappingManifold<dim, spacedim>>(*this);
195}
196
197
198
199template <int dim, int spacedim>
208
209
210
211template <int dim, int spacedim>
215 const Point<dim> &p) const
216{
217 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
218 std::array<double, GeometryInfo<dim>::vertices_per_cell> weights;
219
220 for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
221 {
222 vertices[v] = cell->vertex(v);
224 }
225 return cell->get_manifold().get_new_point(
226 make_array_view(vertices.begin(), vertices.end()),
227 make_array_view(weights.begin(), weights.end()));
228}
229
230
231
232// In the code below, GCC tries to instantiate MappingManifold<3,4> when
233// seeing which of the overloaded versions of
234// do_transform_real_to_unit_cell_internal() to call. This leads to bad
235// error messages and, generally, nothing very good. Avoid this by ensuring
236// that this class exists, but does not have an inner InternalData
237// type, thereby ruling out the codim-1 version of the function
238// below when doing overload resolution.
239template <>
241{};
242
243
244
245template <int dim, int spacedim>
248 const UpdateFlags in) const
249{
250 // add flags if the respective quantities are necessary to compute
251 // what we need. note that some flags appear in both the conditions
252 // and in subsequent set operations. this leads to some circular
253 // logic. the only way to treat this is to iterate. since there are
254 // 5 if-clauses in the loop, it will take at most 5 iterations to
255 // converge. do them:
256 UpdateFlags out = in;
257 for (unsigned int i = 0; i < 5; ++i)
258 {
259 // The following is a little incorrect:
260 // If not applied on a face,
261 // update_boundary_forms does not
262 // make sense. On the other hand,
263 // it is necessary on a
264 // face. Currently,
265 // update_boundary_forms is simply
266 // ignored for the interior of a
267 // cell.
270
275
276 if (out &
281
282 // The contravariant transformation used in the Piola
283 // transformation, which requires the determinant of the Jacobi
284 // matrix of the transformation. Because we have no way of
285 // knowing here whether the finite elements wants to use the
286 // contravariant of the Piola transforms, we add the JxW values
287 // to the list of flags to be updated for each cell.
289 out |= update_JxW_values;
290
291 if (out & update_normal_vectors)
292 out |= update_JxW_values;
293 }
294
295 // Now throw an exception if we stumble upon something that was not
296 // implemented yet
301
302 return out;
303}
304
305
306
307template <int dim, int spacedim>
308std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
309MappingManifold<dim, spacedim>::get_data(const UpdateFlags update_flags,
310 const Quadrature<dim> &q) const
311{
312 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
313 std::make_unique<InternalData>();
314 data_ptr->reinit(this->requires_update_flags(update_flags), q);
315
316 return data_ptr;
317}
318
319
321template <int dim, int spacedim>
322std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
324 const UpdateFlags update_flags,
325 const hp::QCollection<dim - 1> &quadrature) const
326{
327 AssertDimension(quadrature.size(), 1);
328
329 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
330 std::make_unique<InternalData>();
331 auto &data = dynamic_cast<InternalData &>(*data_ptr);
332 data.initialize_face(this->requires_update_flags(update_flags),
334 ReferenceCells::get_hypercube<dim>(), quadrature[0]),
335 quadrature[0].size());
336
337 return data_ptr;
338}
339
341
342template <int dim, int spacedim>
343std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
345 const UpdateFlags update_flags,
346 const Quadrature<dim - 1> &quadrature) const
347{
348 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
349 std::make_unique<InternalData>();
350 auto &data = dynamic_cast<InternalData &>(*data_ptr);
351 data.initialize_face(this->requires_update_flags(update_flags),
354 quadrature.size());
355
356 return data_ptr;
357}
358
359
360
361namespace internal
364 {
365 namespace
366 {
373 template <int dim, int spacedim>
374 void
375 maybe_compute_q_points(
376 const typename QProjector<dim>::DataSetDescriptor data_set,
377 const typename ::MappingManifold<dim, spacedim>::InternalData
378 &data,
379 std::vector<Point<spacedim>> &quadrature_points)
380 {
381 const UpdateFlags update_flags = data.update_each;
382
383 if (update_flags & update_quadrature_points)
384 {
385 for (unsigned int point = 0; point < quadrature_points.size();
386 ++point)
387 {
388 quadrature_points[point] = data.manifold->get_new_point(
389 make_array_view(data.vertices),
391 data.cell_manifold_quadrature_weights[point + data_set]));
392 }
393 }
394 }
395
396
397
404 template <int dim, int spacedim>
405 void
406 maybe_update_Jacobians(
407 const typename ::QProjector<dim>::DataSetDescriptor data_set,
408 const typename ::MappingManifold<dim, spacedim>::InternalData
409 &data)
410 {
411 const UpdateFlags update_flags = data.update_each;
412
413 if (update_flags & update_contravariant_transformation)
414 {
415 const unsigned int n_q_points = data.contravariant.size();
416
417 std::fill(data.contravariant.begin(),
418 data.contravariant.end(),
420
421 for (unsigned int point = 0; point < n_q_points; ++point)
422 {
423 // Start by figuring out how to compute the direction in
424 // the reference space:
425 const Point<dim> &p = data.quad.point(point + data_set);
426
427 // And get its image on the manifold:
428 const Point<spacedim> P = data.manifold->get_new_point(
429 make_array_view(data.vertices),
431 data.cell_manifold_quadrature_weights[point + data_set]));
432
433 // To compute the Jacobian, we choose dim points aligned
434 // with the dim reference axes, which are still in the
435 // given cell, and ask for the tangent vector in these
436 // directions. Choosing the points is somewhat arbitrary,
437 // so we try to be smart and we pick points which are
438 // on the opposite quadrant w.r.t. the evaluation
439 // point.
440 for (unsigned int i = 0; i < dim; ++i)
441 {
443 const double pi = p[i];
444 Assert(pi >= 0 && pi <= 1.0,
446 "Was expecting a quadrature point "
447 "inside the unit reference element."));
448
449 // In the length L, we store also the direction sign,
450 // which is positive, if the coordinate is < .5,
451 const double L = pi > .5 ? -pi : 1 - pi;
452
453 const Point<dim> np(p + L * ei);
454
455 // Get the weights to compute the np point in real space
456 for (const unsigned int j :
458 data.vertex_weights[j] =
460
461 const Point<spacedim> NP = data.manifold->get_new_point(
462 make_array_view(data.vertices),
463 make_array_view(data.vertex_weights));
464
465 const Tensor<1, spacedim> T =
466 data.manifold->get_tangent_vector(P, NP);
467
468 for (unsigned int d = 0; d < spacedim; ++d)
469 data.contravariant[point][d][i] = T[d] / L;
470 }
471 }
472
473 if (update_flags & update_covariant_transformation)
474 {
475 const unsigned int n_q_points = data.contravariant.size();
476 for (unsigned int point = 0; point < n_q_points; ++point)
477 {
478 data.covariant[point] =
479 (data.contravariant[point]).covariant_form();
480 }
481 }
482
483 if (update_flags & update_volume_elements)
484 {
485 const unsigned int n_q_points = data.contravariant.size();
486 for (unsigned int point = 0; point < n_q_points; ++point)
487 data.volume_elements[point] =
488 data.contravariant[point].determinant();
489 }
490 }
491 }
492 } // namespace
493 } // namespace MappingManifoldImplementation
494} // namespace internal
495
496
497
498template <int dim, int spacedim>
503 const Quadrature<dim> &quadrature,
504 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
506 &output_data) const
507{
508 // ensure that the following static_cast is really correct:
509 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
511 const InternalData &data = static_cast<const InternalData &>(internal_data);
512
513 const unsigned int n_q_points = quadrature.size();
514
515 data.store_vertices(cell);
516 data.manifold = &(cell->get_manifold());
517
518 internal::MappingManifoldImplementation::maybe_compute_q_points<dim,
519 spacedim>(
521 data,
522 output_data.quadrature_points);
523
524 internal::MappingManifoldImplementation::maybe_update_Jacobians<dim,
525 spacedim>(
527
528 const UpdateFlags update_flags = data.update_each;
529 const std::vector<double> &weights = quadrature.get_weights();
530
531 // Multiply quadrature weights by absolute value of Jacobian determinants or
532 // the area element g=sqrt(DX^t DX) in case of codim > 0
533
534 if (update_flags & (update_normal_vectors | update_JxW_values))
535 {
536 AssertDimension(output_data.JxW_values.size(), n_q_points);
537
538 Assert(!(update_flags & update_normal_vectors) ||
539 (output_data.normal_vectors.size() == n_q_points),
540 ExcDimensionMismatch(output_data.normal_vectors.size(),
541 n_q_points));
542
543
544 for (unsigned int point = 0; point < n_q_points; ++point)
545 {
546 if (dim == spacedim)
547 {
548 const double det = data.contravariant[point].determinant();
549
550 // check for distorted cells.
551
552 // TODO: this allows for anisotropies of up to 1e6 in 3d and
553 // 1e12 in 2d. might want to find a finer
554 // (dimension-independent) criterion
556 cell->diameter() / std::sqrt(double(dim))),
558 cell->center(), det, point)));
559
560 output_data.JxW_values[point] = weights[point] * det;
561 }
562 // if dim==spacedim, then there is no cell normal to
563 // compute. since this is for FEValues (and not FEFaceValues),
564 // there are also no face normals to compute
565 else // codim>0 case
566 {
567 Tensor<1, spacedim> DX_t[dim];
568 for (unsigned int i = 0; i < spacedim; ++i)
569 for (unsigned int j = 0; j < dim; ++j)
570 DX_t[j][i] = data.contravariant[point][i][j];
571
572 Tensor<2, dim> G; // First fundamental form
573 for (unsigned int i = 0; i < dim; ++i)
574 for (unsigned int j = 0; j < dim; ++j)
575 G[i][j] = DX_t[i] * DX_t[j];
576
577 output_data.JxW_values[point] =
578 std::sqrt(determinant(G)) * weights[point];
579
580 if (update_flags & update_normal_vectors)
581 {
582 Assert(spacedim == dim + 1,
584 "There is no (unique) cell normal for " +
586 "-dimensional cells in " +
587 Utilities::int_to_string(spacedim) +
588 "-dimensional space. This only works if the "
589 "space dimension is one greater than the "
590 "dimensionality of the mesh cells."));
591
592 if (dim == 1)
593 output_data.normal_vectors[point] =
594 cross_product_2d(-DX_t[0]);
595 else // dim == 2
596 output_data.normal_vectors[point] =
597 cross_product_3d(DX_t[0], DX_t[1]);
598
599 output_data.normal_vectors[point] /=
600 output_data.normal_vectors[point].norm();
601
602 if (cell->direction_flag() == false)
603 output_data.normal_vectors[point] *= -1.;
604 }
605 } // codim>0 case
606 }
607 }
608
609
610
611 // copy values from InternalData to vector given by reference
612 if (update_flags & update_jacobians)
613 {
614 AssertDimension(output_data.jacobians.size(), n_q_points);
615 for (unsigned int point = 0; point < n_q_points; ++point)
616 output_data.jacobians[point] = data.contravariant[point];
617 }
618
619 // copy values from InternalData to vector given by reference
620 if (update_flags & update_inverse_jacobians)
621 {
622 AssertDimension(output_data.inverse_jacobians.size(), n_q_points);
623 for (unsigned int point = 0; point < n_q_points; ++point)
624 output_data.inverse_jacobians[point] =
625 data.covariant[point].transpose();
626 }
627
629}
630
631
632
633namespace internal
634{
635 namespace MappingManifoldImplementation
636 {
637 namespace
638 {
648 template <int dim, int spacedim>
649 void
650 maybe_compute_face_data(
651 const ::MappingManifold<dim, spacedim> &mapping,
652 const typename ::Triangulation<dim, spacedim>::cell_iterator
653 &cell,
654 const unsigned int face_no,
655 const unsigned int subface_no,
656 const unsigned int n_q_points,
657 const std::vector<double> &weights,
658 const typename ::MappingManifold<dim, spacedim>::InternalData
659 &data,
661 &output_data)
662 {
663 const UpdateFlags update_flags = data.update_each;
664
665 if (update_flags & update_boundary_forms)
666 {
667 AssertDimension(output_data.boundary_forms.size(), n_q_points);
668 if (update_flags & update_normal_vectors)
669 AssertDimension(output_data.normal_vectors.size(), n_q_points);
670 if (update_flags & update_JxW_values)
671 AssertDimension(output_data.JxW_values.size(), n_q_points);
672
673 // map the unit tangentials to the real cell. checking for d!=dim-1
674 // eliminates compiler warnings regarding unsigned int expressions <
675 // 0.
676 for (unsigned int d = 0; d != dim - 1; ++d)
677 {
679 data.unit_tangentials.size(),
681 Assert(
682 data.aux[d].size() <=
683 data
684 .unit_tangentials[face_no +
686 .size(),
688
689 mapping.transform(
691 data
692 .unit_tangentials[face_no +
695 data,
696 make_array_view(data.aux[d]));
697 }
698
699 // if dim==spacedim, we can use the unit tangentials to compute the
700 // boundary form by simply taking the cross product
701 if (dim == spacedim)
702 {
703 for (unsigned int i = 0; i < n_q_points; ++i)
704 switch (dim)
705 {
706 case 1:
707 // in 1d, we don't have access to any of the data.aux
708 // fields (because it has only dim-1 components), but we
709 // can still compute the boundary form by simply
710 // looking at the number of the face
711 output_data.boundary_forms[i][0] =
712 (face_no == 0 ? -1 : +1);
713 break;
714 case 2:
715 output_data.boundary_forms[i] =
716 cross_product_2d(data.aux[0][i]);
717 break;
718 case 3:
719 output_data.boundary_forms[i] =
720 cross_product_3d(data.aux[0][i], data.aux[1][i]);
721 break;
722 default:
724 }
725 }
726 else //(dim < spacedim)
727 {
728 // in the codim-one case, the boundary form results from the
729 // cross product of all the face tangential vectors and the cell
730 // normal vector
731 //
732 // to compute the cell normal, use the same method used in
733 // fill_fe_values for cells above
734 AssertDimension(data.contravariant.size(), n_q_points);
735
736 for (unsigned int point = 0; point < n_q_points; ++point)
737 {
738 switch (dim)
739 {
740 case 1:
741 {
742 // J is a tangent vector
743 output_data.boundary_forms[point] =
744 data.contravariant[point].transpose()[0];
745 output_data.boundary_forms[point] /=
746 (face_no == 0 ? -1. : +1.) *
747 output_data.boundary_forms[point].norm();
748
749 break;
750 }
751
752 case 2:
753 {
754 const DerivativeForm<1, spacedim, dim> DX_t =
755 data.contravariant[point].transpose();
756
757 Tensor<1, spacedim> cell_normal =
758 cross_product_3d(DX_t[0], DX_t[1]);
759 cell_normal /= cell_normal.norm();
760
761 // then compute the face normal from the face
762 // tangent and the cell normal:
763 output_data.boundary_forms[point] =
764 cross_product_3d(data.aux[0][point], cell_normal);
765
766 break;
767 }
768
769 default:
771 }
772 }
773 }
774
775 if (update_flags & (update_normal_vectors | update_JxW_values))
776 for (unsigned int i = 0; i < output_data.boundary_forms.size();
777 ++i)
778 {
779 if (update_flags & update_JxW_values)
780 {
781 output_data.JxW_values[i] =
782 output_data.boundary_forms[i].norm() * weights[i];
783
784 if (subface_no != numbers::invalid_unsigned_int)
785 {
786 const double area_ratio =
788 cell->subface_case(face_no), subface_no);
789 output_data.JxW_values[i] *= area_ratio;
790 }
791 }
792
793 if (update_flags & update_normal_vectors)
794 output_data.normal_vectors[i] =
795 Point<spacedim>(output_data.boundary_forms[i] /
796 output_data.boundary_forms[i].norm());
797 }
798
799 if (update_flags & update_jacobians)
800 for (unsigned int point = 0; point < n_q_points; ++point)
801 output_data.jacobians[point] = data.contravariant[point];
802
803 if (update_flags & update_inverse_jacobians)
804 for (unsigned int point = 0; point < n_q_points; ++point)
805 output_data.inverse_jacobians[point] =
806 data.covariant[point].transpose();
807 }
808 }
809
810
817 template <int dim, int spacedim>
818 void
820 const ::MappingManifold<dim, spacedim> &mapping,
821 const typename ::Triangulation<dim, spacedim>::cell_iterator
822 &cell,
823 const unsigned int face_no,
824 const unsigned int subface_no,
825 const typename QProjector<dim>::DataSetDescriptor data_set,
826 const Quadrature<dim - 1> &quadrature,
827 const typename ::MappingManifold<dim, spacedim>::InternalData
828 &data,
829 internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
830 &output_data)
831 {
832 data.store_vertices(cell);
833
834 data.manifold = &cell->face(face_no)->get_manifold();
835
836 maybe_compute_q_points<dim, spacedim>(data_set,
837 data,
838 output_data.quadrature_points);
839 maybe_update_Jacobians<dim, spacedim>(data_set, data);
840
842 cell,
843 face_no,
844 subface_no,
845 quadrature.size(),
846 quadrature.get_weights(),
847 data,
848 output_data);
849 }
850
851 template <int dim, int spacedim, int rank>
852 void
854 const ArrayView<const Tensor<rank, dim>> &input,
855 const MappingKind mapping_kind,
856 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
857 const ArrayView<Tensor<rank, spacedim>> &output)
858 {
859 AssertDimension(input.size(), output.size());
860 Assert((dynamic_cast<const typename ::
861 MappingManifold<dim, spacedim>::InternalData *>(
862 &mapping_data) != nullptr),
864 const typename ::MappingManifold<dim, spacedim>::InternalData
865 &data =
866 static_cast<const typename ::MappingManifold<dim, spacedim>::
867 InternalData &>(mapping_data);
868
869 switch (mapping_kind)
870 {
872 {
873 Assert(
874 data.update_each & update_contravariant_transformation,
876 "update_contravariant_transformation"));
877
878 for (unsigned int i = 0; i < output.size(); ++i)
879 output[i] =
880 apply_transformation(data.contravariant[i], input[i]);
881
882 return;
883 }
884
885 case mapping_piola:
886 {
887 Assert(
888 data.update_each & update_contravariant_transformation,
890 "update_contravariant_transformation"));
891 Assert(
892 data.update_each & update_volume_elements,
894 "update_volume_elements"));
895 Assert(rank == 1, ExcMessage("Only for rank 1"));
896 if (rank != 1)
897 return;
898
899 for (unsigned int i = 0; i < output.size(); ++i)
900 {
901 output[i] =
902 apply_transformation(data.contravariant[i], input[i]);
903 output[i] /= data.volume_elements[i];
904 }
905 return;
906 }
907 // We still allow this operation as in the
908 // reference cell Derivatives are Tensor
909 // rather than DerivativeForm
911 {
912 Assert(
913 data.update_each & update_contravariant_transformation,
915 "update_covariant_transformation"));
916
917 for (unsigned int i = 0; i < output.size(); ++i)
918 output[i] = apply_transformation(data.covariant[i], input[i]);
919
920 return;
921 }
922
923 default:
925 }
926 }
927
928
929 template <int dim, int spacedim, int rank>
930 void
932 const ArrayView<const Tensor<rank, dim>> &input,
933 const MappingKind mapping_kind,
934 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
935 const ArrayView<Tensor<rank, spacedim>> &output)
936 {
937 AssertDimension(input.size(), output.size());
938 Assert((dynamic_cast<const typename ::
939 MappingManifold<dim, spacedim>::InternalData *>(
940 &mapping_data) != nullptr),
942 const typename ::MappingManifold<dim, spacedim>::InternalData
943 &data =
944 static_cast<const typename ::MappingManifold<dim, spacedim>::
945 InternalData &>(mapping_data);
946
947 switch (mapping_kind)
948 {
950 {
951 Assert(
952 data.update_each & update_covariant_transformation,
954 "update_covariant_transformation"));
955 Assert(
956 data.update_each & update_contravariant_transformation,
958 "update_contravariant_transformation"));
959 Assert(rank == 2, ExcMessage("Only for rank 2"));
960
961 for (unsigned int i = 0; i < output.size(); ++i)
962 {
963 DerivativeForm<1, spacedim, dim> A =
964 apply_transformation(data.contravariant[i],
965 transpose(input[i]));
966 output[i] =
967 apply_transformation(data.covariant[i], A.transpose());
968 }
969
970 return;
971 }
972
974 {
975 Assert(
976 data.update_each & update_covariant_transformation,
978 "update_covariant_transformation"));
979 Assert(rank == 2, ExcMessage("Only for rank 2"));
980
981 for (unsigned int i = 0; i < output.size(); ++i)
982 {
983 DerivativeForm<1, spacedim, dim> A =
984 apply_transformation(data.covariant[i],
985 transpose(input[i]));
986 output[i] =
987 apply_transformation(data.covariant[i], A.transpose());
988 }
989
990 return;
991 }
992
994 {
995 Assert(
996 data.update_each & update_covariant_transformation,
998 "update_covariant_transformation"));
999 Assert(
1000 data.update_each & update_contravariant_transformation,
1002 "update_contravariant_transformation"));
1003 Assert(
1004 data.update_each & update_volume_elements,
1006 "update_volume_elements"));
1007 Assert(rank == 2, ExcMessage("Only for rank 2"));
1008
1009 for (unsigned int i = 0; i < output.size(); ++i)
1010 {
1011 DerivativeForm<1, spacedim, dim> A =
1012 apply_transformation(data.covariant[i], input[i]);
1013 Tensor<2, spacedim> T =
1014 apply_transformation(data.contravariant[i],
1015 A.transpose());
1016
1017 output[i] = transpose(T);
1018 output[i] /= data.volume_elements[i];
1019 }
1020
1021 return;
1022 }
1023
1024 default:
1026 }
1027 }
1028
1029
1030
1031 template <int dim, int spacedim>
1032 void
1034 const ArrayView<const Tensor<3, dim>> &input,
1035 const MappingKind mapping_kind,
1036 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1037 const ArrayView<Tensor<3, spacedim>> &output)
1038 {
1039 AssertDimension(input.size(), output.size());
1040 Assert((dynamic_cast<const typename ::
1041 MappingManifold<dim, spacedim>::InternalData *>(
1042 &mapping_data) != nullptr),
1044 const typename ::MappingManifold<dim, spacedim>::InternalData
1045 &data =
1046 static_cast<const typename ::MappingManifold<dim, spacedim>::
1047 InternalData &>(mapping_data);
1048
1049 switch (mapping_kind)
1050 {
1052 {
1053 Assert(
1054 data.update_each & update_covariant_transformation,
1056 "update_covariant_transformation"));
1057 Assert(
1058 data.update_each & update_contravariant_transformation,
1060 "update_contravariant_transformation"));
1061
1062 for (unsigned int q = 0; q < output.size(); ++q)
1063 output[q] =
1064 internal::apply_contravariant_hessian(data.covariant[q],
1065 data.contravariant[q],
1066 input[q]);
1067
1068 return;
1069 }
1070
1072 {
1073 Assert(
1074 data.update_each & update_covariant_transformation,
1076 "update_covariant_transformation"));
1077
1078 for (unsigned int q = 0; q < output.size(); ++q)
1079 output[q] =
1080 internal::apply_covariant_hessian(data.covariant[q],
1081 input[q]);
1082
1083 return;
1084 }
1085
1087 {
1088 Assert(
1089 data.update_each & update_covariant_transformation,
1091 "update_covariant_transformation"));
1092 Assert(
1093 data.update_each & update_contravariant_transformation,
1095 "update_contravariant_transformation"));
1096 Assert(
1097 data.update_each & update_volume_elements,
1099 "update_volume_elements"));
1100
1101 for (unsigned int q = 0; q < output.size(); ++q)
1102 output[q] =
1103 internal::apply_piola_hessian(data.covariant[q],
1104 data.contravariant[q],
1105 data.volume_elements[q],
1106 input[q]);
1107
1108 return;
1109 }
1110
1111 default:
1113 }
1114 }
1115
1116
1117
1118 template <int dim, int spacedim, int rank>
1119 void
1121 const ArrayView<const DerivativeForm<rank, dim, spacedim>> &input,
1122 const MappingKind mapping_kind,
1123 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1124 const ArrayView<Tensor<rank + 1, spacedim>> &output)
1125 {
1126 AssertDimension(input.size(), output.size());
1127 Assert((dynamic_cast<const typename ::
1128 MappingManifold<dim, spacedim>::InternalData *>(
1129 &mapping_data) != nullptr),
1131 const typename ::MappingManifold<dim, spacedim>::InternalData
1132 &data =
1133 static_cast<const typename ::MappingManifold<dim, spacedim>::
1134 InternalData &>(mapping_data);
1135
1136 switch (mapping_kind)
1137 {
1138 case mapping_covariant:
1139 {
1140 Assert(
1141 data.update_each & update_contravariant_transformation,
1143 "update_covariant_transformation"));
1144
1145 for (unsigned int i = 0; i < output.size(); ++i)
1146 output[i] = apply_transformation(data.covariant[i], input[i]);
1147
1148 return;
1149 }
1150 default:
1152 }
1153 }
1154 } // namespace
1155 } // namespace MappingManifoldImplementation
1156} // namespace internal
1157
1158
1159
1160template <int dim, int spacedim>
1161void
1164 const unsigned int face_no,
1165 const hp::QCollection<dim - 1> &quadrature,
1166 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1168 &output_data) const
1169{
1170 AssertDimension(quadrature.size(), 1);
1171
1172 // ensure that the following cast is really correct:
1173 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1175 const InternalData &data = static_cast<const InternalData &>(internal_data);
1176
1177 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1178 *this,
1179 cell,
1180 face_no,
1184 face_no,
1185 cell->combined_face_orientation(face_no),
1186 quadrature[0].size()),
1187 quadrature[0],
1188 data,
1189 output_data);
1190}
1191
1192
1193
1194template <int dim, int spacedim>
1195void
1198 const unsigned int face_no,
1199 const unsigned int subface_no,
1200 const Quadrature<dim - 1> &quadrature,
1201 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
1203 &output_data) const
1204{
1205 // ensure that the following cast is really correct:
1206 Assert((dynamic_cast<const InternalData *>(&internal_data) != nullptr),
1208 const InternalData &data = static_cast<const InternalData &>(internal_data);
1209
1210 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1211 *this,
1212 cell,
1213 face_no,
1214 subface_no,
1217 face_no,
1218 subface_no,
1219 cell->combined_face_orientation(face_no),
1220 quadrature.size(),
1221 cell->subface_case(face_no)),
1222 quadrature,
1223 data,
1224 output_data);
1225}
1226
1227
1228
1229template <int dim, int spacedim>
1230void
1232 const ArrayView<const Tensor<1, dim>> &input,
1233 const MappingKind mapping_kind,
1234 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1235 const ArrayView<Tensor<1, spacedim>> &output) const
1236{
1237 internal::MappingManifoldImplementation::transform_fields(input,
1238 mapping_kind,
1239 mapping_data,
1240 output);
1241}
1242
1243
1244
1245template <int dim, int spacedim>
1246void
1248 const ArrayView<const DerivativeForm<1, dim, spacedim>> &input,
1249 const MappingKind mapping_kind,
1250 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1251 const ArrayView<Tensor<2, spacedim>> &output) const
1252{
1253 internal::MappingManifoldImplementation::transform_differential_forms(
1254 input, mapping_kind, mapping_data, output);
1255}
1256
1257
1258
1259template <int dim, int spacedim>
1260void
1262 const ArrayView<const Tensor<2, dim>> &input,
1263 const MappingKind mapping_kind,
1264 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1265 const ArrayView<Tensor<2, spacedim>> &output) const
1266{
1267 switch (mapping_kind)
1268 {
1270 internal::MappingManifoldImplementation::transform_fields(input,
1271 mapping_kind,
1272 mapping_data,
1273 output);
1274 return;
1275
1279 internal::MappingManifoldImplementation::transform_gradients(
1280 input, mapping_kind, mapping_data, output);
1281 return;
1282 default:
1284 }
1285}
1286
1287
1288
1289template <int dim, int spacedim>
1290void
1292 const ArrayView<const DerivativeForm<2, dim, spacedim>> &input,
1293 const MappingKind mapping_kind,
1294 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1295 const ArrayView<Tensor<3, spacedim>> &output) const
1296{
1297 AssertDimension(input.size(), output.size());
1298 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
1300 const InternalData &data = static_cast<const InternalData &>(mapping_data);
1301
1302 switch (mapping_kind)
1303 {
1305 {
1308 "update_covariant_transformation"));
1309
1310 for (unsigned int q = 0; q < output.size(); ++q)
1311 output[q] =
1313
1314 return;
1315 }
1316
1317 default:
1319 }
1320}
1321
1322
1323
1324template <int dim, int spacedim>
1325void
1327 const ArrayView<const Tensor<3, dim>> &input,
1328 const MappingKind mapping_kind,
1329 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1330 const ArrayView<Tensor<3, spacedim>> &output) const
1331{
1332 switch (mapping_kind)
1333 {
1337 internal::MappingManifoldImplementation::transform_hessians(
1338 input, mapping_kind, mapping_data, output);
1339 return;
1340 default:
1342 }
1343}
1344
1345//--------------------------- Explicit instantiations -----------------------
1346#include "fe/mapping_manifold.inst"
1347
1348
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
double diameter(const Mapping< dim, spacedim > &mapping) const
bool direction_flag() const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) 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::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > vertices
std::vector< std::array< double, GeometryInfo< dim >::vertices_per_cell > > cell_manifold_quadrature_weights
std::vector< std::vector< Tensor< 1, spacedim > > > aux
virtual std::size_t memory_consumption() const override
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
ObserverPointer< const Manifold< dim, spacedim > > manifold
std::array< double, GeometryInfo< dim >::vertices_per_cell > 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
UpdateFlags update_each
Definition mapping.h:704
virtual std::size_t memory_consumption() const
Definition point.h:113
static constexpr Point< dim, Number > unit_vector(const unsigned int i)
Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to...
Definition qprojector.h:334
static DataSetDescriptor subface(const ReferenceCell &reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static DataSetDescriptor cell()
Definition qprojector.h:516
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
static Quadrature< dim > project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
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
const Manifold< dim, spacedim > & get_manifold() const
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
Point< spacedim > & vertex(const unsigned int i) const
unsigned int size() const
Definition collection.h:316
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Tensor< 1, spacedim, typename ProductType< Number1, Number2 >::type > apply_transformation(const DerivativeForm< 1, dim, spacedim, Number1 > &grad_F, const Tensor< 1, dim, Number2 > &d_x)
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
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)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
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:81
@ mapping_piola
Definition mapping.h:116
@ mapping_covariant_gradient
Definition mapping.h:102
@ mapping_covariant
Definition mapping.h:91
@ mapping_contravariant
Definition mapping.h:96
@ mapping_contravariant_hessian
Definition mapping.h:158
@ mapping_covariant_hessian
Definition mapping.h:152
@ mapping_contravariant_gradient
Definition mapping.h:108
@ mapping_piola_gradient
Definition mapping.h:122
@ mapping_piola_hessian
Definition mapping.h:164
types::geometric_orientation combined_face_orientation(const unsigned int face) const
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
constexpr char T
constexpr char A
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:193
constexpr const ReferenceCell & get_hypercube()
constexpr T fixed_power(const T t)
Definition utilities.h:943
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:466
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)
Tensor< 3, spacedim, Number > apply_contravariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_piola_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Number &volume_element, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_gradient(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 2, dim, spacedim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const Tensor< 3, dim, Number > &input)
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
static constexpr ndarray< Tensor< 1, dim >, faces_per_cell, dim - 1 > unit_tangential_vectors
static constexpr unsigned int vertices_per_cell
static constexpr unsigned int faces_per_cell
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 > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition tensor.h:2647
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition tensor.h:2672