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
manifold.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2014 - 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#include <deal.II/base/table.h>
16#include <deal.II/base/tensor.h>
18
19#include <deal.II/fe/fe_q.h>
20
23#include <deal.II/grid/tria.h>
26
27#include <boost/container/small_vector.hpp>
28
29#include <algorithm>
30#include <cmath>
31#include <limits>
32#include <memory>
33#include <numeric>
34
35
37
38/* -------------------------- Manifold --------------------- */
39#ifndef DOXYGEN
40
41template <int dim, int spacedim>
44 const ArrayView<const Point<spacedim>> &,
45 const Point<spacedim> &) const
46{
48 return {};
49}
50
51
52
53template <int dim, int spacedim>
56 const Point<spacedim> &p2,
57 const double w) const
58{
59 const std::array<Point<spacedim>, 2> vertices{{p1, p2}};
60 return project_to_manifold(make_array_view(vertices), w * p2 + (1 - w) * p1);
61}
62
63
64
65template <int dim, int spacedim>
68 const ArrayView<const Point<spacedim>> &surrounding_points,
69 const ArrayView<const double> &weights) const
70{
71 const double tol = 1e-10;
72 const unsigned int n_points = surrounding_points.size();
73
74 Assert(n_points > 0, ExcMessage("There should be at least one point."));
75
76 Assert(n_points == weights.size(),
78 "There should be as many surrounding points as weights given."));
79
80 Assert(std::abs(std::accumulate(weights.begin(), weights.end(), 0.0) - 1.0) <
81 tol,
82 ExcMessage("The weights for the individual points should sum to 1!"));
83
84 // First sort points in the order of their weights. This is done to
85 // produce unique points even if get_intermediate_points is not
86 // associative (as for the SphericalManifold).
87 boost::container::small_vector<unsigned int, 100> permutation(n_points);
88 std::iota(permutation.begin(), permutation.end(), 0u);
89 std::sort(permutation.begin(),
90 permutation.end(),
91 [&weights](const std::size_t a, const std::size_t b) {
92 return weights[a] < weights[b];
93 });
94
95 // Now loop over points in the order of their associated weight
97 Point<spacedim> p = surrounding_points[permutation[0]];
99 double w = weights[permutation[0]];
100 for (unsigned int i = 1; i < n_points; ++i)
101 {
102 double weight = 0.0;
103 if (std::abs(weights[permutation[i]] + w) < tol)
104 weight = 0.0;
105 else
106 weight = w / (weights[permutation[i]] + w);
107
108 if (std::abs(weight) > 1e-14)
109 {
110 p = get_intermediate_point(p,
111 surrounding_points[permutation[i]],
112 1.0 - weight);
113 }
114 else
115 {
116 p = surrounding_points[permutation[i]];
117 }
118 w += weights[permutation[i]];
119 }
120
121 return p;
122}
123
124
125
126template <int dim, int spacedim>
127void
129 const ArrayView<const Point<spacedim>> &surrounding_points,
130 const Table<2, double> &weights,
131 ArrayView<Point<spacedim>> new_points) const
132{
133 AssertDimension(surrounding_points.size(), weights.size(1));
134
135 for (unsigned int row = 0; row < weights.size(0); ++row)
136 {
137 new_points[row] =
138 get_new_point(surrounding_points, make_array_view(weights, row));
139 }
140}
141
142
143
144template <>
147 const Point<2> &p) const
148{
149 const int spacedim = 2;
150
151 // get the tangent vector from the point 'p' in the direction of the further
152 // one of the two vertices that make up the face of this 2d cell
153 const Tensor<1, spacedim> tangent =
154 ((p - face->vertex(0)).norm_square() > (p - face->vertex(1)).norm_square() ?
155 -get_tangent_vector(p, face->vertex(0)) :
156 get_tangent_vector(p, face->vertex(1)));
157
158 // then rotate it by 90 degrees
159 const Tensor<1, spacedim> normal = cross_product_2d(tangent);
160 return normal / normal.norm();
161}
162
163
164
165template <>
168 const Point<3> &p) const
169{
170 const int spacedim = 3;
171
172 const std::array<Point<spacedim>, 4> vertices{
173 {face->vertex(0), face->vertex(1), face->vertex(2), face->vertex(3)}};
174 const std::array<double, 4> distances{{vertices[0].distance(p),
175 vertices[1].distance(p),
176 vertices[2].distance(p),
177 vertices[3].distance(p)}};
178 const double max_distance =
179 std::max({distances[0], distances[1], distances[2], distances[3]});
180
181 // We need to find two tangential vectors to the given point p, but we do
182 // not know how the point is oriented against the face. We guess the two
183 // directions by assuming a flat topology and take the two directions that
184 // indicate the angle closest to a perpendicular one (i.e., cos(theta) close
185 // to zero). We start with an invalid value but the loops should always find
186 // a value.
187 double abs_cos_angle = std::numeric_limits<double>::max();
188 unsigned int first_index = numbers::invalid_unsigned_int,
189 second_index = numbers::invalid_unsigned_int;
190 for (unsigned int i = 0; i < 3; ++i)
191 if (distances[i] > 1e-8 * max_distance)
192 for (unsigned int j = i + 1; j < 4; ++j)
193 if (distances[j] > 1e-8 * max_distance)
194 {
195 const double new_angle = (p - vertices[i]) * (p - vertices[j]) /
196 (distances[i] * distances[j]);
197 // multiply by factor 0.999 to bias the search in a way that
198 // avoids trouble with roundoff
199 if (std::abs(new_angle) < 0.999 * abs_cos_angle)
200 {
201 abs_cos_angle = std::abs(new_angle);
202 first_index = i;
203 second_index = j;
204 }
205 }
207 ExcMessage("The search for possible directions did not succeed."));
208
209 // Compute tangents and normal for selected vertices
212 Tensor<1, spacedim> normal;
213
214 bool done = false;
215 std::vector<bool> tested_vertices(vertices.size(), false);
216 tested_vertices[first_index] = true;
217 tested_vertices[second_index] = true;
218
219 do
220 {
221 // Compute tangents and normal for selected vertices
222 t1 = get_tangent_vector(p, vertices[first_index]);
223 t2 = get_tangent_vector(p, vertices[second_index]);
224 normal = cross_product_3d(t1, t2);
225
226 // if normal is zero, try some other combination of vertices
227 if (normal.norm_square() < 1e4 * std::numeric_limits<double>::epsilon() *
228 t1.norm_square() * t2.norm_square())
229 {
230 // See if we have tested all vertices already
231 auto first_false =
232 std::find(tested_vertices.begin(), tested_vertices.end(), false);
233 if (first_false == tested_vertices.end())
234 {
235 done = true;
236 }
237 else
238 {
239 *first_false = true;
240 second_index = first_false - tested_vertices.begin();
241 }
242 }
243 else
244 {
245 done = true;
246 }
247 }
248 while (!done);
249
250 Assert(
251 normal.norm_square() > 1e4 * std::numeric_limits<double>::epsilon() *
252 t1.norm_square() * t2.norm_square(),
254 "Manifold::normal_vector was unable to find a suitable combination "
255 "of vertices to compute a normal on this face. We chose the secants "
256 "that are as orthogonal as possible, but tangents appear to be "
257 "linearly dependent. Check for distorted faces in your triangulation."));
258
259 // Now figure out if we need to flip the direction, we do this by comparing
260 // to a reference normal that would be the correct result if all vertices
261 // would lie in a plane
262 const Tensor<1, spacedim> rt1 = vertices[3] - vertices[0];
263 const Tensor<1, spacedim> rt2 = vertices[2] - vertices[1];
264 const Tensor<1, spacedim> reference_normal = cross_product_3d(rt1, rt2);
265
266 if (reference_normal * normal < 0.0)
267 normal *= -1.0;
268
269 return normal / normal.norm();
270}
271
272
273
274template <int dim, int spacedim>
277 const typename Triangulation<dim, spacedim>::face_iterator & /*face*/,
278 const Point<spacedim> & /*p*/) const
279{
281 return Tensor<1, spacedim>();
282}
283
284
285
286template <>
287void
290 FaceVertexNormals &n) const
291{
292 n[0] = cross_product_2d(get_tangent_vector(face->vertex(0), face->vertex(1)));
293 n[1] =
294 -cross_product_2d(get_tangent_vector(face->vertex(1), face->vertex(0)));
295
296 for (unsigned int i = 0; i < 2; ++i)
297 {
298 Assert(n[i].norm() != 0,
299 ExcInternalError("Something went wrong. The "
300 "computed normals have "
301 "zero length."));
302 n[i] /= n[i].norm();
303 }
304}
305
306
307
308template <>
309void
312 FaceVertexNormals &n) const
313{
314 n[0] = cross_product_3d(get_tangent_vector(face->vertex(0), face->vertex(1)),
315 get_tangent_vector(face->vertex(0), face->vertex(2)));
316
317 n[1] = cross_product_3d(get_tangent_vector(face->vertex(1), face->vertex(3)),
318 get_tangent_vector(face->vertex(1), face->vertex(0)));
319
320 n[2] = cross_product_3d(get_tangent_vector(face->vertex(2), face->vertex(0)),
321 get_tangent_vector(face->vertex(2), face->vertex(3)));
322
323 n[3] = cross_product_3d(get_tangent_vector(face->vertex(3), face->vertex(2)),
324 get_tangent_vector(face->vertex(3), face->vertex(1)));
325
326 for (unsigned int i = 0; i < 4; ++i)
327 {
328 Assert(n[i].norm() != 0,
329 ExcInternalError("Something went wrong. The "
330 "computed normals have "
331 "zero length."));
332 n[i] /= n[i].norm();
333 }
334}
335
336
337
338template <int dim, int spacedim>
339void
342 FaceVertexNormals &n) const
343{
344 for (unsigned int v = 0; v < face->reference_cell().n_vertices(); ++v)
345 {
346 n[v] = normal_vector(face, face->vertex(v));
347 n[v] /= n[v].norm();
348 }
349}
350
351
352
353template <int dim, int spacedim>
356 const typename Triangulation<dim, spacedim>::line_iterator &line) const
357{
358 const auto points_weights = Manifolds::get_default_points_and_weights(line);
359 return get_new_point(make_array_view(points_weights.first),
360 make_array_view(points_weights.second));
361}
362
363
364
365template <int dim, int spacedim>
368 const typename Triangulation<dim, spacedim>::quad_iterator &quad) const
369{
370 const auto points_weights = Manifolds::get_default_points_and_weights(quad);
371 return get_new_point(make_array_view(points_weights.first),
372 make_array_view(points_weights.second));
373}
374
375
376
377template <int dim, int spacedim>
380 const typename Triangulation<dim, spacedim>::face_iterator &face) const
381{
382 Assert(dim > 1, ExcImpossibleInDim(dim));
383
384 switch (dim)
385 {
386 case 2:
387 return get_new_point_on_line(face);
388 case 3:
389 return get_new_point_on_quad(face);
390 }
391
392 return {};
393}
394
395
396
397template <int dim, int spacedim>
400 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
401{
402 switch (dim)
403 {
404 case 1:
405 return get_new_point_on_line(cell);
406 case 2:
407 return get_new_point_on_quad(cell);
408 case 3:
409 return get_new_point_on_hex(cell);
410 }
411
412 return {};
413}
414
415
416
417template <>
421{
422 Assert(false, ExcImpossibleInDim(1));
423 return {};
424}
425
426
427
428template <>
432{
433 Assert(false, ExcImpossibleInDim(1));
434 return {};
435}
436
437
438
439template <>
443{
444 Assert(false, ExcImpossibleInDim(1));
445 return {};
446}
447
448
449
450template <>
454{
455 Assert(false, ExcImpossibleInDim(1));
456 return {};
457}
458
459
460
461template <>
465{
466 Assert(false, ExcImpossibleInDim(1));
467 return {};
468}
469
470
471
472template <>
476{
477 Assert(false, ExcImpossibleInDim(1));
478 return {};
479}
480
481
482
483template <int dim, int spacedim>
486 const typename Triangulation<dim, spacedim>::hex_iterator & /*hex*/) const
487{
488 Assert(false, ExcImpossibleInDim(dim));
489 return {};
490}
491
492
493
494template <>
497 const Triangulation<3, 3>::hex_iterator &hex) const
498{
499 const auto points_weights =
501 return get_new_point(make_array_view(points_weights.first),
502 make_array_view(points_weights.second));
503}
504
505
506
507template <int dim, int spacedim>
510 const Point<spacedim> &x2) const
511{
512 const double epsilon = 1e-8;
513
514 const std::array<Point<spacedim>, 2> points{{x1, x2}};
515 const std::array<double, 2> weights{{epsilon, 1.0 - epsilon}};
516 const Point<spacedim> neighbor_point =
517 get_new_point(make_array_view(points), make_array_view(weights));
518 return (neighbor_point - x1) / epsilon;
519}
520#endif
521/* -------------------------- FlatManifold --------------------- */
522
523namespace internal
524{
525 namespace
526 {
527 Tensor<1, 3>
528 normalized_alternating_product(const Tensor<1, 3> (&)[1])
529 {
530 // we get here from FlatManifold<2,3>::normal_vector, but
531 // the implementation below is bogus for this case anyway
532 // (see the assert at the beginning of that function).
534 return {};
535 }
536
537
538
539 Tensor<1, 3>
540 normalized_alternating_product(const Tensor<1, 3> (&basis_vectors)[2])
541 {
542 Tensor<1, 3> tmp = cross_product_3d(basis_vectors[0], basis_vectors[1]);
543 return tmp / tmp.norm();
544 }
545
546 } // namespace
547} // namespace internal
548
549#ifndef DOXYGEN
550
551template <int dim, int spacedim>
553 const Tensor<1, spacedim> &periodicity,
554 const double tolerance)
555 : periodicity(periodicity)
556 , tolerance(tolerance)
557{}
558
559
560
561template <int dim, int spacedim>
562std::unique_ptr<Manifold<dim, spacedim>>
564{
565 return std::make_unique<FlatManifold<dim, spacedim>>(periodicity, tolerance);
566}
567
568
569
570template <int dim, int spacedim>
573 const ArrayView<const Point<spacedim>> &surrounding_points,
574 const ArrayView<const double> &weights) const
575{
576 Assert(std::abs(std::accumulate(weights.begin(), weights.end(), 0.0) - 1.0) <
577 1e-10,
578 ExcMessage("The weights for the new point should sum to 1!"));
579
581
582 // if there is no periodicity, use a shortcut
583 if (periodicity == Tensor<1, spacedim>())
584 {
585 for (unsigned int i = 0; i < surrounding_points.size(); ++i)
586 p += surrounding_points[i] * weights[i];
587 }
588 else
589 {
590 Tensor<1, spacedim> minP = periodicity;
591
592 for (unsigned int d = 0; d < spacedim; ++d)
593 if (periodicity[d] > 0)
594 for (unsigned int i = 0; i < surrounding_points.size(); ++i)
595 {
596 minP[d] = std::min(minP[d], surrounding_points[i][d]);
597 Assert((surrounding_points[i][d] <
598 periodicity[d] + tolerance * periodicity[d]) ||
599 (surrounding_points[i][d] >=
600 -tolerance * periodicity[d]),
601 ExcPeriodicBox(d, surrounding_points[i], periodicity[d]));
602 }
603
604 // compute the weighted average point, possibly taking into account
605 // periodicity
606 for (unsigned int i = 0; i < surrounding_points.size(); ++i)
607 {
609 for (unsigned int d = 0; d < spacedim; ++d)
610 if (periodicity[d] > 0)
611 dp[d] =
612 ((surrounding_points[i][d] - minP[d]) > periodicity[d] / 2.0 ?
613 -periodicity[d] :
614 0.0);
615
616 p += (surrounding_points[i] + dp) * weights[i];
617 }
618
619 // if necessary, also adjust the weighted point by the periodicity
620 for (unsigned int d = 0; d < spacedim; ++d)
621 if (periodicity[d] > 0)
622 if (p[d] < 0)
623 p[d] += periodicity[d];
624 }
625
626 return project_to_manifold(surrounding_points, p);
627}
628
629
630
631template <int dim, int spacedim>
632void
634 const ArrayView<const Point<spacedim>> &surrounding_points,
635 const Table<2, double> &weights,
636 ArrayView<Point<spacedim>> new_points) const
637{
638 AssertDimension(surrounding_points.size(), weights.size(1));
639 if (weights.size(0) == 0)
640 return;
641 AssertDimension(new_points.size(), weights.size(0));
642
643 const std::size_t n_points = surrounding_points.size();
644
645 // if there is no periodicity, use an optimized implementation with
646 // VectorizedArray, otherwise go to the get_new_point function for adjusting
647 // the domain
648 if (periodicity == Tensor<1, spacedim>())
649 {
650 for (unsigned int row = 0; row < weights.size(0); ++row)
651 Assert(std::abs(std::accumulate(&weights(row, 0),
652 &weights(row, 0) + n_points,
653 0.0) -
654 1.0) < 1e-10,
655 ExcMessage("The weights for each of the points should sum to "
656 "1!"));
657
658 constexpr std::size_t n_lanes =
660 using VectorizedArrayType = VectorizedArray<double, n_lanes>;
661 const std::size_t n_regular_cols = (n_points / n_lanes) * n_lanes;
662 for (unsigned int row = 0; row < weights.size(0); row += n_lanes)
663 {
664 std::array<unsigned int, n_lanes> offsets;
665 // ensure to not access out of bounds, possibly duplicating some
666 // entries
667 for (std::size_t i = 0; i < n_lanes; ++i)
668 offsets[i] =
669 std::min<unsigned int>((row + i) * n_points,
670 (weights.size(0) - 1) * n_points);
672 for (std::size_t col = 0; col < n_regular_cols; col += n_lanes)
673 {
674 std::array<VectorizedArrayType, n_lanes> vectorized_weights;
676 &weights(0, 0) + col,
677 offsets.data(),
678 vectorized_weights.data());
679 for (std::size_t i = 0; i < n_lanes; ++i)
680 point += vectorized_weights[i] * surrounding_points[col + i];
681 }
682 for (std::size_t col = n_regular_cols; col < n_points; ++col)
683 {
684 VectorizedArrayType vectorized_weights;
685 vectorized_weights.gather(&weights(0, 0) + col, offsets.data());
686 point += vectorized_weights * surrounding_points[col];
687 }
688 for (unsigned int r = row;
689 r < std::min<unsigned int>(weights.size(0), row + n_lanes);
690 ++r)
691 {
692 // unpack and project to manifold
693 for (unsigned int d = 0; d < spacedim; ++d)
694 new_points[r][d] = point[d][r - row];
695 new_points[r] =
696 project_to_manifold(surrounding_points, new_points[r]);
697 }
698 }
699 }
700 else
701 for (unsigned int row = 0; row < weights.size(0); ++row)
702 new_points[row] =
703 get_new_point(surrounding_points,
704 ArrayView<const double>(&weights(row, 0), n_points));
705}
706
707
708
709template <int dim, int spacedim>
712 const ArrayView<const Point<spacedim>> & /*vertices*/,
713 const Point<spacedim> &candidate) const
714{
715 return candidate;
716}
717
718
719
720template <int dim, int spacedim>
723{
724 return periodicity;
725}
726
727
728
729template <int dim, int spacedim>
732 const Point<spacedim> &x2) const
733{
734 Tensor<1, spacedim> direction = x2 - x1;
735
736 // see if we have to take into account periodicity. if so, we need
737 // to make sure that if a distance in one coordinate direction
738 // is larger than half of the box length, then go the other way
739 // around (i.e., via the periodic box)
740 for (unsigned int d = 0; d < spacedim; ++d)
741 if (periodicity[d] > tolerance)
742 {
743 if (direction[d] < -periodicity[d] / 2)
744 direction[d] += periodicity[d];
745 else if (direction[d] > periodicity[d] / 2)
746 direction[d] -= periodicity[d];
747 }
748
749 return direction;
750}
751
752
753
754template <>
755void
759{
760 Assert(false, ExcImpossibleInDim(1));
761}
762
763
764
765template <>
766void
770{
772}
773
774
775
776template <>
777void
781{
783}
784
785
786
787template <>
788void
791 Manifold<2, 2>::FaceVertexNormals &face_vertex_normals) const
792{
793 const Tensor<1, 2> tangent = face->vertex(1) - face->vertex(0);
794 // We're in 2d. Faces are edges:
795 for (const unsigned int vertex : ReferenceCells::Line.vertex_indices())
796 // compute normals from tangent
797 face_vertex_normals[vertex] = Tensor<1, 2>({tangent[1], -tangent[0]});
798}
799
800
801
802template <>
803void
805 const Triangulation<2, 3>::face_iterator & /*face*/,
806 Manifold<2, 3>::FaceVertexNormals & /*face_vertex_normals*/) const
807{
809}
810
811
812
813template <>
814void
817 Manifold<3, 3>::FaceVertexNormals &face_vertex_normals) const
818{
820
821 static const unsigned int neighboring_vertices[4][2] = {{1, 2},
822 {3, 0},
823 {0, 3},
824 {2, 1}};
825 for (unsigned int vertex = 0; vertex < vertices_per_face; ++vertex)
826 {
827 // first define the two tangent vectors at the vertex by using the
828 // two lines radiating away from this vertex
829 const Tensor<1, 3> tangents[2] = {
830 face->vertex(neighboring_vertices[vertex][0]) - face->vertex(vertex),
831 face->vertex(neighboring_vertices[vertex][1]) - face->vertex(vertex)};
832
833 // then compute the normal by taking the cross product. since the
834 // normal is not required to be normalized, no problem here
835 face_vertex_normals[vertex] = cross_product_3d(tangents[0], tangents[1]);
836 }
837}
838
839
840
841template <>
844 const Point<1> &) const
845{
847 return {};
848}
849
850
851
852template <>
855 const Point<2> &) const
856{
858 return {};
859}
860
861
862
863template <>
866 const Point<3> &) const
867{
869 return {};
870}
871
872
873
874template <>
878 const Point<2> &p) const
879{
880 // In 2d, a face is just a straight line and
881 // we can use the 'standard' implementation.
882 return Manifold<2, 2>::normal_vector(face, p);
883}
884
885
886
887template <int dim, int spacedim>
891 const Point<spacedim> &p) const
892{
893 // I don't think the implementation below will work when dim!=spacedim;
894 // in fact, I believe that we don't even have enough information here,
895 // because we would need to know not only about the tangent vectors
896 // of the face, but also of the cell, to compute the normal vector.
897 // Someone will have to think about this some more.
898 Assert(dim == spacedim, ExcNotImplemented());
899
900 // in order to find out what the normal vector is, we first need to
901 // find the reference coordinates of the point p on the given face,
902 // or at least the reference coordinates of the closest point on the
903 // face
904 //
905 // in other words, we need to find a point xi so that f(xi)=||F(xi)-p||^2->min
906 // where F(xi) is the mapping. this algorithm is implemented in
907 // MappingQ1<dim,spacedim>::transform_real_to_unit_cell but only for cells,
908 // while we need it for faces here. it's also implemented in somewhat
909 // more generality there using the machinery of the MappingQ1 class
910 // while we really only need it for a specific case here
911 //
912 // in any case, the iteration we use here is a Gauss-Newton's iteration with
913 // xi^{n+1} = xi^n - H(xi^n)^{-1} J(xi^n)
914 // where
915 // J(xi) = (grad F(xi))^T (F(xi)-p)
916 // and
917 // H(xi) = [grad F(xi)]^T [grad F(xi)]
918 // In all this,
919 // F(xi) = sum_v vertex[v] phi_v(xi)
920 // We get the shape functions phi_v from an object of type FE_Q<dim-1>(1)
921
922 // We start at the center of the cell. If the face is a line or
923 // square, then the center is at 0.5 or (0.5,0.5). If the face is
924 // a triangle, then we start at the point (1/3,1/3).
925 const unsigned int facedim = dim - 1;
926
928
929 const auto face_reference_cell = face->reference_cell();
930
931 if ((dim <= 2) ||
932 (face_reference_cell == ReferenceCells::get_hypercube<facedim>()))
933 {
934 for (unsigned int i = 0; i < facedim; ++i)
935 xi[i] = 1. / 2;
936 }
937 else
938 {
939 for (unsigned int i = 0; i < facedim; ++i)
940 xi[i] = 1. / 3;
941 }
943 const double eps = 1e-12;
944 Tensor<1, spacedim> grad_F[facedim];
945 unsigned int iteration = 0;
946 while (true)
947 {
949 for (const unsigned int v : face->vertex_indices())
950 F +=
951 face->vertex(v) * face_reference_cell.d_linear_shape_function(xi, v);
952
953 for (unsigned int i = 0; i < facedim; ++i)
954 {
955 grad_F[i] = 0;
956 for (const unsigned int v : face->vertex_indices())
957 grad_F[i] +=
958 face->vertex(v) *
959 face_reference_cell.d_linear_shape_function_gradient(xi, v)[i];
960 }
961
963 for (unsigned int i = 0; i < facedim; ++i)
964 for (unsigned int j = 0; j < spacedim; ++j)
965 J[i] += grad_F[i][j] * (F - p)[j];
966
968 for (unsigned int i = 0; i < facedim; ++i)
969 for (unsigned int j = 0; j < facedim; ++j)
970 for (unsigned int k = 0; k < spacedim; ++k)
971 H[i][j] += grad_F[i][k] * grad_F[j][k];
972
973 const Tensor<1, facedim> delta_xi = -invert(H) * J;
974 xi += delta_xi;
975 ++iteration;
977 AssertThrow(iteration < 10,
979 "The Newton iteration to find the reference point "
980 "did not converge in 10 iterations. Do you have a "
981 "deformed cell? (See the glossary for a definition "
982 "of what a deformed cell is. You may want to output "
983 "the vertices of your cell."));
984
985 // It turns out that the check in reference coordinates with an absolute
986 // tolerance can cause a convergence failure of the Newton method as
987 // seen in tests/manifold/flat_manifold_09.cc. To work around this, also
988 // use a convergence check in world coordinates. This check has to be
989 // relative to the size of the face of course. Here we decided to use
990 // diameter because it works for non-planar faces and is cheap to
991 // compute:
992 const double normalized_delta_world = (F - p).norm() / face->diameter();
993
994 if (delta_xi.norm() < eps || normalized_delta_world < eps)
995 break;
996 }
997
998 // so now we have the reference coordinates xi of the point p.
999 // we then have to compute the normal vector, which we can do
1000 // by taking the (normalize) alternating product of all the tangent
1001 // vectors given by grad_F
1002 return internal::normalized_alternating_product(grad_F);
1003}
1004#endif
1005
1006/* -------------------------- ChartManifold --------------------- */
1007template <int dim, int spacedim, int chartdim>
1012
1013
1015template <int dim, int spacedim, int chartdim>
1018 const Point<spacedim> &p1,
1019 const Point<spacedim> &p2,
1020 const double w) const
1021{
1022 const std::array<Point<spacedim>, 2> points{{p1, p2}};
1023 const std::array<double, 2> weights{{1. - w, w}};
1024 return get_new_point(make_array_view(points), make_array_view(weights));
1025}
1026
1027
1028
1029template <int dim, int spacedim, int chartdim>
1032 const ArrayView<const Point<spacedim>> &surrounding_points,
1033 const ArrayView<const double> &weights) const
1034{
1035 const std::size_t n_points = surrounding_points.size();
1036
1037 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1038
1039 for (unsigned int i = 0; i < n_points; ++i)
1040 chart_points[i] = pull_back(surrounding_points[i]);
1041
1042 const Point<chartdim> p_chart =
1043 sub_manifold.get_new_point(chart_points, weights);
1044
1045 return push_forward(p_chart);
1046}
1047
1048
1049
1050template <int dim, int spacedim, int chartdim>
1051void
1053 const ArrayView<const Point<spacedim>> &surrounding_points,
1054 const Table<2, double> &weights,
1055 ArrayView<Point<spacedim>> new_points) const
1056{
1057 Assert(weights.size(0) > 0, ExcEmptyObject());
1058 AssertDimension(surrounding_points.size(), weights.size(1));
1059
1060 const std::size_t n_points = surrounding_points.size();
1061
1062 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1063 for (std::size_t i = 0; i < n_points; ++i)
1064 chart_points[i] = pull_back(surrounding_points[i]);
1065
1066 boost::container::small_vector<Point<chartdim>, 200> new_points_on_chart(
1067 weights.size(0));
1068 sub_manifold.get_new_points(chart_points, weights, new_points_on_chart);
1069
1070 for (std::size_t row = 0; row < weights.size(0); ++row)
1071 new_points[row] = push_forward(new_points_on_chart[row]);
1073
1074
1075
1076template <int dim, int spacedim, int chartdim>
1079 const Point<chartdim> &) const
1080{
1081 // function must be implemented in a derived class to be usable,
1082 // as discussed in this function's documentation
1083 Assert(false, ExcPureFunctionCalled());
1085}
1086
1087
1088
1089template <int dim, int spacedim, int chartdim>
1092 const Point<spacedim> &x1,
1093 const Point<spacedim> &x2) const
1094{
1097
1098 // ensure that the chart is not singular by asserting that its
1099 // derivative has a positive determinant. we need to make this
1100 // comparison relative to the size of the derivative. since the
1101 // determinant is the product of chartdim factors, take the
1102 // chartdim-th root of it in comparing against the size of the
1103 // derivative
1104 Assert(std::pow(std::abs(F_prime.determinant()), 1. / chartdim) >=
1105 1e-12 * F_prime.norm(),
1106 ExcMessage(
1107 "The derivative of a chart function must not be singular."));
1108
1109 const Tensor<1, chartdim> delta =
1110 sub_manifold.get_tangent_vector(pull_back(x1), pull_back(x2));
1111
1112 Tensor<1, spacedim> result;
1113 for (unsigned int i = 0; i < spacedim; ++i)
1114 result[i] += F_prime[i] * delta;
1115
1116 return result;
1117}
1118
1119
1120
1121template <int dim, int spacedim, int chartdim>
1122const Tensor<1, chartdim> &
1124{
1125 return sub_manifold.get_periodicity();
1126}
1127
1128// explicit instantiations
1129#include "grid/manifold.inst"
1130
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
iterator begin() const
Definition array_view.h:707
iterator end() const
Definition array_view.h:716
std::size_t size() const
Definition array_view.h:689
virtual Point< spacedim > push_forward(const Point< chartdim > &chart_point) const =0
const FlatManifold< chartdim, chartdim > sub_manifold
Definition manifold.h:1093
ChartManifold(const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
Definition manifold.cc:1008
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const override
Definition manifold.cc:1052
virtual Point< chartdim > pull_back(const Point< spacedim > &space_point) const =0
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
Definition manifold.cc:1091
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const
Definition manifold.cc:1078
const Tensor< 1, chartdim > & get_periodicity() const
Definition manifold.cc:1123
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
Definition manifold.cc:1017
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
Definition manifold.cc:1031
Number determinant() const
numbers::NumberTraits< Number >::real_type norm() const
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Manifold< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim > > &points, const Point< spacedim > &candidate) const override
const Tensor< 1, spacedim > & get_periodicity() const
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
FlatManifold(const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const override
virtual Point< spacedim > get_new_point_on_hex(const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) 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 void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const
Point< spacedim > get_new_point_on_face(const typename Triangulation< dim, spacedim >::face_iterator &face) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face > FaceVertexNormals
Definition manifold.h:305
Point< spacedim > get_new_point_on_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const
Definition point.h:113
numbers::NumberTraits< Number >::real_type norm() const
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
#define DEAL_II_NOT_IMPLEMENTED()
constexpr unsigned int GeometryInfo< dim >::vertices_per_face
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcPureFunctionCalled()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1692
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1668
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1644
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face_iterator
Definition tria.h:1596
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
std::pair< std::array< Point< MeshIteratorType::AccessorType::space_dimension >, n_default_points_per_cell< MeshIteratorType >()>, std::array< double, n_default_points_per_cell< MeshIteratorType >()> > get_default_points_and_weights(const MeshIteratorType &iterator, const bool with_interpolation=false)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:193
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
constexpr const ReferenceCell & get_hypercube()
constexpr ReferenceCell Line
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static constexpr unsigned int vertices_per_face
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition tensor.h:2647
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition tensor.h:2672
void vectorized_load_and_transpose(const unsigned int n_entries, const Number *in, const unsigned int *offsets, VectorizedArray< Number, width > *out)