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