Reference documentation for deal.II version 9.2.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\}}\)
manifold_lib.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2013 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 #include <deal.II/base/table.h>
18 #include <deal.II/base/tensor.h>
19 
20 #include <deal.II/fe/mapping.h>
21 
23 #include <deal.II/grid/tria.h>
26 
27 #include <deal.II/lac/vector.h>
28 
29 #include <cmath>
30 
32 
33 
34 namespace internal
35 {
36  // The pull_back function fails regularly in the compute_chart_points
37  // method, and, instead of throwing an exception, returns a point outside
38  // the unit cell. The individual coordinates of that point are given by the
39  // value below.
40  static constexpr double invalid_pull_back_coordinate = 20.0;
41 
42  // Rotate a given unit vector u around the axis dir
43  // where the angle is given by the length of dir.
44  // This is the exponential map for a sphere.
47  {
48  const double theta = dir.norm();
49  if (theta < 1.e-10)
50  {
51  return u;
52  }
53  else
54  {
55  const Tensor<1, 3> dirUnit = dir / theta;
56  const Tensor<1, 3> tmp =
57  std::cos(theta) * u + std::sin(theta) * dirUnit;
58  return tmp / tmp.norm();
59  }
60  }
61 
62  // Returns the direction to go from v to u
63  // projected to the plane perpendicular to the unit vector v.
64  // This one is more stable when u and v are nearly equal.
67  {
68  Tensor<1, 3> ans = u - v;
69  ans -= (ans * v) * v;
70  return ans; // ans = (u-v) - ((u-v)*v)*v
71  }
72 
73  // helper function to compute a vector orthogonal to a given one.
74  // does nothing unless spacedim == 3.
75  template <int spacedim>
77  compute_normal(const Tensor<1, spacedim> & /*vector*/,
78  bool /*normalize*/ = false)
79  {
80  return {};
81  }
82 
83  Point<3>
84  compute_normal(const Tensor<1, 3> &vector, bool normalize = false)
85  {
86  Assert(vector.norm_square() != 0.,
87  ExcMessage("The direction parameter must not be zero!"));
88  Point<3> normal;
89  if (std::abs(vector[0]) >= std::abs(vector[1]) &&
90  std::abs(vector[0]) >= std::abs(vector[2]))
91  {
92  normal[1] = -1.;
93  normal[2] = -1.;
94  normal[0] = (vector[1] + vector[2]) / vector[0];
95  }
96  else if (std::abs(vector[1]) >= std::abs(vector[0]) &&
97  std::abs(vector[1]) >= std::abs(vector[2]))
98  {
99  normal[0] = -1.;
100  normal[2] = -1.;
101  normal[1] = (vector[0] + vector[2]) / vector[1];
102  }
103  else
104  {
105  normal[0] = -1.;
106  normal[1] = -1.;
107  normal[2] = (vector[0] + vector[1]) / vector[2];
108  }
109  if (normalize)
110  normal /= normal.norm();
111  return normal;
112  }
113 } // namespace internal
114 
115 
116 
117 // ============================================================
118 // PolarManifold
119 // ============================================================
120 
121 template <int dim, int spacedim>
123  : ChartManifold<dim, spacedim, spacedim>(
124  PolarManifold<dim, spacedim>::get_periodicity())
125  , center(center)
126 {}
127 
128 
129 
130 template <int dim, int spacedim>
131 std::unique_ptr<Manifold<dim, spacedim>>
133 {
134  return std_cxx14::make_unique<PolarManifold<dim, spacedim>>(center);
135 }
136 
137 
138 
139 template <int dim, int spacedim>
142 {
143  Tensor<1, spacedim> periodicity;
144  // In two dimensions, theta is periodic.
145  // In three dimensions things are a little more complicated, since the only
146  // variable that is truly periodic is phi, while theta should be bounded
147  // between 0 and pi. There is currently no way to enforce this, so here we
148  // only fix periodicity for the last variable, corresponding to theta in 2d
149  // and phi in 3d.
150  periodicity[spacedim - 1] = 2 * numbers::PI;
151  return periodicity;
152 }
153 
154 
155 
156 template <int dim, int spacedim>
159  const Point<spacedim> &spherical_point) const
160 {
161  Assert(spherical_point[0] >= 0.0,
162  ExcMessage("Negative radius for given point."));
163  const double rho = spherical_point[0];
164  const double theta = spherical_point[1];
165 
166  Point<spacedim> p;
167  if (rho > 1e-10)
168  switch (spacedim)
169  {
170  case 2:
171  p[0] = rho * std::cos(theta);
172  p[1] = rho * std::sin(theta);
173  break;
174  case 3:
175  {
176  const double phi = spherical_point[2];
177  p[0] = rho * std::sin(theta) * std::cos(phi);
178  p[1] = rho * std::sin(theta) * std::sin(phi);
179  p[2] = rho * std::cos(theta);
180  break;
181  }
182  default:
183  Assert(false, ExcNotImplemented());
184  }
185  return p + center;
186 }
187 
188 
189 
190 template <int dim, int spacedim>
193  const Point<spacedim> &space_point) const
194 {
195  const Tensor<1, spacedim> R = space_point - center;
196  const double rho = R.norm();
197 
198  Point<spacedim> p;
199  p[0] = rho;
200 
201  switch (spacedim)
202  {
203  case 2:
204  {
205  p[1] = std::atan2(R[1], R[0]);
206  if (p[1] < 0)
207  p[1] += 2 * numbers::PI;
208  break;
209  }
210 
211  case 3:
212  {
213  const double z = R[2];
214  p[2] = std::atan2(R[1], R[0]); // phi
215  if (p[2] < 0)
216  p[2] += 2 * numbers::PI; // phi is periodic
217  p[1] = std::atan2(std::sqrt(R[0] * R[0] + R[1] * R[1]), z); // theta
218  break;
219  }
220 
221  default:
222  Assert(false, ExcNotImplemented());
223  }
224  return p;
225 }
226 
227 
228 
229 template <int dim, int spacedim>
232  const Point<spacedim> &spherical_point) const
233 {
234  Assert(spherical_point[0] >= 0.0,
235  ExcMessage("Negative radius for given point."));
236  const double rho = spherical_point[0];
237  const double theta = spherical_point[1];
238 
240  if (rho > 1e-10)
241  switch (spacedim)
242  {
243  case 2:
244  {
245  DX[0][0] = std::cos(theta);
246  DX[0][1] = -rho * std::sin(theta);
247  DX[1][0] = std::sin(theta);
248  DX[1][1] = rho * std::cos(theta);
249  break;
250  }
251 
252  case 3:
253  {
254  const double phi = spherical_point[2];
255  DX[0][0] = std::sin(theta) * std::cos(phi);
256  DX[0][1] = rho * std::cos(theta) * std::cos(phi);
257  DX[0][2] = -rho * std::sin(theta) * std::sin(phi);
258 
259  DX[1][0] = std::sin(theta) * std::sin(phi);
260  DX[1][1] = rho * std::cos(theta) * std::sin(phi);
261  DX[1][2] = rho * std::sin(theta) * std::cos(phi);
262 
263  DX[2][0] = std::cos(theta);
264  DX[2][1] = -rho * std::sin(theta);
265  DX[2][2] = 0;
266  break;
267  }
268 
269  default:
270  Assert(false, ExcNotImplemented());
271  }
272  return DX;
273 }
274 
275 
276 
277 namespace
278 {
279  template <int dim, int spacedim>
280  bool
281  spherical_face_is_horizontal(
282  const typename Triangulation<dim, spacedim>::face_iterator &face,
283  const Point<spacedim> & manifold_center)
284  {
285  // We test whether a face is horizontal by checking that the vertices
286  // all have roughly the same distance from the center: If the
287  // maximum deviation for the distances from the vertices to the
288  // center is less than 1.e-5 of the distance between vertices (as
289  // measured by the minimum distance from any of the other vertices
290  // to the first vertex), then we call this a horizontal face.
291  constexpr unsigned int n_vertices =
293  std::array<double, n_vertices> sqr_distances_to_center;
294  std::array<double, n_vertices - 1> sqr_distances_to_first_vertex;
295  sqr_distances_to_center[0] =
296  (face->vertex(0) - manifold_center).norm_square();
297  for (unsigned int i = 1; i < n_vertices; ++i)
298  {
299  sqr_distances_to_center[i] =
300  (face->vertex(i) - manifold_center).norm_square();
301  sqr_distances_to_first_vertex[i - 1] =
302  (face->vertex(i) - face->vertex(0)).norm_square();
303  }
304  const auto minmax_sqr_distance =
305  std::minmax_element(sqr_distances_to_center.begin(),
306  sqr_distances_to_center.end());
307  const auto min_sqr_distance_to_first_vertex =
308  std::min_element(sqr_distances_to_first_vertex.begin(),
309  sqr_distances_to_first_vertex.end());
310 
311  return (*minmax_sqr_distance.second - *minmax_sqr_distance.first <
312  1.e-10 * *min_sqr_distance_to_first_vertex);
313  }
314 } // namespace
315 
316 
317 
318 template <int dim, int spacedim>
321  const typename Triangulation<dim, spacedim>::face_iterator &face,
322  const Point<spacedim> & p) const
323 {
324  // Let us first test whether we are on a "horizontal" face
325  // (tangential to the sphere). In this case, the normal vector is
326  // easy to compute since it is proportional to the vector from the
327  // center to the point 'p'.
328  if (spherical_face_is_horizontal<dim, spacedim>(face, center))
329  {
330  // So, if this is a "horizontal" face, then just compute the normal
331  // vector as the one from the center to the point 'p', adequately
332  // scaled.
333  const Tensor<1, spacedim> unnormalized_spherical_normal = p - center;
334  const Tensor<1, spacedim> normalized_spherical_normal =
335  unnormalized_spherical_normal / unnormalized_spherical_normal.norm();
336  return normalized_spherical_normal;
337  }
338  else
339  // If it is not a horizontal face, just use the machinery of the
340  // base class.
342 
343  return Tensor<1, spacedim>();
344 }
345 
346 
347 
348 // ============================================================
349 // SphericalManifold
350 // ============================================================
351 
352 template <int dim, int spacedim>
354  const Point<spacedim> center)
355  : center(center)
356  , polar_manifold(center)
357 {}
358 
359 
360 
361 template <int dim, int spacedim>
362 std::unique_ptr<Manifold<dim, spacedim>>
364 {
365  return std_cxx14::make_unique<SphericalManifold<dim, spacedim>>(center);
366 }
367 
368 
369 
370 template <int dim, int spacedim>
373  const Point<spacedim> &p1,
374  const Point<spacedim> &p2,
375  const double w) const
376 {
377  const double tol = 1e-10;
378 
379  if ((p1 - p2).norm_square() < tol * tol || std::abs(w) < tol)
380  return p1;
381  else if (std::abs(w - 1.0) < tol)
382  return p2;
383 
384  // If the points are one dimensional then there is no need for anything but
385  // a linear combination.
386  if (spacedim == 1)
387  return Point<spacedim>(w * p2 + (1 - w) * p1);
388 
389  const Tensor<1, spacedim> v1 = p1 - center;
390  const Tensor<1, spacedim> v2 = p2 - center;
391  const double r1 = v1.norm();
392  const double r2 = v2.norm();
393 
394  Assert(r1 > tol && r2 > tol,
395  ExcMessage("p1 and p2 cannot coincide with the center."));
396 
397  const Tensor<1, spacedim> e1 = v1 / r1;
398  const Tensor<1, spacedim> e2 = v2 / r2;
399 
400  // Find the cosine of the angle gamma described by v1 and v2.
401  const double cosgamma = e1 * e2;
402 
403  // Points are collinear with the center (allow for 8*eps as a tolerance)
404  if (cosgamma < -1 + 8. * std::numeric_limits<double>::epsilon())
405  return center;
406 
407  // Points are along a line, in which case e1 and e2 are essentially the same.
408  if (cosgamma > 1 - 8. * std::numeric_limits<double>::epsilon())
409  return Point<spacedim>(center + w * v2 + (1 - w) * v1);
410 
411  // Find the angle sigma that corresponds to arclength equal to w. acos
412  // should never be undefined because we have ruled out the two special cases
413  // above.
414  const double sigma = w * std::acos(cosgamma);
415 
416  // Normal to v1 in the plane described by v1,v2,and the origin.
417  // Since p1 and p2 do not coincide n is not zero and well defined.
418  Tensor<1, spacedim> n = v2 - (v2 * e1) * e1;
419  const double n_norm = n.norm();
420  Assert(n_norm > 0,
421  ExcInternalError("n should be different from the null vector. "
422  "Probably, this means v1==v2 or v2==0."));
423 
424  n /= n_norm;
425 
426  // Find the point Q along O,v1 such that
427  // P1,V,P2 has measure sigma.
428  const Tensor<1, spacedim> P = std::cos(sigma) * e1 + std::sin(sigma) * n;
429 
430  // Project this point on the manifold.
431  return Point<spacedim>(center + (w * r2 + (1.0 - w) * r1) * P);
432 }
433 
434 
435 
436 template <int dim, int spacedim>
439  const Point<spacedim> &p1,
440  const Point<spacedim> &p2) const
441 {
442  const double tol = 1e-10;
443  (void)tol;
444 
445  Assert(p1 != p2, ExcMessage("p1 and p2 should not concide."));
446 
447  const Tensor<1, spacedim> v1 = p1 - center;
448  const Tensor<1, spacedim> v2 = p2 - center;
449  const double r1 = v1.norm();
450  const double r2 = v2.norm();
451 
452  Assert(r1 > tol, ExcMessage("p1 cannot coincide with the center."));
453 
454  Assert(r2 > tol, ExcMessage("p2 cannot coincide with the center."));
455 
456  const Tensor<1, spacedim> e1 = v1 / r1;
457  const Tensor<1, spacedim> e2 = v2 / r2;
458 
459  // Find the cosine of the angle gamma described by v1 and v2.
460  const double cosgamma = e1 * e2;
461 
462  Assert(cosgamma > -1 + 8. * std::numeric_limits<double>::epsilon(),
463  ExcMessage("p1 and p2 cannot lie on the same diameter and be opposite "
464  "respect to the center."));
465 
466  if (cosgamma > 1 - 8. * std::numeric_limits<double>::epsilon())
467  return v2 - v1;
468 
469  // Normal to v1 in the plane described by v1,v2,and the origin.
470  // Since p1 and p2 do not coincide n is not zero and well defined.
471  Tensor<1, spacedim> n = v2 - (v2 * e1) * e1;
472  const double n_norm = n.norm();
473  Assert(n_norm > 0,
474  ExcInternalError("n should be different from the null vector. "
475  "Probably, this means v1==v2 or v2==0."));
476 
477  n /= n_norm;
478 
479  // this is the derivative of the geodesic in get_intermediate_point
480  // derived with respect to w and inserting w=0.
481  const double gamma = std::acos(cosgamma);
482  return (r2 - r1) * e1 + r1 * gamma * n;
483 }
484 
485 
486 
487 template <int dim, int spacedim>
490  const typename Triangulation<dim, spacedim>::face_iterator &face,
491  const Point<spacedim> & p) const
492 {
493  // Let us first test whether we are on a "horizontal" face
494  // (tangential to the sphere). In this case, the normal vector is
495  // easy to compute since it is proportional to the vector from the
496  // center to the point 'p'.
497  if (spherical_face_is_horizontal<dim, spacedim>(face, center))
498  {
499  // So, if this is a "horizontal" face, then just compute the normal
500  // vector as the one from the center to the point 'p', adequately
501  // scaled.
502  const Tensor<1, spacedim> unnormalized_spherical_normal = p - center;
503  const Tensor<1, spacedim> normalized_spherical_normal =
504  unnormalized_spherical_normal / unnormalized_spherical_normal.norm();
505  return normalized_spherical_normal;
506  }
507  else
508  // If it is not a horizontal face, just use the machinery of the
509  // base class.
511 
512  return Tensor<1, spacedim>();
513 }
514 
515 
516 
517 template <>
518 void
522 {
523  Assert(false, ExcImpossibleInDim(1));
524 }
525 
526 
527 
528 template <>
529 void
533 {
534  Assert(false, ExcImpossibleInDim(1));
535 }
536 
537 
538 
539 template <int dim, int spacedim>
540 void
542  const typename Triangulation<dim, spacedim>::face_iterator &face,
543  typename Manifold<dim, spacedim>::FaceVertexNormals &face_vertex_normals)
544  const
545 {
546  // Let us first test whether we are on a "horizontal" face
547  // (tangential to the sphere). In this case, the normal vector is
548  // easy to compute since it is proportional to the vector from the
549  // center to the point 'p'.
550  if (spherical_face_is_horizontal<dim, spacedim>(face, center))
551  {
552  // So, if this is a "horizontal" face, then just compute the normal
553  // vector as the one from the center to the point 'p', adequately
554  // scaled.
555  for (unsigned int vertex = 0;
556  vertex < GeometryInfo<spacedim>::vertices_per_face;
557  ++vertex)
558  face_vertex_normals[vertex] = face->vertex(vertex) - center;
559  }
560  else
561  Manifold<dim, spacedim>::get_normals_at_vertices(face, face_vertex_normals);
562 }
563 
564 
565 
566 template <int dim, int spacedim>
567 void
569  const ArrayView<const Point<spacedim>> &surrounding_points,
570  const Table<2, double> & weights,
571  ArrayView<Point<spacedim>> new_points) const
572 {
573  AssertDimension(new_points.size(), weights.size(0));
574  AssertDimension(surrounding_points.size(), weights.size(1));
575 
576  get_new_points(surrounding_points, make_array_view(weights), new_points);
577 
578  return;
579 }
580 
581 
582 
583 template <int dim, int spacedim>
586  const ArrayView<const Point<spacedim>> &vertices,
587  const ArrayView<const double> & weights) const
588 {
589  // To avoid duplicating all of the logic in get_new_points, simply call it
590  // for one position.
591  Point<spacedim> new_point;
592  get_new_points(vertices,
593  weights,
594  make_array_view(&new_point, &new_point + 1));
595 
596  return new_point;
597 }
598 
599 
600 
601 template <int dim, int spacedim>
602 void
604  const ArrayView<const Point<spacedim>> &surrounding_points,
605  const ArrayView<const double> & weights,
606  ArrayView<Point<spacedim>> new_points) const
607 {
608  AssertDimension(weights.size(),
609  new_points.size() * surrounding_points.size());
610  const unsigned int weight_rows = new_points.size();
611  const unsigned int weight_columns = surrounding_points.size();
612 
613  if (surrounding_points.size() == 2)
614  {
615  for (unsigned int row = 0; row < weight_rows; ++row)
616  new_points[row] =
617  get_intermediate_point(surrounding_points[0],
618  surrounding_points[1],
619  weights[row * weight_columns + 1]);
620  return;
621  }
622 
623  boost::container::small_vector<std::pair<double, Tensor<1, spacedim>>, 100>
624  new_candidates(new_points.size());
625  boost::container::small_vector<Tensor<1, spacedim>, 100> directions(
626  surrounding_points.size(), Point<spacedim>());
627  boost::container::small_vector<double, 100> distances(
628  surrounding_points.size(), 0.0);
629  double max_distance = 0.;
630  for (unsigned int i = 0; i < surrounding_points.size(); ++i)
631  {
632  directions[i] = surrounding_points[i] - center;
633  distances[i] = directions[i].norm();
634 
635  if (distances[i] != 0.)
636  directions[i] /= distances[i];
637  else
638  Assert(false,
639  ExcMessage("One of the vertices coincides with the center. "
640  "This is not allowed!"));
641 
642  // Check if an estimate is good enough,
643  // this is often the case for sufficiently refined meshes.
644  for (unsigned int k = 0; k < i; ++k)
645  {
646  const double squared_distance =
647  (directions[i] - directions[k]).norm_square();
648  max_distance = std::max(max_distance, squared_distance);
649  }
650  }
651 
652  // Step 1: Check for some special cases, create simple linear guesses
653  // otherwise.
654  const double tolerance = 1e-10;
655  boost::container::small_vector<bool, 100> accurate_point_was_found(
656  new_points.size(), false);
657  const ArrayView<const Tensor<1, spacedim>> array_directions =
658  make_array_view(directions.begin(), directions.end());
659  const ArrayView<const double> array_distances =
660  make_array_view(distances.begin(), distances.end());
661  for (unsigned int row = 0; row < weight_rows; ++row)
662  {
663  new_candidates[row] =
664  guess_new_point(array_directions,
665  array_distances,
666  ArrayView<const double>(&weights[row * weight_columns],
667  weight_columns));
668 
669  // If the candidate is the center, mark it as found to avoid entering
670  // the Newton iteration in step 2, which would crash.
671  if (new_candidates[row].first == 0.0)
672  {
673  new_points[row] = center;
674  accurate_point_was_found[row] = true;
675  continue;
676  }
677 
678  // If not in 3D, just use the implementation from PolarManifold
679  // after we verified that the candidate is not the center.
680  if (spacedim < 3)
681  new_points[row] = polar_manifold.get_new_point(
682  surrounding_points,
683  ArrayView<const double>(&weights[row * weight_columns],
684  weight_columns));
685  }
686 
687  // In this case, we treated the case that the candidate is the center and
688  // obtained the new locations from the PolarManifold object otherwise.
689  if (spacedim < 3)
690  return;
691 
692  // If all the points are close to each other, we expect the estimate to
693  // be good enough. This tolerance was chosen such that the first iteration
694  // for a at least three time refined HyperShell mesh with radii .5 and 1.
695  // doesn't already succeed.
696  if (max_distance < 2e-2)
697  {
698  for (unsigned int row = 0; row < weight_rows; ++row)
699  new_points[row] =
700  center + new_candidates[row].first * new_candidates[row].second;
701 
702  return;
703  }
704 
705  // Step 2:
706  // Do more expensive Newton-style iterations to improve the estimate.
707 
708  // Search for duplicate directions and merge them to minimize the cost of
709  // the get_new_point function call below.
710  boost::container::small_vector<double, 1000> merged_weights(weights.size());
711  boost::container::small_vector<Tensor<1, spacedim>, 100> merged_directions(
712  surrounding_points.size(), Point<spacedim>());
713  boost::container::small_vector<double, 100> merged_distances(
714  surrounding_points.size(), 0.0);
715 
716  unsigned int n_unique_directions = 0;
717  for (unsigned int i = 0; i < surrounding_points.size(); ++i)
718  {
719  bool found_duplicate = false;
720 
721  // This inner loop is of @f$O(N^2)@f$ complexity, but
722  // surrounding_points.size() is usually at most 8 points large.
723  for (unsigned int j = 0; j < n_unique_directions; ++j)
724  {
725  const double squared_distance =
726  (directions[i] - directions[j]).norm_square();
727  if (!found_duplicate && squared_distance < 1e-28)
728  {
729  found_duplicate = true;
730  for (unsigned int row = 0; row < weight_rows; ++row)
731  merged_weights[row * weight_columns + j] +=
732  weights[row * weight_columns + i];
733  }
734  }
735 
736  if (found_duplicate == false)
737  {
738  merged_directions[n_unique_directions] = directions[i];
739  merged_distances[n_unique_directions] = distances[i];
740  for (unsigned int row = 0; row < weight_rows; ++row)
741  merged_weights[row * weight_columns + n_unique_directions] =
742  weights[row * weight_columns + i];
743 
744  ++n_unique_directions;
745  }
746  }
747 
748  // Search for duplicate weight rows and merge them to minimize the cost of
749  // the get_new_point function call below.
750  boost::container::small_vector<unsigned int, 100> merged_weights_index(
751  new_points.size(), numbers::invalid_unsigned_int);
752  for (unsigned int row = 0; row < weight_rows; ++row)
753  {
754  for (unsigned int existing_row = 0; existing_row < row; ++existing_row)
755  {
756  bool identical_weights = true;
757 
758  for (unsigned int weight_index = 0;
759  weight_index < n_unique_directions;
760  ++weight_index)
761  if (std::abs(merged_weights[row * weight_columns + weight_index] -
762  merged_weights[existing_row * weight_columns +
763  weight_index]) > tolerance)
764  {
765  identical_weights = false;
766  break;
767  }
768 
769  if (identical_weights)
770  {
771  merged_weights_index[row] = existing_row;
772  break;
773  }
774  }
775  }
776 
777  // Note that we only use the n_unique_directions first entries in the
778  // ArrayView
779  const ArrayView<const Tensor<1, spacedim>> array_merged_directions =
780  make_array_view(merged_directions.begin(),
781  merged_directions.begin() + n_unique_directions);
782  const ArrayView<const double> array_merged_distances =
783  make_array_view(merged_distances.begin(),
784  merged_distances.begin() + n_unique_directions);
785 
786  for (unsigned int row = 0; row < weight_rows; ++row)
787  if (!accurate_point_was_found[row])
788  {
789  if (merged_weights_index[row] == numbers::invalid_unsigned_int)
790  {
791  const ArrayView<const double> array_merged_weights(
792  &merged_weights[row * weight_columns], n_unique_directions);
793  new_candidates[row].second =
794  get_new_point(array_merged_directions,
795  array_merged_distances,
796  array_merged_weights,
797  Point<spacedim>(new_candidates[row].second));
798  }
799  else
800  new_candidates[row].second =
801  new_candidates[merged_weights_index[row]].second;
802 
803  new_points[row] =
804  center + new_candidates[row].first * new_candidates[row].second;
805  }
806 }
807 
808 
809 
810 template <int dim, int spacedim>
811 std::pair<double, Tensor<1, spacedim>>
813  const ArrayView<const Tensor<1, spacedim>> &directions,
814  const ArrayView<const double> & distances,
815  const ArrayView<const double> & weights) const
816 {
817  const double tolerance = 1e-10;
818  double rho = 0.;
819  Tensor<1, spacedim> candidate;
820 
821  // Perform a simple average ...
822  double total_weights = 0.;
823  for (unsigned int i = 0; i < directions.size(); i++)
824  {
825  // if one weight is one, return its direction
826  if (std::abs(1 - weights[i]) < tolerance)
827  return std::make_pair(distances[i], directions[i]);
828 
829  rho += distances[i] * weights[i];
830  candidate += directions[i] * weights[i];
831  total_weights += weights[i];
832  }
833 
834  // ... and normalize if the candidate is different from the origin.
835  const double norm = candidate.norm();
836  if (norm == 0.)
837  return std::make_pair(0.0, Point<spacedim>());
838  candidate /= norm;
839  rho /= total_weights;
840 
841  return std::make_pair(rho, candidate);
842 }
843 
844 
845 namespace
846 {
847  template <int spacedim>
849  do_get_new_point(const ArrayView<const Tensor<1, spacedim>> & /*directions*/,
850  const ArrayView<const double> & /*distances*/,
851  const ArrayView<const double> & /*weights*/,
852  const Point<spacedim> & /*candidate_point*/)
853  {
854  Assert(false, ExcNotImplemented());
855  return Point<spacedim>();
856  }
857 
858  template <>
859  Point<3>
860  do_get_new_point(const ArrayView<const Tensor<1, 3>> &directions,
861  const ArrayView<const double> & distances,
862  const ArrayView<const double> & weights,
863  const Point<3> & candidate_point)
864  {
865  (void)distances;
866 
867  AssertDimension(directions.size(), distances.size());
868  AssertDimension(directions.size(), weights.size());
869 
870  Point<3> candidate = candidate_point;
871  const unsigned int n_merged_points = directions.size();
872  const double tolerance = 1e-10;
873  const int max_iterations = 10;
874 
875  {
876  // If the candidate happens to coincide with a normalized
877  // direction, we return it. Otherwise, the Hessian would be singular.
878  for (unsigned int i = 0; i < n_merged_points; ++i)
879  {
880  const double squared_distance =
881  (candidate - directions[i]).norm_square();
882  if (squared_distance < tolerance * tolerance)
883  return candidate;
884  }
885 
886  // check if we only have two points now, in which case we can use the
887  // get_intermediate_point function
888  if (n_merged_points == 2)
889  {
890  SphericalManifold<3, 3> unit_manifold;
891  Assert(std::abs(weights[0] + weights[1] - 1.0) < 1e-13,
892  ExcMessage("Weights do not sum up to 1"));
893  Point<3> intermediate =
894  unit_manifold.get_intermediate_point(Point<3>(directions[0]),
895  Point<3>(directions[1]),
896  weights[1]);
897  return intermediate;
898  }
899 
900  Tensor<1, 3> vPerp;
901  Tensor<2, 2> Hessian;
902  Tensor<1, 2> gradient;
903  Tensor<1, 2> gradlocal;
904 
905  // On success we exit the loop early.
906  // Otherwise, we just take the result after max_iterations steps.
907  for (unsigned int i = 0; i < max_iterations; ++i)
908  {
909  // Step 2a: Find new descent direction
910 
911  // Get local basis for the estimate candidate
912  const Tensor<1, 3> Clocalx = internal::compute_normal(candidate);
913  const Tensor<1, 3> Clocaly = cross_product_3d(candidate, Clocalx);
914 
915  // For each vertices vector, compute the tangent vector from candidate
916  // towards the vertices vector -- its length is the spherical length
917  // from candidate to the vertices vector.
918  // Then compute its contribution to the Hessian.
919  gradient = 0.;
920  Hessian = 0.;
921  for (unsigned int i = 0; i < n_merged_points; ++i)
922  if (std::abs(weights[i]) > 1.e-15)
923  {
924  vPerp = internal::projected_direction(directions[i], candidate);
925  const double sinthetaSq = vPerp.norm_square();
926  const double sintheta = std::sqrt(sinthetaSq);
927  if (sintheta < tolerance)
928  {
929  Hessian[0][0] += weights[i];
930  Hessian[1][1] += weights[i];
931  }
932  else
933  {
934  const double costheta = (directions[i]) * candidate;
935  const double theta = std::atan2(sintheta, costheta);
936  const double sincthetaInv = theta / sintheta;
937 
938  const double cosphi = vPerp * Clocalx;
939  const double sinphi = vPerp * Clocaly;
940 
941  gradlocal[0] = cosphi;
942  gradlocal[1] = sinphi;
943  gradient += (weights[i] * sincthetaInv) * gradlocal;
944 
945  const double wt = weights[i] / sinthetaSq;
946  const double sinphiSq = sinphi * sinphi;
947  const double cosphiSq = cosphi * cosphi;
948  const double tt = sincthetaInv * costheta;
949  const double offdiag = cosphi * sinphi * wt * (1.0 - tt);
950  Hessian[0][0] += wt * (cosphiSq + tt * sinphiSq);
951  Hessian[0][1] += offdiag;
952  Hessian[1][0] += offdiag;
953  Hessian[1][1] += wt * (sinphiSq + tt * cosphiSq);
954  }
955  }
956 
957  Assert(determinant(Hessian) > tolerance, ExcInternalError());
958 
959  const Tensor<2, 2> inverse_Hessian = invert(Hessian);
960 
961  const Tensor<1, 2> xDisplocal = inverse_Hessian * gradient;
962  const Tensor<1, 3> xDisp =
963  xDisplocal[0] * Clocalx + xDisplocal[1] * Clocaly;
964 
965  // Step 2b: rotate candidate in direction xDisp for a new candidate.
966  const Point<3> candidateOld = candidate;
967  candidate =
968  Point<3>(internal::apply_exponential_map(candidate, xDisp));
969 
970  // Step 2c: return the new candidate if we didn't move
971  if ((candidate - candidateOld).norm_square() < tolerance * tolerance)
972  break;
973  }
974  }
975  return candidate;
976  }
977 } // namespace
978 
979 
980 
981 template <int dim, int spacedim>
984  const ArrayView<const Tensor<1, spacedim>> &,
985  const ArrayView<const double> &,
986  const ArrayView<const double> &,
987  const Point<spacedim> &) const
988 {
989  Assert(false, ExcNotImplemented());
990  return Point<spacedim>();
991 }
992 
993 
994 
995 template <>
996 Point<3>
998  const ArrayView<const Tensor<1, 3>> &directions,
999  const ArrayView<const double> & distances,
1000  const ArrayView<const double> & weights,
1001  const Point<3> & candidate_point) const
1002 {
1003  return do_get_new_point(directions, distances, weights, candidate_point);
1004 }
1005 
1006 
1007 
1008 template <>
1009 Point<3>
1011  const ArrayView<const Tensor<1, 3>> &directions,
1012  const ArrayView<const double> & distances,
1013  const ArrayView<const double> & weights,
1014  const Point<3> & candidate_point) const
1015 {
1016  return do_get_new_point(directions, distances, weights, candidate_point);
1017 }
1018 
1019 
1020 
1021 template <>
1022 Point<3>
1024  const ArrayView<const Tensor<1, 3>> &directions,
1025  const ArrayView<const double> & distances,
1026  const ArrayView<const double> & weights,
1027  const Point<3> & candidate_point) const
1028 {
1029  return do_get_new_point(directions, distances, weights, candidate_point);
1030 }
1031 
1032 
1033 
1034 // ============================================================
1035 // CylindricalManifold
1036 // ============================================================
1037 template <int dim, int spacedim>
1039  const double tolerance)
1040  : CylindricalManifold<dim, spacedim>(Point<spacedim>::unit_vector(axis),
1041  Point<spacedim>(),
1042  tolerance)
1043 {
1044  // do not use static_assert to make dimension-independent programming
1045  // easier.
1046  Assert(spacedim == 3,
1047  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1048 }
1049 
1050 
1051 
1052 template <int dim, int spacedim>
1054  const Tensor<1, spacedim> &direction,
1055  const Point<spacedim> & point_on_axis,
1056  const double tolerance)
1057  : ChartManifold<dim, spacedim, 3>(Tensor<1, 3>({0, 2. * numbers::PI, 0}))
1058  , normal_direction(internal::compute_normal(direction, true))
1059  , direction(direction / direction.norm())
1060  , point_on_axis(point_on_axis)
1061  , tolerance(tolerance)
1062 {
1063  // do not use static_assert to make dimension-independent programming
1064  // easier.
1065  Assert(spacedim == 3,
1066  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1067 }
1068 
1069 
1070 
1071 template <int dim, int spacedim>
1072 std::unique_ptr<Manifold<dim, spacedim>>
1074 {
1075  return std_cxx14::make_unique<CylindricalManifold<dim, spacedim>>(
1076  direction, point_on_axis, tolerance);
1077 }
1078 
1079 
1080 
1081 template <int dim, int spacedim>
1084  const ArrayView<const Point<spacedim>> &surrounding_points,
1085  const ArrayView<const double> & weights) const
1086 {
1087  Assert(spacedim == 3,
1088  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1089 
1090  // First check if the average in space lies on the axis.
1091  Point<spacedim> middle;
1092  double average_length = 0.;
1093  for (unsigned int i = 0; i < surrounding_points.size(); ++i)
1094  {
1095  middle += surrounding_points[i] * weights[i];
1096  average_length += surrounding_points[i].square() * weights[i];
1097  }
1098  middle -= point_on_axis;
1099  const double lambda = middle * direction;
1100 
1101  if ((middle - direction * lambda).square() < tolerance * average_length)
1102  return point_on_axis + direction * lambda;
1103  else // If not, using the ChartManifold should yield valid results.
1104  return ChartManifold<dim, spacedim, 3>::get_new_point(surrounding_points,
1105  weights);
1106 }
1107 
1108 
1109 
1110 template <int dim, int spacedim>
1111 Point<3>
1113  const Point<spacedim> &space_point) const
1114 {
1115  Assert(spacedim == 3,
1116  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1117 
1118  // First find the projection of the given point to the axis.
1119  const Tensor<1, spacedim> normalized_point = space_point - point_on_axis;
1120  const double lambda = normalized_point * direction;
1121  const Point<spacedim> projection = point_on_axis + direction * lambda;
1122  const Tensor<1, spacedim> p_diff = space_point - projection;
1123 
1124  // Then compute the angle between the projection direction and
1125  // another vector orthogonal to the direction vector.
1126  const double dot = normal_direction * p_diff;
1127  const double det = direction * cross_product_3d(normal_direction, p_diff);
1128  const double phi = std::atan2(det, dot);
1129 
1130  // Return distance from the axis, angle and signed distance on the axis.
1131  return Point<3>(p_diff.norm(), phi, lambda);
1132 }
1133 
1134 
1135 
1136 template <int dim, int spacedim>
1139  const Point<3> &chart_point) const
1140 {
1141  Assert(spacedim == 3,
1142  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1143 
1144  // Rotate the orthogonal direction by the given angle
1145  const double sine_r = std::sin(chart_point(1)) * chart_point(0);
1146  const double cosine_r = std::cos(chart_point(1)) * chart_point(0);
1147  const Tensor<1, spacedim> dxn = cross_product_3d(direction, normal_direction);
1148  const Tensor<1, spacedim> intermediate =
1149  normal_direction * cosine_r + dxn * sine_r;
1150 
1151  // Finally, put everything together.
1152  return point_on_axis + direction * chart_point(2) + intermediate;
1153 }
1154 
1155 
1156 
1157 template <int dim, int spacedim>
1160  const Point<3> &chart_point) const
1161 {
1162  Assert(spacedim == 3,
1163  ExcMessage("CylindricalManifold can only be used for spacedim==3!"));
1164 
1165  Tensor<2, 3> derivatives;
1166 
1167  // Rotate the orthogonal direction by the given angle
1168  const double sine = std::sin(chart_point(1));
1169  const double cosine = std::cos(chart_point(1));
1170  const Tensor<1, spacedim> dxn = cross_product_3d(direction, normal_direction);
1171  const Tensor<1, spacedim> intermediate =
1172  normal_direction * cosine + dxn * sine;
1173 
1174  // avoid compiler warnings
1175  constexpr int s0 = 0 % spacedim;
1176  constexpr int s1 = 1 % spacedim;
1177  constexpr int s2 = 2 % spacedim;
1178 
1179  // derivative w.r.t the radius
1180  derivatives[s0][s0] = intermediate[s0];
1181  derivatives[s1][s0] = intermediate[s1];
1182  derivatives[s2][s0] = intermediate[s2];
1183 
1184  // derivatives w.r.t the angle
1185  derivatives[s0][s1] = -normal_direction[s0] * sine + dxn[s0] * cosine;
1186  derivatives[s1][s1] = -normal_direction[s1] * sine + dxn[s1] * cosine;
1187  derivatives[s2][s1] = -normal_direction[s2] * sine + dxn[s2] * cosine;
1188 
1189  // derivatives w.r.t the direction of the axis
1190  derivatives[s0][s2] = direction[s0];
1191  derivatives[s1][s2] = direction[s1];
1192  derivatives[s2][s2] = direction[s2];
1193 
1194  return derivatives;
1195 }
1196 
1197 
1198 
1199 // ============================================================
1200 // EllipticalManifold
1201 // ============================================================
1202 template <int dim, int spacedim>
1204  const Point<spacedim> & center,
1205  const Tensor<1, spacedim> &major_axis_direction,
1206  const double eccentricity)
1207  : ChartManifold<dim, spacedim, spacedim>(
1208  EllipticalManifold<dim, spacedim>::get_periodicity())
1209  , direction(major_axis_direction)
1210  , center(center)
1211  , cosh_u(1.0 / eccentricity)
1212  , sinh_u(std::sqrt(cosh_u * cosh_u - 1.0))
1213 {
1214  // Throws an exception if dim!=2 || spacedim!=2.
1215  Assert(dim == 2 && spacedim == 2, ExcNotImplemented());
1216  // Throws an exception if eccentricity is not in range.
1217  Assert(std::signbit(cosh_u * cosh_u - 1.0) == false,
1218  ExcMessage(
1219  "Invalid eccentricity: It must satisfy 0 < eccentricity < 1."));
1220  const double direction_norm = direction.norm();
1221  Assert(direction_norm != 0.0,
1222  ExcMessage(
1223  "Invalid major axis direction vector: Null vector not allowed."));
1224  direction /= direction_norm;
1225 }
1226 
1227 
1228 
1229 template <int dim, int spacedim>
1230 std::unique_ptr<Manifold<dim, spacedim>>
1232 {
1233  const double eccentricity = 1.0 / cosh_u;
1234  return std_cxx14::make_unique<EllipticalManifold<dim, spacedim>>(
1235  center, direction, eccentricity);
1236 }
1237 
1238 
1239 
1240 template <int dim, int spacedim>
1243 {
1244  Tensor<1, spacedim> periodicity;
1245  // The second elliptical coordinate is periodic, while the first is not.
1246  // Enforce periodicity on the last variable.
1247  periodicity[spacedim - 1] = 2.0 * numbers::PI;
1248  return periodicity;
1249 }
1250 
1251 
1252 
1253 template <int dim, int spacedim>
1256 {
1257  Assert(false, ExcNotImplemented());
1258  return Point<spacedim>();
1259 }
1260 
1261 
1262 
1263 template <>
1264 Point<2>
1266 {
1267  const double cs = std::cos(chart_point[1]);
1268  const double sn = std::sin(chart_point[1]);
1269  // Coordinates in the reference frame (i.e. major axis direction is
1270  // x-axis)
1271  const double x = chart_point[0] * cosh_u * cs;
1272  const double y = chart_point[0] * sinh_u * sn;
1273  // Rotate them according to the major axis direction
1274  const Point<2> p(direction[0] * x - direction[1] * y,
1275  direction[1] * x + direction[0] * y);
1276  return p + center;
1277 }
1278 
1279 
1280 
1281 template <int dim, int spacedim>
1284 {
1285  Assert(false, ExcNotImplemented());
1286  return Point<spacedim>();
1287 }
1288 
1289 
1290 
1291 template <>
1292 Point<2>
1294 {
1295  // Moving space_point in the reference coordinate system.
1296  const double x0 = space_point[0] - center[0];
1297  const double y0 = space_point[1] - center[1];
1298  const double x = direction[0] * x0 + direction[1] * y0;
1299  const double y = -direction[1] * x0 + direction[0] * y0;
1300  const double pt0 =
1301  std::sqrt((x * x) / (cosh_u * cosh_u) + (y * y) / (sinh_u * sinh_u));
1302  // If the radius is exactly zero, the point coincides with the origin.
1303  if (pt0 == 0.0)
1304  {
1305  return center;
1306  }
1307  double cos_eta = x / (pt0 * cosh_u);
1308  if (cos_eta < -1.0)
1309  {
1310  cos_eta = -1.0;
1311  }
1312  if (cos_eta > 1.0)
1313  {
1314  cos_eta = 1.0;
1315  }
1316  const double eta = std::acos(cos_eta);
1317  const double pt1 = (std::signbit(y) ? 2.0 * numbers::PI - eta : eta);
1318  return {pt0, pt1};
1319 }
1320 
1321 
1322 
1323 template <int dim, int spacedim>
1326  const Point<spacedim> &) const
1327 {
1328  Assert(false, ExcNotImplemented());
1329  return {};
1330 }
1331 
1332 
1333 
1334 template <>
1337  const Point<2> &chart_point) const
1338 {
1339  const double cs = std::cos(chart_point[1]);
1340  const double sn = std::sin(chart_point[1]);
1341  Tensor<2, 2> dX;
1342  dX[0][0] = cosh_u * cs;
1343  dX[0][1] = -chart_point[0] * cosh_u * sn;
1344  dX[1][0] = sinh_u * sn;
1345  dX[1][1] = chart_point[0] * sinh_u * cs;
1346 
1347  // rotate according to the major axis direction
1349  {{+direction[0], -direction[1]}, {direction[1], direction[0]}}};
1350 
1351  return rot * dX;
1352 }
1353 
1354 
1355 
1356 // ============================================================
1357 // FunctionManifold
1358 // ============================================================
1359 template <int dim, int spacedim, int chartdim>
1361  const Function<chartdim> & push_forward_function,
1362  const Function<spacedim> & pull_back_function,
1363  const Tensor<1, chartdim> &periodicity,
1364  const double tolerance)
1365  : ChartManifold<dim, spacedim, chartdim>(periodicity)
1366  , const_map()
1367  , push_forward_function(&push_forward_function)
1368  , pull_back_function(&pull_back_function)
1369  , tolerance(tolerance)
1370  , owns_pointers(false)
1371  , finite_difference_step(0)
1372 {
1373  AssertDimension(push_forward_function.n_components, spacedim);
1374  AssertDimension(pull_back_function.n_components, chartdim);
1375 }
1376 
1377 
1378 
1379 template <int dim, int spacedim, int chartdim>
1381  std::unique_ptr<Function<chartdim>> push_forward,
1382  std::unique_ptr<Function<spacedim>> pull_back,
1383  const Tensor<1, chartdim> & periodicity,
1384  const double tolerance)
1385  : ChartManifold<dim, spacedim, chartdim>(periodicity)
1386  , const_map()
1387  , push_forward_function(push_forward.release())
1388  , pull_back_function(pull_back.release())
1389  , tolerance(tolerance)
1390  , owns_pointers(true)
1391  , finite_difference_step(0)
1392 {
1393  AssertDimension(push_forward_function->n_components, spacedim);
1394  AssertDimension(pull_back_function->n_components, chartdim);
1395 }
1396 
1397 
1398 
1399 template <int dim, int spacedim, int chartdim>
1401  const std::string push_forward_expression,
1402  const std::string pull_back_expression,
1403  const Tensor<1, chartdim> & periodicity,
1404  const typename FunctionParser<spacedim>::ConstMap const_map,
1405  const std::string chart_vars,
1406  const std::string space_vars,
1407  const double tolerance,
1408  const double h)
1409  : ChartManifold<dim, spacedim, chartdim>(periodicity)
1410  , const_map(const_map)
1411  , tolerance(tolerance)
1412  , owns_pointers(true)
1413  , push_forward_expression(push_forward_expression)
1414  , pull_back_expression(pull_back_expression)
1415  , chart_vars(chart_vars)
1416  , space_vars(space_vars)
1417  , finite_difference_step(h)
1418 {
1419  FunctionParser<chartdim> *pf = new FunctionParser<chartdim>(spacedim, 0.0, h);
1420  FunctionParser<spacedim> *pb = new FunctionParser<spacedim>(chartdim, 0.0, h);
1423  push_forward_function = pf;
1424  pull_back_function = pb;
1425 }
1426 
1427 
1428 
1429 template <int dim, int spacedim, int chartdim>
1431 {
1432  if (owns_pointers == true)
1433  {
1434  const Function<chartdim> *pf = push_forward_function;
1435  push_forward_function = nullptr;
1436  delete pf;
1437 
1438  const Function<spacedim> *pb = pull_back_function;
1439  pull_back_function = nullptr;
1440  delete pb;
1441  }
1442 }
1443 
1444 
1445 
1446 template <int dim, int spacedim, int chartdim>
1447 std::unique_ptr<Manifold<dim, spacedim>>
1449 {
1450  // This manifold can be constructed either by providing an expression for the
1451  // push forward and the pull back charts, or by providing two Function
1452  // objects. In the first case, the push_forward and pull_back functions are
1453  // created internally in FunctionManifold, and destroyed when this object is
1454  // deleted. In the second case, the function objects are destroyed if they
1455  // are passed as pointers upon construction.
1456  // We need to make sure that our cloned object is constructed in the
1457  // same way this class was constructed, and that its internal Function
1458  // pointers point either to the same Function objects used to construct this
1459  // function or that the newly generated manifold creates internally the
1460  // push_forward and pull_back functions using the same expressions that were
1461  // used to construct this class.
1462  if (!(push_forward_expression.empty() && pull_back_expression.empty()))
1463  {
1464  return std_cxx14::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1465  push_forward_expression,
1466  pull_back_expression,
1467  this->get_periodicity(),
1468  const_map,
1469  chart_vars,
1470  space_vars,
1471  tolerance,
1472  finite_difference_step);
1473  }
1474  else
1475  {
1476  return std_cxx14::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1477  *push_forward_function,
1478  *pull_back_function,
1479  this->get_periodicity(),
1480  tolerance);
1481  }
1482 }
1483 
1484 
1485 
1486 template <int dim, int spacedim, int chartdim>
1489  const Point<chartdim> &chart_point) const
1490 {
1491  Vector<double> pf(spacedim);
1492  Point<spacedim> result;
1493  push_forward_function->vector_value(chart_point, pf);
1494  for (unsigned int i = 0; i < spacedim; ++i)
1495  result[i] = pf[i];
1496 
1497 #ifdef DEBUG
1498  Vector<double> pb(chartdim);
1499  pull_back_function->vector_value(result, pb);
1500  for (unsigned int i = 0; i < chartdim; ++i)
1501  Assert(
1502  (chart_point.norm() > tolerance &&
1503  (std::abs(pb[i] - chart_point[i]) < tolerance * chart_point.norm())) ||
1504  (std::abs(pb[i] - chart_point[i]) < tolerance),
1505  ExcMessage(
1506  "The push forward is not the inverse of the pull back! Bailing out."));
1507 #endif
1508 
1509  return result;
1510 }
1511 
1512 
1513 
1514 template <int dim, int spacedim, int chartdim>
1517  const Point<chartdim> &chart_point) const
1518 {
1520  for (unsigned int i = 0; i < spacedim; ++i)
1521  {
1522  const auto gradient = push_forward_function->gradient(chart_point, i);
1523  for (unsigned int j = 0; j < chartdim; ++j)
1524  DF[i][j] = gradient[j];
1525  }
1526  return DF;
1527 }
1528 
1529 
1530 
1531 template <int dim, int spacedim, int chartdim>
1534  const Point<spacedim> &space_point) const
1535 {
1536  Vector<double> pb(chartdim);
1537  Point<chartdim> result;
1538  pull_back_function->vector_value(space_point, pb);
1539  for (unsigned int i = 0; i < chartdim; ++i)
1540  result[i] = pb[i];
1541  return result;
1542 }
1543 
1544 
1545 
1546 // ============================================================
1547 // TorusManifold
1548 // ============================================================
1549 template <int dim>
1550 Point<3>
1552 {
1553  double x = p(0);
1554  double z = p(1);
1555  double y = p(2);
1556  double phi = std::atan2(y, x);
1557  double theta = std::atan2(z, std::sqrt(x * x + y * y) - R);
1558  double w = std::sqrt(std::pow(y - std::sin(phi) * R, 2.0) +
1559  std::pow(x - std::cos(phi) * R, 2.0) + z * z) /
1560  r;
1561  return {phi, theta, w};
1562 }
1563 
1564 
1565 
1566 template <int dim>
1567 Point<3>
1569 {
1570  double phi = chart_point(0);
1571  double theta = chart_point(1);
1572  double w = chart_point(2);
1573 
1574  return {std::cos(phi) * R + r * w * std::cos(theta) * std::cos(phi),
1575  r * w * std::sin(theta),
1576  std::sin(phi) * R + r * w * std::cos(theta) * std::sin(phi)};
1577 }
1578 
1579 
1580 
1581 template <int dim>
1582 TorusManifold<dim>::TorusManifold(const double R, const double r)
1583  : ChartManifold<dim, 3, 3>(Point<3>(2 * numbers::PI, 2 * numbers::PI, 0.0))
1584  , r(r)
1585  , R(R)
1586 {
1587  Assert(R > r,
1588  ExcMessage("Outer radius R must be greater than the inner "
1589  "radius r."));
1590  Assert(r > 0.0, ExcMessage("inner radius must be positive."));
1591 }
1592 
1593 
1594 
1595 template <int dim>
1596 std::unique_ptr<Manifold<dim, 3>>
1598 {
1599  return std_cxx14::make_unique<TorusManifold<dim>>(R, r);
1600 }
1601 
1602 
1603 
1604 template <int dim>
1607 {
1609 
1610  double phi = chart_point(0);
1611  double theta = chart_point(1);
1612  double w = chart_point(2);
1613 
1614  DX[0][0] = -std::sin(phi) * R - r * w * std::cos(theta) * std::sin(phi);
1615  DX[0][1] = -r * w * std::sin(theta) * std::cos(phi);
1616  DX[0][2] = r * std::cos(theta) * std::cos(phi);
1617 
1618  DX[1][0] = 0;
1619  DX[1][1] = r * w * std::cos(theta);
1620  DX[1][2] = r * std::sin(theta);
1621 
1622  DX[2][0] = std::cos(phi) * R + r * w * std::cos(theta) * std::cos(phi);
1623  DX[2][1] = -r * w * std::sin(theta) * std::sin(phi);
1624  DX[2][2] = r * std::cos(theta) * std::sin(phi);
1625 
1626  return DX;
1627 }
1628 
1629 
1630 
1631 // ============================================================
1632 // TransfiniteInterpolationManifold
1633 // ============================================================
1634 template <int dim, int spacedim>
1637  : triangulation(nullptr)
1638  , level_coarse(-1)
1639 {
1640  AssertThrow(dim > 1, ExcNotImplemented());
1641 }
1642 
1643 
1644 
1645 template <int dim, int spacedim>
1647  spacedim>::~TransfiniteInterpolationManifold()
1648 {
1649  if (clear_signal.connected())
1650  clear_signal.disconnect();
1651 }
1652 
1653 
1654 
1655 template <int dim, int spacedim>
1656 std::unique_ptr<Manifold<dim, spacedim>>
1658 {
1660  if (triangulation)
1661  ptr->initialize(*triangulation);
1662  return std::unique_ptr<Manifold<dim, spacedim>>(ptr);
1663 }
1664 
1665 
1666 
1667 template <int dim, int spacedim>
1668 void
1671 {
1672  this->triangulation = &triangulation;
1673  // in case the triangulatoin is cleared, remove the pointers by a signal
1674  clear_signal = triangulation.signals.clear.connect([&]() -> void {
1675  this->triangulation = nullptr;
1676  this->level_coarse = -1;
1677  });
1678  level_coarse = triangulation.last()->level();
1679  coarse_cell_is_flat.resize(triangulation.n_cells(level_coarse), false);
1681  cell = triangulation.begin(level_coarse),
1682  endc = triangulation.end(level_coarse);
1683  for (; cell != endc; ++cell)
1684  {
1685  bool cell_is_flat = true;
1686  for (unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
1687  if (cell->line(l)->manifold_id() != cell->manifold_id() &&
1688  cell->line(l)->manifold_id() != numbers::flat_manifold_id)
1689  cell_is_flat = false;
1690  if (dim > 2)
1691  for (unsigned int q = 0; q < GeometryInfo<dim>::quads_per_cell; ++q)
1692  if (cell->quad(q)->manifold_id() != cell->manifold_id() &&
1693  cell->quad(q)->manifold_id() != numbers::flat_manifold_id)
1694  cell_is_flat = false;
1695  AssertIndexRange(static_cast<unsigned int>(cell->index()),
1696  coarse_cell_is_flat.size());
1697  coarse_cell_is_flat[cell->index()] = cell_is_flat;
1698  }
1699 }
1700 
1701 
1702 
1703 namespace
1704 {
1705  // version for 1D
1706  template <typename AccessorType>
1708  compute_transfinite_interpolation(const AccessorType &cell,
1709  const Point<1> & chart_point,
1710  const bool /*cell_is_flat*/)
1711  {
1712  return cell.vertex(0) * (1. - chart_point[0]) +
1713  cell.vertex(1) * chart_point[0];
1714  }
1715 
1716  // version for 2D
1717  template <typename AccessorType>
1719  compute_transfinite_interpolation(const AccessorType &cell,
1720  const Point<2> & chart_point,
1721  const bool cell_is_flat)
1722  {
1723  const unsigned int dim = AccessorType::dimension;
1724  const unsigned int spacedim = AccessorType::space_dimension;
1725  const types::manifold_id my_manifold_id = cell.manifold_id();
1726  const Triangulation<dim, spacedim> &tria = cell.get_triangulation();
1727 
1728  // formula see wikipedia
1729  // https://en.wikipedia.org/wiki/Transfinite_interpolation
1730  // S(u,v) = (1-v)c_1(u)+v c_3(u) + (1-u)c_2(v) + u c_4(v) -
1731  // [(1-u)(1-v)P_0 + u(1-v) P_1 + (1-u)v P_2 + uv P_3]
1732  const std::array<Point<spacedim>, 4> vertices{
1733  {cell.vertex(0), cell.vertex(1), cell.vertex(2), cell.vertex(3)}};
1734 
1735  // this evaluates all bilinear shape functions because we need them
1736  // repeatedly. we will update this values in the complicated case with
1737  // curved lines below
1738  std::array<double, 4> weights_vertices{
1739  {(1. - chart_point[0]) * (1. - chart_point[1]),
1740  chart_point[0] * (1. - chart_point[1]),
1741  (1. - chart_point[0]) * chart_point[1],
1742  chart_point[0] * chart_point[1]}};
1743 
1744  Point<spacedim> new_point;
1745  if (cell_is_flat)
1746  for (const unsigned int v : GeometryInfo<2>::vertex_indices())
1747  new_point += weights_vertices[v] * vertices[v];
1748  else
1749  {
1750  // The second line in the formula tells us to subtract the
1751  // contribution of the vertices. If a line employs the same manifold
1752  // as the cell, we can merge the weights of the line with the weights
1753  // of the vertex with a negative sign while going through the faces
1754  // (this is a bit artificial in 2D but it becomes clear in 3D where we
1755  // avoid looking at the faces' orientation and other complications).
1756 
1757  // add the contribution from the lines around the cell (first line in
1758  // formula)
1759  std::array<double, GeometryInfo<2>::vertices_per_face> weights;
1760  std::array<Point<spacedim>, GeometryInfo<2>::vertices_per_face> points;
1761  // note that the views are immutable, but the arrays are not
1762  const auto weights_view =
1763  make_array_view(weights.begin(), weights.end());
1764  const auto points_view = make_array_view(points.begin(), points.end());
1765 
1766  for (unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1767  ++line)
1768  {
1769  const double my_weight =
1770  (line % 2) ? chart_point[line / 2] : 1 - chart_point[line / 2];
1771  const double line_point = chart_point[1 - line / 2];
1772 
1773  // Same manifold or invalid id which will go back to the same
1774  // class -> contribution should be added for the final point,
1775  // which means that we subtract the current weight from the
1776  // negative weight applied to the vertex
1777  const types::manifold_id line_manifold_id =
1778  cell.line(line)->manifold_id();
1779  if (line_manifold_id == my_manifold_id ||
1780  line_manifold_id == numbers::flat_manifold_id)
1781  {
1782  weights_vertices[GeometryInfo<2>::line_to_cell_vertices(line,
1783  0)] -=
1784  my_weight * (1. - line_point);
1785  weights_vertices[GeometryInfo<2>::line_to_cell_vertices(line,
1786  1)] -=
1787  my_weight * line_point;
1788  }
1789  else
1790  {
1791  points[0] =
1793  points[1] =
1795  weights[0] = 1. - line_point;
1796  weights[1] = line_point;
1797  new_point +=
1798  my_weight * tria.get_manifold(line_manifold_id)
1799  .get_new_point(points_view, weights_view);
1800  }
1801  }
1802 
1803  // subtract contribution from the vertices (second line in formula)
1804  for (const unsigned int v : GeometryInfo<2>::vertex_indices())
1805  new_point -= weights_vertices[v] * vertices[v];
1806  }
1807 
1808  return new_point;
1809  }
1810 
1811  // this is replicated from GeometryInfo::face_to_cell_vertices since we need
1812  // it very often in compute_transfinite_interpolation and the function is
1813  // performance critical
1814  static constexpr unsigned int face_to_cell_vertices_3d[6][4] = {{0, 2, 4, 6},
1815  {1, 3, 5, 7},
1816  {0, 4, 1, 5},
1817  {2, 6, 3, 7},
1818  {0, 1, 2, 3},
1819  {4, 5, 6, 7}};
1820 
1821  // this is replicated from GeometryInfo::face_to_cell_lines since we need it
1822  // very often in compute_transfinite_interpolation and the function is
1823  // performance critical
1824  static constexpr unsigned int face_to_cell_lines_3d[6][4] = {{8, 10, 0, 4},
1825  {9, 11, 1, 5},
1826  {2, 6, 8, 9},
1827  {3, 7, 10, 11},
1828  {0, 1, 2, 3},
1829  {4, 5, 6, 7}};
1830 
1831  // version for 3D
1832  template <typename AccessorType>
1834  compute_transfinite_interpolation(const AccessorType &cell,
1835  const Point<3> & chart_point,
1836  const bool cell_is_flat)
1837  {
1838  const unsigned int dim = AccessorType::dimension;
1839  const unsigned int spacedim = AccessorType::space_dimension;
1840  const types::manifold_id my_manifold_id = cell.manifold_id();
1841  const Triangulation<dim, spacedim> &tria = cell.get_triangulation();
1842 
1843  // Same approach as in 2D, but adding the faces, subtracting the edges, and
1844  // adding the vertices
1845  const std::array<Point<spacedim>, 8> vertices{{cell.vertex(0),
1846  cell.vertex(1),
1847  cell.vertex(2),
1848  cell.vertex(3),
1849  cell.vertex(4),
1850  cell.vertex(5),
1851  cell.vertex(6),
1852  cell.vertex(7)}};
1853 
1854  // store the components of the linear shape functions because we need them
1855  // repeatedly. we allow for 10 such shape functions to wrap around the
1856  // first four once again for easier face access.
1857  double linear_shapes[10];
1858  for (unsigned int d = 0; d < 3; ++d)
1859  {
1860  linear_shapes[2 * d] = 1. - chart_point[d];
1861  linear_shapes[2 * d + 1] = chart_point[d];
1862  }
1863 
1864  // wrap linear shape functions around for access in face loop
1865  for (unsigned int d = 6; d < 10; ++d)
1866  linear_shapes[d] = linear_shapes[d - 6];
1867 
1868  std::array<double, 8> weights_vertices;
1869  for (unsigned int i2 = 0, v = 0; i2 < 2; ++i2)
1870  for (unsigned int i1 = 0; i1 < 2; ++i1)
1871  for (unsigned int i0 = 0; i0 < 2; ++i0, ++v)
1872  weights_vertices[v] =
1873  (linear_shapes[4 + i2] * linear_shapes[2 + i1]) * linear_shapes[i0];
1874 
1875  Point<spacedim> new_point;
1876  if (cell_is_flat)
1877  for (unsigned int v = 0; v < 8; ++v)
1878  new_point += weights_vertices[v] * vertices[v];
1879  else
1880  {
1881  // identify the weights for the lines to be accumulated (vertex
1882  // weights are set outside and coincide with the flat manifold case)
1883 
1884  std::array<double, GeometryInfo<3>::lines_per_cell> weights_lines;
1885  std::fill(weights_lines.begin(), weights_lines.end(), 0.0);
1886 
1887  // start with the contributions of the faces
1888  std::array<double, GeometryInfo<2>::vertices_per_cell> weights;
1889  std::array<Point<spacedim>, GeometryInfo<2>::vertices_per_cell> points;
1890  // note that the views are immutable, but the arrays are not
1891  const auto weights_view =
1892  make_array_view(weights.begin(), weights.end());
1893  const auto points_view = make_array_view(points.begin(), points.end());
1894 
1895  for (const unsigned int face : GeometryInfo<3>::face_indices())
1896  {
1897  const double my_weight = linear_shapes[face];
1898  const unsigned int face_even = face - face % 2;
1899 
1900  if (std::abs(my_weight) < 1e-13)
1901  continue;
1902 
1903  // same manifold or invalid id which will go back to the same class
1904  // -> face will interpolate from the surrounding lines and vertices
1905  const types::manifold_id face_manifold_id =
1906  cell.face(face)->manifold_id();
1907  if (face_manifold_id == my_manifold_id ||
1908  face_manifold_id == numbers::flat_manifold_id)
1909  {
1910  for (unsigned int line = 0;
1911  line < GeometryInfo<2>::lines_per_cell;
1912  ++line)
1913  {
1914  const double line_weight =
1915  linear_shapes[face_even + 2 + line];
1916  weights_lines[face_to_cell_lines_3d[face][line]] +=
1917  my_weight * line_weight;
1918  }
1919  // as to the indices inside linear_shapes: we use the index
1920  // wrapped around at 2*d, ensuring the correct orientation of
1921  // the face's coordinate system with respect to the
1922  // lexicographic indices
1923  weights_vertices[face_to_cell_vertices_3d[face][0]] -=
1924  linear_shapes[face_even + 2] *
1925  (linear_shapes[face_even + 4] * my_weight);
1926  weights_vertices[face_to_cell_vertices_3d[face][1]] -=
1927  linear_shapes[face_even + 3] *
1928  (linear_shapes[face_even + 4] * my_weight);
1929  weights_vertices[face_to_cell_vertices_3d[face][2]] -=
1930  linear_shapes[face_even + 2] *
1931  (linear_shapes[face_even + 5] * my_weight);
1932  weights_vertices[face_to_cell_vertices_3d[face][3]] -=
1933  linear_shapes[face_even + 3] *
1934  (linear_shapes[face_even + 5] * my_weight);
1935  }
1936  else
1937  {
1938  for (const unsigned int v : GeometryInfo<2>::vertex_indices())
1939  points[v] = vertices[face_to_cell_vertices_3d[face][v]];
1940  weights[0] =
1941  linear_shapes[face_even + 2] * linear_shapes[face_even + 4];
1942  weights[1] =
1943  linear_shapes[face_even + 3] * linear_shapes[face_even + 4];
1944  weights[2] =
1945  linear_shapes[face_even + 2] * linear_shapes[face_even + 5];
1946  weights[3] =
1947  linear_shapes[face_even + 3] * linear_shapes[face_even + 5];
1948  new_point +=
1949  my_weight * tria.get_manifold(face_manifold_id)
1950  .get_new_point(points_view, weights_view);
1951  }
1952  }
1953 
1954  // next subtract the contributions of the lines
1955  const auto weights_view_line =
1956  make_array_view(weights.begin(), weights.begin() + 2);
1957  const auto points_view_line =
1958  make_array_view(points.begin(), points.begin() + 2);
1959  for (unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
1960  ++line)
1961  {
1962  const double line_point =
1963  (line < 8 ? chart_point[1 - (line % 4) / 2] : chart_point[2]);
1964  double my_weight = 0.;
1965  if (line < 8)
1966  my_weight = linear_shapes[line % 4] * linear_shapes[4 + line / 4];
1967  else
1968  {
1969  const unsigned int subline = line - 8;
1970  my_weight =
1971  linear_shapes[subline % 2] * linear_shapes[2 + subline / 2];
1972  }
1973  my_weight -= weights_lines[line];
1974 
1975  if (std::abs(my_weight) < 1e-13)
1976  continue;
1977 
1978  const types::manifold_id line_manifold_id =
1979  cell.line(line)->manifold_id();
1980  if (line_manifold_id == my_manifold_id ||
1981  line_manifold_id == numbers::flat_manifold_id)
1982  {
1983  weights_vertices[GeometryInfo<3>::line_to_cell_vertices(line,
1984  0)] -=
1985  my_weight * (1. - line_point);
1986  weights_vertices[GeometryInfo<3>::line_to_cell_vertices(line,
1987  1)] -=
1988  my_weight * (line_point);
1989  }
1990  else
1991  {
1992  points[0] =
1994  points[1] =
1996  weights[0] = 1. - line_point;
1997  weights[1] = line_point;
1998  new_point -= my_weight * tria.get_manifold(line_manifold_id)
1999  .get_new_point(points_view_line,
2000  weights_view_line);
2001  }
2002  }
2003 
2004  // finally add the contribution of the
2005  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
2006  new_point += weights_vertices[v] * vertices[v];
2007  }
2008  return new_point;
2009  }
2010 } // namespace
2011 
2012 
2013 
2014 template <int dim, int spacedim>
2017  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
2018  const Point<dim> & chart_point) const
2019 {
2020  AssertDimension(cell->level(), level_coarse);
2021 
2022  // check that the point is in the unit cell which is the current chart
2023  // Tolerance 5e-4 chosen that the method also works with manifolds
2024  // that have some discretization error like SphericalManifold
2026  ExcMessage("chart_point is not in unit interval"));
2027 
2028  return compute_transfinite_interpolation(*cell,
2029  chart_point,
2030  coarse_cell_is_flat[cell->index()]);
2031 }
2032 
2033 
2034 
2035 template <int dim, int spacedim>
2038  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
2039  const Point<dim> & chart_point,
2040  const Point<spacedim> &pushed_forward_chart_point) const
2041 {
2042  // compute the derivative with the help of finite differences
2044  for (unsigned int d = 0; d < dim; ++d)
2045  {
2046  Point<dim> modified = chart_point;
2047  const double step = chart_point[d] > 0.5 ? -1e-8 : 1e-8;
2048 
2049  // avoid checking outside of the unit interval
2050  modified[d] += step;
2051  Tensor<1, spacedim> difference =
2052  compute_transfinite_interpolation(*cell,
2053  modified,
2054  coarse_cell_is_flat[cell->index()]) -
2055  pushed_forward_chart_point;
2056  for (unsigned int e = 0; e < spacedim; ++e)
2057  grad[e][d] = difference[e] / step;
2058  }
2059  return grad;
2060 }
2061 
2062 
2063 
2064 template <int dim, int spacedim>
2065 Point<dim>
2067  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
2068  const Point<spacedim> & point,
2069  const Point<dim> &initial_guess) const
2070 {
2071  Point<dim> outside;
2072  for (unsigned int d = 0; d < dim; ++d)
2074 
2075  // project the user-given input to unit cell
2076  Point<dim> chart_point =
2078 
2079  // run quasi-Newton iteration with a combination of finite differences for
2080  // the exact Jacobian and "Broyden's good method". As opposed to the various
2081  // mapping implementations, this class does not throw exception upon failure
2082  // as those are relatively expensive and failure occurs quite regularly in
2083  // the implementation of the compute_chart_points method.
2084  Tensor<1, spacedim> residual =
2085  point -
2086  compute_transfinite_interpolation(*cell,
2087  chart_point,
2088  coarse_cell_is_flat[cell->index()]);
2089  const double tolerance = 1e-21 * Utilities::fixed_power<2>(cell->diameter());
2090  double residual_norm_square = residual.norm_square();
2092  bool must_recompute_jacobian = true;
2093  for (unsigned int i = 0; i < 100; ++i)
2094  {
2095  if (residual_norm_square < tolerance)
2096  {
2097  // do a final update of the point with the last available Jacobian
2098  // information. The residual is close to zero due to the check
2099  // above, but me might improve some of the last digits by a final
2100  // Newton-like step with step length 1
2101  Tensor<1, dim> update;
2102  for (unsigned int d = 0; d < spacedim; ++d)
2103  for (unsigned int e = 0; e < dim; ++e)
2104  update[e] += inv_grad[d][e] * residual[d];
2105  return chart_point + update;
2106  }
2107 
2108  // every 9 iterations, including the first time around, we create an
2109  // approximation of the Jacobian with finite differences. Broyden's
2110  // method usually does not need more than 5-8 iterations, but sometimes
2111  // we might have had a bad initial guess and then we can accelerate
2112  // convergence considerably with getting the actual Jacobian rather than
2113  // using secant-like methods (one gradient calculation in 3D costs as
2114  // much as 3 more iterations). this usually happens close to convergence
2115  // and one more step with the finite-differenced Jacobian leads to
2116  // convergence
2117  if (must_recompute_jacobian || i % 9 == 0)
2118  {
2119  // if the determinant is zero or negative, the mapping is either not
2120  // invertible or already has inverted and we are outside the valid
2121  // chart region. Note that the Jacobian here represents the
2122  // derivative of the forward map and should have a positive
2123  // determinant since we use properly oriented meshes.
2125  push_forward_gradient(cell,
2126  chart_point,
2127  Point<spacedim>(point - residual));
2128  if (grad.determinant() <= 0.0)
2129  return outside;
2130  inv_grad = grad.covariant_form();
2131  must_recompute_jacobian = false;
2132  }
2133  Tensor<1, dim> update;
2134  for (unsigned int d = 0; d < spacedim; ++d)
2135  for (unsigned int e = 0; e < dim; ++e)
2136  update[e] += inv_grad[d][e] * residual[d];
2137 
2138  // Line search, accept step if the residual has decreased
2139  double alpha = 1.;
2140 
2141  // check if point is inside 1.2 times the unit cell to avoid
2142  // hitting points very far away from valid ones in the manifolds
2143  while (
2144  !GeometryInfo<dim>::is_inside_unit_cell(chart_point + alpha * update,
2145  0.2) &&
2146  alpha > 1e-7)
2147  alpha *= 0.5;
2148 
2149  const Tensor<1, spacedim> old_residual = residual;
2150  while (alpha > 1e-4)
2151  {
2152  Point<dim> guess = chart_point + alpha * update;
2153  residual =
2154  point - compute_transfinite_interpolation(
2155  *cell, guess, coarse_cell_is_flat[cell->index()]);
2156  const double residual_norm_new = residual.norm_square();
2157  if (residual_norm_new < residual_norm_square)
2158  {
2159  residual_norm_square = residual_norm_new;
2160  chart_point += alpha * update;
2161  break;
2162  }
2163  else
2164  alpha *= 0.5;
2165  }
2166  // If alpha got very small, it is likely due to a bad Jacobian
2167  // approximation with Broyden's method (relatively far away from the
2168  // zero), which can be corrected by the outer loop when a Newton update
2169  // is recomputed. The second case is when the Jacobian is actually bad
2170  // and we should fail as early as possible. Since we cannot really
2171  // distinguish the two, we must continue here in any case.
2172  if (alpha <= 1e-4)
2173  must_recompute_jacobian = true;
2174 
2175  // update the inverse Jacobian with "Broyden's good method" and
2176  // Sherman-Morrison formula for the update of the inverse, see
2177  // https://en.wikipedia.org/wiki/Broyden%27s_method
2178  // J^{-1}_n = J^{-1}_{n-1} + (delta x_n - J^{-1}_{n-1} delta f_n) /
2179  // (delta x_n^T J_{-1}_{n-1} delta f_n) delta x_n^T J^{-1}_{n-1}
2180 
2181  // switch sign in residual as compared to the formula above because we
2182  // use a negative definition of the residual with respect to the
2183  // Jacobian
2184  const Tensor<1, spacedim> delta_f = old_residual - residual;
2185 
2186  Tensor<1, dim> Jinv_deltaf;
2187  for (unsigned int d = 0; d < spacedim; ++d)
2188  for (unsigned int e = 0; e < dim; ++e)
2189  Jinv_deltaf[e] += inv_grad[d][e] * delta_f[d];
2190 
2191  const Tensor<1, dim> delta_x = alpha * update;
2192 
2193  // prevent division by zero. This number should be scale-invariant
2194  // because Jinv_deltaf carries no units and x is in reference
2195  // coordinates.
2196  if (std::abs(delta_x * Jinv_deltaf) > 1e-12)
2197  {
2198  const Tensor<1, dim> factor =
2199  (delta_x - Jinv_deltaf) / (delta_x * Jinv_deltaf);
2200  Tensor<1, spacedim> jac_update;
2201  for (unsigned int d = 0; d < spacedim; ++d)
2202  for (unsigned int e = 0; e < dim; ++e)
2203  jac_update[d] += delta_x[e] * inv_grad[d][e];
2204  for (unsigned int d = 0; d < spacedim; ++d)
2205  for (unsigned int e = 0; e < dim; ++e)
2206  inv_grad[d][e] += factor[e] * jac_update[d];
2207  }
2208  }
2209  return outside;
2210 }
2211 
2212 
2213 
2214 template <int dim, int spacedim>
2215 std::array<unsigned int, 20>
2218  const ArrayView<const Point<spacedim>> &points) const
2219 {
2220  // The methods to identify cells around points in GridTools are all written
2221  // for the active cells, but we are here looking at some cells at the coarse
2222  // level.
2223  Assert(triangulation != nullptr, ExcNotInitialized());
2224  Assert(triangulation->begin_active()->level() >= level_coarse,
2225  ExcMessage("The manifold was initialized with level " +
2226  std::to_string(level_coarse) + " but there are now" +
2227  "active cells on a lower level. Coarsening the mesh is " +
2228  "currently not supported"));
2229 
2230  // This computes the distance of the surrounding points transformed to the
2231  // unit cell from the unit cell.
2233  triangulation->begin(
2234  level_coarse),
2235  endc =
2236  triangulation->end(
2237  level_coarse);
2238  boost::container::small_vector<std::pair<double, unsigned int>, 200>
2239  distances_and_cells;
2240  for (; cell != endc; ++cell)
2241  {
2242  // only consider cells where the current manifold is attached
2243  if (&cell->get_manifold() != this)
2244  continue;
2245 
2246  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
2247  vertices;
2248  for (const unsigned int vertex_n : GeometryInfo<dim>::vertex_indices())
2249  {
2250  vertices[vertex_n] = cell->vertex(vertex_n);
2251  }
2252 
2253  // cheap check: if any of the points is not inside a circle around the
2254  // center of the loop, we can skip the expensive part below (this assumes
2255  // that the manifold does not deform the grid too much)
2257  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
2258  center += vertices[v];
2260  double radius_square = 0.;
2261  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
2262  radius_square =
2263  std::max(radius_square, (center - vertices[v]).norm_square());
2264  bool inside_circle = true;
2265  for (unsigned int i = 0; i < points.size(); ++i)
2266  if ((center - points[i]).norm_square() > radius_square * 1.5)
2267  {
2268  inside_circle = false;
2269  break;
2270  }
2271  if (inside_circle == false)
2272  continue;
2273 
2274  // slightly more expensive search
2275  double current_distance = 0;
2276  for (unsigned int i = 0; i < points.size(); ++i)
2277  {
2278  Point<dim> point =
2279  cell->real_to_unit_cell_affine_approximation(points[i]);
2280  current_distance += GeometryInfo<dim>::distance_to_unit_cell(point);
2281  }
2282  distances_and_cells.push_back(
2283  std::make_pair(current_distance, cell->index()));
2284  }
2285  // no coarse cell could be found -> transformation failed
2286  AssertThrow(distances_and_cells.size() > 0,
2288  std::sort(distances_and_cells.begin(), distances_and_cells.end());
2289  std::array<unsigned int, 20> cells;
2290  cells.fill(numbers::invalid_unsigned_int);
2291  for (unsigned int i = 0; i < distances_and_cells.size() && i < cells.size();
2292  ++i)
2293  cells[i] = distances_and_cells[i].second;
2294 
2295  return cells;
2296 }
2297 
2298 
2299 
2300 template <int dim, int spacedim>
2303  const ArrayView<const Point<spacedim>> &surrounding_points,
2304  ArrayView<Point<dim>> chart_points) const
2305 {
2306  Assert(surrounding_points.size() == chart_points.size(),
2307  ExcMessage("The chart points array view must be as large as the "
2308  "surrounding points array view."));
2309 
2310  std::array<unsigned int, 20> nearby_cells =
2311  get_possible_cells_around_points(surrounding_points);
2312 
2313  // This function is nearly always called to place new points on a cell or
2314  // cell face. In this case, the general structure of the surrounding points
2315  // is known (i.e., if there are eight surrounding points, then they will
2316  // almost surely be either eight points around a quadrilateral or the eight
2317  // vertices of a cube). Hence, making this assumption, we use two
2318  // optimizations (one for structdim == 2 and one for structdim == 3) that
2319  // guess the locations of some of the chart points more efficiently than the
2320  // affine map approximation. The affine map approximation is used whenever
2321  // we don't have a cheaper guess available.
2322 
2323  // Function that can guess the location of a chart point by assuming that
2324  // the eight surrounding points are points on a two-dimensional object
2325  // (either a cell in 2D or the face of a hexahedron in 3D), arranged like
2326  //
2327  // 2 - 7 - 3
2328  // | |
2329  // 4 5
2330  // | |
2331  // 0 - 6 - 1
2332  //
2333  // This function assumes that the first three chart points have been
2334  // computed since there is no effective way to guess them.
2335  auto guess_chart_point_structdim_2 = [&](const unsigned int i) -> Point<dim> {
2336  Assert(surrounding_points.size() == 8 && 2 < i && i < 8,
2337  ExcMessage("This function assumes that there are eight surrounding "
2338  "points around a two-dimensional object. It also assumes "
2339  "that the first three chart points have already been "
2340  "computed."));
2341  switch (i)
2342  {
2343  case 0:
2344  case 1:
2345  case 2:
2346  Assert(false, ExcInternalError());
2347  break;
2348  case 3:
2349  return chart_points[1] + (chart_points[2] - chart_points[0]);
2350  case 4:
2351  return 0.5 * (chart_points[0] + chart_points[2]);
2352  case 5:
2353  return 0.5 * (chart_points[1] + chart_points[3]);
2354  case 6:
2355  return 0.5 * (chart_points[0] + chart_points[1]);
2356  case 7:
2357  return 0.5 * (chart_points[2] + chart_points[3]);
2358  default:
2359  Assert(false, ExcInternalError());
2360  }
2361 
2362  return Point<dim>();
2363  };
2364 
2365  // Function that can guess the location of a chart point by assuming that
2366  // the eight surrounding points form the vertices of a hexahedron, arranged
2367  // like
2368  //
2369  // 6-------7
2370  // /| /|
2371  // / / |
2372  // / | / |
2373  // 4-------5 |
2374  // | 2- -|- -3
2375  // | / | /
2376  // | | /
2377  // |/ |/
2378  // 0-------1
2379  //
2380  // (where vertex 2 is the back left vertex) we can estimate where chart
2381  // points 5 - 7 are by computing the height (in chart coordinates) as c4 -
2382  // c0 and then adding that onto the appropriate bottom vertex.
2383  //
2384  // This function assumes that the first five chart points have been computed
2385  // since there is no effective way to guess them.
2386  auto guess_chart_point_structdim_3 = [&](const unsigned int i) -> Point<dim> {
2387  Assert(surrounding_points.size() == 8 && 4 < i && i < 8,
2388  ExcMessage("This function assumes that there are eight surrounding "
2389  "points around a three-dimensional object. It also "
2390  "assumes that the first five chart points have already "
2391  "been computed."));
2392  return chart_points[i - 4] + (chart_points[4] - chart_points[0]);
2393  };
2394 
2395  // Check if we can use the two chart point shortcuts above before we start:
2396  bool use_structdim_2_guesses = false;
2397  bool use_structdim_3_guesses = false;
2398  // note that in the structdim 2 case: 0 - 6 and 2 - 7 should be roughly
2399  // parallel, while in the structdim 3 case, 0 - 6 and 2 - 7 should be roughly
2400  // orthogonal. Use the angle between these two vectors to figure out if we
2401  // should turn on either structdim optimization.
2402  if (surrounding_points.size() == 8)
2403  {
2404  const Tensor<1, spacedim> v06 =
2405  surrounding_points[6] - surrounding_points[0];
2406  const Tensor<1, spacedim> v27 =
2407  surrounding_points[7] - surrounding_points[2];
2408 
2409  // note that we can save a call to sqrt() by rearranging
2410  const double cosine = scalar_product(v06, v27) /
2411  std::sqrt(v06.norm_square() * v27.norm_square());
2412  if (0.707 < cosine)
2413  // the angle is less than pi/4, so these vectors are roughly parallel:
2414  // enable the structdim 2 optimization
2415  use_structdim_2_guesses = true;
2416  else if (spacedim == 3)
2417  // otherwise these vectors are roughly orthogonal: enable the
2418  // structdim 3 optimization if we are in 3D
2419  use_structdim_3_guesses = true;
2420  }
2421  // we should enable at most one of the optimizations
2422  Assert((!use_structdim_2_guesses && !use_structdim_3_guesses) ||
2423  (use_structdim_2_guesses ^ use_structdim_3_guesses),
2424  ExcInternalError());
2425 
2426 
2427 
2428  auto compute_chart_point =
2429  [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
2430  const unsigned int point_index) {
2431  Point<dim> guess;
2432  // an optimization: keep track of whether or not we used the affine
2433  // approximation so that we don't call pull_back with the same
2434  // initial guess twice (i.e., if pull_back fails the first time,
2435  // don't try again with the same function arguments).
2436  bool used_affine_approximation = false;
2437  // if we have already computed three points, we can guess the fourth
2438  // to be the missing corner point of a rectangle
2439  if (point_index == 3 && surrounding_points.size() >= 8)
2440  guess = chart_points[1] + (chart_points[2] - chart_points[0]);
2441  else if (use_structdim_2_guesses && 3 < point_index)
2442  guess = guess_chart_point_structdim_2(point_index);
2443  else if (use_structdim_3_guesses && 4 < point_index)
2444  guess = guess_chart_point_structdim_3(point_index);
2445  else if (dim == 3 && point_index > 7 && surrounding_points.size() == 26)
2446  {
2447  if (point_index < 20)
2448  guess =
2449  0.5 * (chart_points[GeometryInfo<dim>::line_to_cell_vertices(
2450  point_index - 8, 0)] +
2452  point_index - 8, 1)]);
2453  else
2454  guess =
2455  0.25 * (chart_points[GeometryInfo<dim>::face_to_cell_vertices(
2456  point_index - 20, 0)] +
2458  point_index - 20, 1)] +
2460  point_index - 20, 2)] +
2462  point_index - 20, 3)]);
2463  }
2464  else
2465  {
2466  guess = cell->real_to_unit_cell_affine_approximation(
2467  surrounding_points[point_index]);
2468  used_affine_approximation = true;
2469  }
2470  chart_points[point_index] =
2471  pull_back(cell, surrounding_points[point_index], guess);
2472 
2473  // the initial guess may not have been good enough: if applicable,
2474  // try again with the affine approximation (which is more accurate
2475  // than the cheap methods used above)
2476  if (chart_points[point_index][0] ==
2478  !used_affine_approximation)
2479  {
2480  guess = cell->real_to_unit_cell_affine_approximation(
2481  surrounding_points[point_index]);
2482  chart_points[point_index] =
2483  pull_back(cell, surrounding_points[point_index], guess);
2484  }
2485 
2486  if (chart_points[point_index][0] ==
2488  {
2489  for (unsigned int d = 0; d < dim; ++d)
2490  guess[d] = 0.5;
2491  chart_points[point_index] =
2492  pull_back(cell, surrounding_points[point_index], guess);
2493  }
2494  };
2495 
2496  // check whether all points are inside the unit cell of the current chart
2497  for (unsigned int c = 0; c < nearby_cells.size(); ++c)
2498  {
2500  triangulation, level_coarse, nearby_cells[c]);
2501  bool inside_unit_cell = true;
2502  for (unsigned int i = 0; i < surrounding_points.size(); ++i)
2503  {
2504  compute_chart_point(cell, i);
2505 
2506  // Tolerance 5e-4 chosen that the method also works with manifolds
2507  // that have some discretization error like SphericalManifold
2508  if (GeometryInfo<dim>::is_inside_unit_cell(chart_points[i], 5e-4) ==
2509  false)
2510  {
2511  inside_unit_cell = false;
2512  break;
2513  }
2514  }
2515  if (inside_unit_cell == true)
2516  {
2517  return cell;
2518  }
2519 
2520  // if we did not find a point and this was the last valid cell (the next
2521  // iterate being the end of the array or an invalid tag), we must stop
2522  if (c == nearby_cells.size() - 1 ||
2523  nearby_cells[c + 1] == numbers::invalid_unsigned_int)
2524  {
2525  // generate additional information to help debugging why we did not
2526  // get a point
2527  std::ostringstream message;
2528  for (unsigned int b = 0; b <= c; ++b)
2529  {
2531  triangulation, level_coarse, nearby_cells[b]);
2532  message << "Looking at cell " << cell->id()
2533  << " with vertices: " << std::endl;
2534  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
2535  message << cell->vertex(v) << " ";
2536  message << std::endl;
2537  message << "Transformation to chart coordinates: " << std::endl;
2538  for (unsigned int i = 0; i < surrounding_points.size(); ++i)
2539  {
2540  compute_chart_point(cell, i);
2541  message << surrounding_points[i] << " -> " << chart_points[i]
2542  << std::endl;
2543  }
2544  }
2545 
2546  AssertThrow(false,
2548  message.str())));
2549  }
2550  }
2551 
2552  // a valid inversion should have returned a point above. an invalid
2553  // inversion should have triggered the assertion, so we should never end up
2554  // here
2555  Assert(false, ExcInternalError());
2557 }
2558 
2559 
2560 
2561 template <int dim, int spacedim>
2564  const ArrayView<const Point<spacedim>> &surrounding_points,
2565  const ArrayView<const double> & weights) const
2566 {
2567  boost::container::small_vector<Point<dim>, 100> chart_points(
2568  surrounding_points.size());
2569  ArrayView<Point<dim>> chart_points_view =
2570  make_array_view(chart_points.begin(), chart_points.end());
2571  const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2572 
2573  const Point<dim> p_chart =
2574  chart_manifold.get_new_point(chart_points_view, weights);
2575 
2576  return push_forward(cell, p_chart);
2577 }
2578 
2579 
2580 
2581 template <int dim, int spacedim>
2582 void
2584  const ArrayView<const Point<spacedim>> &surrounding_points,
2585  const Table<2, double> & weights,
2586  ArrayView<Point<spacedim>> new_points) const
2587 {
2588  Assert(weights.size(0) > 0, ExcEmptyObject());
2589  AssertDimension(surrounding_points.size(), weights.size(1));
2590 
2591  boost::container::small_vector<Point<dim>, 100> chart_points(
2592  surrounding_points.size());
2593  ArrayView<Point<dim>> chart_points_view =
2594  make_array_view(chart_points.begin(), chart_points.end());
2595  const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2596 
2597  boost::container::small_vector<Point<dim>, 100> new_points_on_chart(
2598  weights.size(0));
2599  chart_manifold.get_new_points(chart_points_view,
2600  weights,
2601  make_array_view(new_points_on_chart.begin(),
2602  new_points_on_chart.end()));
2603 
2604  for (unsigned int row = 0; row < weights.size(0); ++row)
2605  new_points[row] = push_forward(cell, new_points_on_chart[row]);
2606 }
2607 
2608 
2609 
2610 // explicit instantiations
2611 #include "manifold_lib.inst"
2612 
Physics::Elasticity::Kinematics::epsilon
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
CylindricalManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:1073
internal::invalid_pull_back_coordinate
static constexpr double invalid_pull_back_coordinate
Definition: manifold_lib.cc:40
TransfiniteInterpolationManifold::compute_chart_points
Triangulation< dim, spacedim >::cell_iterator compute_chart_points(const ArrayView< const Point< spacedim >> &surrounding_points, ArrayView< Point< dim >> chart_points) const
Definition: manifold_lib.cc:2302
internal::QGaussLobatto::gamma
long double gamma(const unsigned int n)
Definition: quadrature_lib.cc:96
tria_accessor.h
PolarManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:132
ChartManifold
Definition: manifold.h:951
SphericalManifold::normal_vector
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
Definition: manifold_lib.cc:489
TorusManifold::push_forward_gradient
virtual DerivativeForm< 1, 3, 3 > push_forward_gradient(const Point< 3 > &chart_point) const override
Definition: manifold_lib.cc:1606
PolarManifold::get_periodicity
static Tensor< 1, spacedim > get_periodicity()
Definition: manifold_lib.cc:141
EllipticalManifold::push_forward_gradient
virtual DerivativeForm< 1, spacedim, spacedim > push_forward_gradient(const Point< spacedim > &chart_point) const override
Definition: manifold_lib.cc:1325
TransfiniteInterpolationManifold::get_possible_cells_around_points
std::array< unsigned int, 20 > get_possible_cells_around_points(const ArrayView< const Point< spacedim >> &surrounding_points) const
Definition: manifold_lib.cc:2217
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation::get_triangulation
Triangulation< dim, spacedim > & get_triangulation()
Definition: tria.cc:13338
EllipticalManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:1231
Triangulation
Definition: tria.h:1109
tria.h
Manifold::get_normals_at_vertices
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
Definition: manifold.cc:301
ArrayView
Definition: array_view.h:77
FunctionManifold::const_map
const FunctionParser< spacedim >::ConstMap const_map
Definition: manifold_lib.h:713
EllipticalManifold::pull_back
virtual Point< spacedim > pull_back(const Point< spacedim > &space_point) const override
Definition: manifold_lib.cc:1283
SymmetricTensor::scalar_product
constexpr ProductType< Number, OtherNumber >::type scalar_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, OtherNumber > &t2)
Definition: symmetric_tensor.h:3749
FunctionManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:1448
PolarManifold
Definition: manifold_lib.h:63
tria_iterator.h
TorusManifold::r
double r
Definition: manifold_lib.h:826
SphericalManifold::get_new_point
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &vertices, const ArrayView< const double > &weights) const override
Definition: manifold_lib.cc:585
GeometryInfo::project_to_unit_cell
static Point< dim > project_to_unit_cell(const Point< dim > &p)
FunctionManifold::pull_back_function
SmartPointer< const Function< spacedim >, FunctionManifold< dim, spacedim, chartdim > > pull_back_function
Definition: manifold_lib.h:727
GeometryInfo
Definition: geometry_info.h:1224
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Patterns::Tools::to_string
std::string to_string(const T &t)
Definition: patterns.h:2360
SymmetricTensor::invert
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &t)
Definition: symmetric_tensor.h:3467
SphericalManifold
Definition: manifold_lib.h:231
Physics::Elasticity::Kinematics::e
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
ArrayView::size
std::size_t size() const
Definition: array_view.h:484
PolarManifold::normal_vector
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
Definition: manifold_lib.cc:320
FunctionManifold::chart_vars
const std::string chart_vars
Definition: manifold_lib.h:758
numbers::flat_manifold_id
const types::manifold_id flat_manifold_id
Definition: types.h:273
FunctionManifold::push_forward_function
SmartPointer< const Function< chartdim >, FunctionManifold< dim, spacedim, chartdim > > push_forward_function
Definition: manifold_lib.h:720
Manifold::normal_vector
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
Definition: manifold.cc:237
second
Point< 2 > second
Definition: grid_out.cc:4353
ChartManifold::get_new_point
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
Definition: manifold.cc:987
TransfiniteInterpolationManifold::initialize
void initialize(const Triangulation< dim, spacedim > &triangulation)
Definition: manifold_lib.cc:1669
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor::cross_product_3d
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:2407
Physics::Transformations::Contravariant::pull_back
Tensor< 1, dim, Number > pull_back(const Tensor< 1, dim, Number > &v, const Tensor< 2, dim, Number > &F)
FunctionManifold::push_forward_gradient
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const override
Definition: manifold_lib.cc:1516
Table< 2, double >
SphericalManifold::get_new_points
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_lib.cc:568
EllipticalManifold
Definition: manifold_lib.h:508
OpenCASCADE::point
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
CylindricalManifold
Definition: manifold_lib.h:387
DerivativeForm
Definition: derivative_form.h:60
FunctionManifold::push_forward_expression
const std::string push_forward_expression
Definition: manifold_lib.h:748
Differentiation::SD::atan2
Expression atan2(const Expression &y, const Expression &x)
Definition: symengine_math.cc:154
TorusManifold::R
double R
Definition: manifold_lib.h:826
Physics::Elasticity::Kinematics::w
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
tensor.h
OpenCASCADE::push_forward
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
Definition: utilities.cc:828
internal::compute_normal
Point< spacedim > compute_normal(const Tensor< 1, spacedim > &, bool=false)
Definition: manifold_lib.cc:77
FunctionManifold::space_vars
const std::string space_vars
Definition: manifold_lib.h:763
GeometryInfo::face_to_cell_vertices
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
TransfiniteInterpolationManifold
Definition: manifold_lib.h:971
TorusManifold::push_forward
virtual Point< 3 > push_forward(const Point< 3 > &chart_point) const override
Definition: manifold_lib.cc:1568
Tensor::norm
numbers::NumberTraits< Number >::real_type norm() const
FunctionParser::ConstMap
std::map< std::string, double > ConstMap
Definition: function_parser.h:285
SphericalManifold::get_tangent_vector
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
Definition: manifold_lib.cc:438
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
double
StandardExceptions::ExcEmptyObject
static ::ExceptionBase & ExcEmptyObject()
TransfiniteInterpolationManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:1657
DerivativeForm::determinant
Number determinant() const
SphericalManifold::clone
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Definition: manifold_lib.cc:363
TransfiniteInterpolationManifold::pull_back
Point< dim > pull_back(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_guess) const
Definition: manifold_lib.cc:2066
Tensor
Definition: tensor.h:450
LocalIntegrators::Divergence::norm
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim >>> &Du)
Definition: divergence.h:548
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Physics::Elasticity::Kinematics::l
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Manifold::FaceVertexNormals
std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face > FaceVertexNormals
Definition: manifold.h:355
TransfiniteInterpolationManifold::get_new_points
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_lib.cc:2583
Tensor::norm_square
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
FunctionManifold::FunctionManifold
FunctionManifold(const Function< chartdim > &push_forward_function, const Function< spacedim > &pull_back_function, const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >(), const double tolerance=1e-10)
Definition: manifold_lib.cc:1360
PolarManifold::push_forward
virtual Point< spacedim > push_forward(const Point< spacedim > &chart_point) const override
Definition: manifold_lib.cc:158
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
StandardExceptions::ExcNotInitialized
static ::ExceptionBase & ExcNotInitialized()
FunctionManifold::pull_back_expression
const std::string pull_back_expression
Definition: manifold_lib.h:753
GeometryInfo::line_to_cell_vertices
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
numbers
Definition: numbers.h:207
EllipticalManifold::direction
Tensor< 1, spacedim > direction
Definition: manifold_lib.h:553
EllipticalManifold::cosh_u
const double cosh_u
Definition: manifold_lib.h:561
TransfiniteInterpolationManifold::push_forward
Point< spacedim > push_forward(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &chart_point) const
Definition: manifold_lib.cc:2016
TriaActiveIterator
Definition: tria_iterator.h:759
TableBase::size
size_type size(const unsigned int i) const
CylindricalManifold::pull_back
virtual Point< 3 > pull_back(const Point< spacedim > &space_point) const override
Definition: manifold_lib.cc:1112
ArrayView::make_array_view
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:607
FunctionParser::initialize
void initialize(const std::string &vars, const std::vector< std::string > &expressions, const ConstMap &constants, const bool time_dependent=false)
Definition: function_parser.cc:94
manifold_lib.h
EllipticalManifold::get_periodicity
static Tensor< 1, spacedim > get_periodicity()
Definition: manifold_lib.cc:1242
CylindricalManifold::push_forward_gradient
virtual DerivativeForm< 1, 3, spacedim > push_forward_gradient(const Point< 3 > &chart_point) const override
Definition: manifold_lib.cc:1159
unsigned int
Triangulation::cell_iterator
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1343
EllipticalManifold::EllipticalManifold
EllipticalManifold(const Point< spacedim > &center, const Tensor< 1, spacedim > &major_axis_direction, const double eccentricity)
Definition: manifold_lib.cc:1203
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
CylindricalManifold::get_new_point
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
Definition: manifold_lib.cc:1083
TransfiniteInterpolationManifold::get_new_point
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
Definition: manifold_lib.cc:2563
mapping.h
std::sqrt
inline ::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &x)
Definition: vectorization.h:5412
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
CylindricalManifold::push_forward
virtual Point< spacedim > push_forward(const Point< 3 > &chart_point) const override
Definition: manifold_lib.cc:1138
SymmetricTensor::determinant
constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &t)
Definition: symmetric_tensor.h:2645
DerivativeForm::covariant_form
DerivativeForm< 1, dim, spacedim, Number > covariant_form() const
VectorizedArray::sqrt
VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &x)
Definition: vectorization.h:5412
vector.h
StandardExceptions::ExcImpossibleInDim
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
PolarManifold::pull_back
virtual Point< spacedim > pull_back(const Point< spacedim > &space_point) const override
Definition: manifold_lib.cc:192
FunctionManifold::~FunctionManifold
virtual ~FunctionManifold() override
Definition: manifold_lib.cc:1430
CylindricalManifold::CylindricalManifold
CylindricalManifold(const unsigned int axis=0, const double tolerance=1e-10)
Definition: manifold_lib.cc:1038
SphericalManifold::get_normals_at_vertices
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Manifold< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
Definition: manifold_lib.cc:541
internal::apply_exponential_map
Tensor< 1, 3 > apply_exponential_map(const Tensor< 1, 3 > &u, const Tensor< 1, 3 > &dir)
Definition: manifold_lib.cc:46
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
PolarManifold::PolarManifold
PolarManifold(const Point< spacedim > center=Point< spacedim >())
Definition: manifold_lib.cc:122
SphericalManifold::get_intermediate_point
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
Definition: manifold_lib.cc:372
internal::projected_direction
Tensor< 1, 3 > projected_direction(const Tensor< 1, 3 > &u, const Tensor< 1, 3 > &v)
Definition: manifold_lib.cc:66
Point< spacedim >
FunctionManifold::pull_back
virtual Point< chartdim > pull_back(const Point< spacedim > &space_point) const override
Definition: manifold_lib.cc:1533
Function
Definition: function.h:151
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
internal
Definition: aligned_vector.h:369
TorusManifold::clone
virtual std::unique_ptr< Manifold< dim, 3 > > clone() const override
Definition: manifold_lib.cc:1597
memory.h
Differentiation::SD::acos
Expression acos(const Expression &x)
Definition: symengine_math.cc:140
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TransfiniteInterpolationManifold::push_forward_gradient
DerivativeForm< 1, dim, spacedim > push_forward_gradient(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &chart_point, const Point< spacedim > &pushed_forward_chart_point) const
Definition: manifold_lib.cc:2037
EllipticalManifold::push_forward
virtual Point< spacedim > push_forward(const Point< spacedim > &chart_point) const override
Definition: manifold_lib.cc:1255
first
Point< 2 > first
Definition: grid_out.cc:4352
numbers::PI
static constexpr double PI
Definition: numbers.h:237
GeometryInfo::distance_to_unit_cell
static double distance_to_unit_cell(const Point< dim > &p)
Vector< double >
PolarManifold::push_forward_gradient
virtual DerivativeForm< 1, spacedim, spacedim > push_forward_gradient(const Point< spacedim > &chart_point) const override
Definition: manifold_lib.cc:231
TorusManifold::pull_back
virtual Point< 3 > pull_back(const Point< 3 > &p) const override
Definition: manifold_lib.cc:1551
TriaIterator
Definition: tria_iterator.h:578
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
center
Point< 3 > center
Definition: data_out_base.cc:169
table.h
FunctionParser
Definition: function_parser.h:224
SphericalManifold::SphericalManifold
SphericalManifold(const Point< spacedim > center=Point< spacedim >())
Definition: manifold_lib.cc:353
Triangulation::get_manifold
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
Definition: tria.cc:10361
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
TorusManifold::TorusManifold
TorusManifold(const double R, const double r)
Definition: manifold_lib.cc:1582
FunctionManifold::push_forward
virtual Point< spacedim > push_forward(const Point< chartdim > &chart_point) const override
Definition: manifold_lib.cc:1488
SphericalManifold::guess_new_point
std::pair< double, Tensor< 1, spacedim > > guess_new_point(const ArrayView< const Tensor< 1, spacedim >> &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights) const
Definition: manifold_lib.cc:812