31#include <boost/container/small_vector.hpp>
54 const double theta = dir.
norm();
64 return tmp / tmp.
norm();
81 template <
int spacedim>
93 ExcMessage(
"The direction parameter must not be zero!"));
100 normal[0] = (vector[1] + vector[2]) / vector[0];
107 normal[1] = (vector[0] + vector[2]) / vector[1];
113 normal[2] = (vector[0] + vector[1]) / vector[2];
116 normal /= normal.
norm();
127template <
int dim,
int spacedim>
136template <
int dim,
int spacedim>
137std::unique_ptr<Manifold<dim, spacedim>>
140 return std::make_unique<PolarManifold<dim, spacedim>>(*this);
145template <
int dim,
int spacedim>
162template <
int dim,
int spacedim>
167 Assert(spherical_point[0] >= 0.0,
168 ExcMessage(
"Negative radius for given point."));
169 const double rho = spherical_point[0];
170 const double theta = spherical_point[1];
182 const double phi = spherical_point[2];
196template <
int dim,
int spacedim>
202 const double rho = R.
norm();
211 p[1] = std::atan2(R[1], R[0]);
219 const double z = R[2];
220 p[2] = std::atan2(R[1], R[0]);
223 p[1] = std::atan2(
std::sqrt(R[0] * R[0] + R[1] * R[1]), z);
235template <
int dim,
int spacedim>
240 Assert(spherical_point[0] >= 0.0,
241 ExcMessage(
"Negative radius for given point."));
242 const double rho = spherical_point[0];
243 const double theta = spherical_point[1];
260 const double phi = spherical_point[2];
285 template <
int dim,
int spacedim>
287 spherical_face_is_horizontal(
297 constexpr unsigned int n_vertices =
299 std::array<double, n_vertices> sqr_distances_to_center;
300 std::array<double, n_vertices - 1> sqr_distances_to_first_vertex;
301 sqr_distances_to_center[0] =
302 (face->vertex(0) - manifold_center).norm_square();
303 for (
unsigned int i = 1; i < n_vertices; ++i)
305 sqr_distances_to_center[i] =
306 (face->vertex(i) - manifold_center).norm_square();
307 sqr_distances_to_first_vertex[i - 1] =
308 (face->vertex(i) - face->vertex(0)).norm_square();
310 const auto minmax_sqr_distance =
311 std::minmax_element(sqr_distances_to_center.begin(),
312 sqr_distances_to_center.end());
313 const auto min_sqr_distance_to_first_vertex =
314 std::min_element(sqr_distances_to_first_vertex.begin(),
315 sqr_distances_to_first_vertex.end());
317 return (*minmax_sqr_distance.second - *minmax_sqr_distance.first <
318 1.e-10 * *min_sqr_distance_to_first_vertex);
324template <
int dim,
int spacedim>
334 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
341 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
342 return normalized_spherical_normal;
358template <
int dim,
int spacedim>
367template <
int dim,
int spacedim>
368std::unique_ptr<Manifold<dim, spacedim>>
371 return std::make_unique<SphericalManifold<dim, spacedim>>(*this);
376template <
int dim,
int spacedim>
381 const double w)
const
383 const double tol = 1e-10;
385 if ((p1 - p2).norm_square() < tol * tol ||
std::abs(w) < tol)
397 const double r1 =
v1.norm();
398 const double r2 = v2.
norm();
400 Assert(r1 > tol && r2 > tol,
401 ExcMessage(
"p1 and p2 cannot coincide with the center."));
407 const double cosgamma = e1 * e2;
410 if (cosgamma < -1 + 8. * std::numeric_limits<double>::epsilon())
414 if (cosgamma > 1 - 8. * std::numeric_limits<double>::epsilon())
420 const double sigma = w * std::acos(cosgamma);
425 const double n_norm = n.
norm();
428 "Probably, this means v1==v2 or v2==0."));
442template <
int dim,
int spacedim>
448 const double tol = 1e-10;
455 const double r1 =
v1.norm();
456 const double r2 = v2.
norm();
466 const double cosgamma = e1 * e2;
468 Assert(cosgamma > -1 + 8. * std::numeric_limits<double>::epsilon(),
469 ExcMessage(
"p1 and p2 cannot lie on the same diameter and be opposite "
470 "respect to the center."));
472 if (cosgamma > 1 - 8. * std::numeric_limits<double>::epsilon())
478 const double n_norm = n.
norm();
481 "Probably, this means v1==v2 or v2==0."));
487 const double gamma = std::acos(cosgamma);
488 return (r2 - r1) * e1 + r1 * gamma * n;
493template <
int dim,
int spacedim>
503 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
510 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
511 return normalized_spherical_normal;
545template <
int dim,
int spacedim>
556 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
561 for (
unsigned int vertex = 0;
562 vertex < GeometryInfo<spacedim>::vertices_per_face;
564 face_vertex_normals[vertex] = face->vertex(vertex) -
center;
572template <
int dim,
int spacedim>
582 get_new_points(surrounding_points,
make_array_view(weights), new_points);
589template <
int dim,
int spacedim>
607template <
int dim,
int spacedim>
615 new_points.size() * surrounding_points.size());
616 const unsigned int weight_rows = new_points.size();
617 const unsigned int weight_columns = surrounding_points.size();
619 if (surrounding_points.size() == 2)
621 for (
unsigned int row = 0; row < weight_rows; ++row)
623 get_intermediate_point(surrounding_points[0],
624 surrounding_points[1],
625 weights[row * weight_columns + 1]);
629 boost::container::small_vector<std::pair<double, Tensor<1, spacedim>>, 100>
630 new_candidates(new_points.size());
631 boost::container::small_vector<Tensor<1, spacedim>, 100> directions(
633 boost::container::small_vector<double, 100> distances(
634 surrounding_points.size(), 0.0);
635 double max_distance = 0.;
636 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
638 directions[i] = surrounding_points[i] -
center;
639 distances[i] = directions[i].
norm();
641 if (distances[i] != 0.)
642 directions[i] /= distances[i];
645 ExcMessage(
"One of the vertices coincides with the center. "
646 "This is not allowed!"));
650 for (
unsigned int k = 0; k < i; ++k)
652 const double squared_distance =
653 (directions[i] - directions[k]).norm_square();
654 max_distance =
std::max(max_distance, squared_distance);
660 const double tolerance = 1e-10;
661 boost::container::small_vector<bool, 100> accurate_point_was_found(
662 new_points.size(),
false);
667 for (
unsigned int row = 0; row < weight_rows; ++row)
669 new_candidates[row] =
670 guess_new_point(array_directions,
677 if (new_candidates[row].
first == 0.0)
680 accurate_point_was_found[row] =
true;
687 new_points[row] = polar_manifold.get_new_point(
702 if (max_distance < 2e-2)
704 for (
unsigned int row = 0; row < weight_rows; ++row)
716 boost::container::small_vector<double, 1000> merged_weights(weights.
size());
717 boost::container::small_vector<Tensor<1, spacedim>, 100> merged_directions(
719 boost::container::small_vector<double, 100> merged_distances(
720 surrounding_points.size(), 0.0);
722 unsigned int n_unique_directions = 0;
723 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
725 bool found_duplicate =
false;
729 for (
unsigned int j = 0; j < n_unique_directions; ++j)
731 const double squared_distance =
732 (directions[i] - directions[j]).norm_square();
733 if (!found_duplicate && squared_distance < 1e-28)
735 found_duplicate =
true;
736 for (
unsigned int row = 0; row < weight_rows; ++row)
737 merged_weights[row * weight_columns + j] +=
738 weights[row * weight_columns + i];
742 if (found_duplicate ==
false)
744 merged_directions[n_unique_directions] = directions[i];
745 merged_distances[n_unique_directions] = distances[i];
746 for (
unsigned int row = 0; row < weight_rows; ++row)
747 merged_weights[row * weight_columns + n_unique_directions] =
748 weights[row * weight_columns + i];
750 ++n_unique_directions;
756 boost::container::small_vector<unsigned int, 100> merged_weights_index(
758 for (
unsigned int row = 0; row < weight_rows; ++row)
760 for (
unsigned int existing_row = 0; existing_row < row; ++existing_row)
762 bool identical_weights =
true;
764 for (
unsigned int weight_index = 0;
765 weight_index < n_unique_directions;
767 if (
std::abs(merged_weights[row * weight_columns + weight_index] -
768 merged_weights[existing_row * weight_columns +
769 weight_index]) > tolerance)
771 identical_weights =
false;
775 if (identical_weights)
777 merged_weights_index[row] = existing_row;
787 merged_directions.begin() + n_unique_directions);
790 merged_distances.begin() + n_unique_directions);
792 for (
unsigned int row = 0; row < weight_rows; ++row)
793 if (!accurate_point_was_found[row])
798 &merged_weights[row * weight_columns], n_unique_directions);
799 new_candidates[row].second =
800 get_new_point(array_merged_directions,
801 array_merged_distances,
802 array_merged_weights,
806 new_candidates[row].second =
807 new_candidates[merged_weights_index[row]].second;
810 center + new_candidates[row].first * new_candidates[row].second;
816template <
int dim,
int spacedim>
817std::pair<double, Tensor<1, spacedim>>
823 const double tolerance = 1e-10;
828 double total_weights = 0.;
829 for (
unsigned int i = 0; i < directions.size(); ++i)
832 if (
std::abs(1 - weights[i]) < tolerance)
833 return std::make_pair(distances[i], directions[i]);
835 rho += distances[i] * weights[i];
836 candidate += directions[i] * weights[i];
837 total_weights += weights[i];
841 const double norm = candidate.
norm();
845 rho /= total_weights;
847 return std::make_pair(rho, candidate);
853 template <
int spacedim>
876 Point<3> candidate = candidate_point;
877 const unsigned int n_merged_points = directions.size();
878 const double tolerance = 1
e-10;
879 const int max_iterations = 10;
884 for (
unsigned int i = 0; i < n_merged_points; ++i)
886 const double squared_distance =
887 (candidate - directions[i]).norm_square();
888 if (squared_distance < tolerance * tolerance)
894 if (n_merged_points == 2)
913 for (
unsigned int i = 0; i < max_iterations; ++i)
919 const Tensor<1, 3> Clocaly = cross_product_3d(candidate, Clocalx);
927 for (
unsigned int i = 0; i < n_merged_points; ++i)
932 const double sintheta =
std::sqrt(sinthetaSq);
933 if (sintheta < tolerance)
935 Hessian[0][0] += weights[i];
936 Hessian[1][1] += weights[i];
940 const double costheta = (directions[i]) * candidate;
941 const double theta = std::atan2(sintheta, costheta);
942 const double sincthetaInv = theta / sintheta;
944 const double cosphi = vPerp * Clocalx;
945 const double sinphi = vPerp * Clocaly;
947 gradlocal[0] = cosphi;
948 gradlocal[1] = sinphi;
949 gradient += (weights[i] * sincthetaInv) * gradlocal;
951 const double wt = weights[i] / sinthetaSq;
952 const double sinphiSq = sinphi * sinphi;
953 const double cosphiSq = cosphi * cosphi;
954 const double tt = sincthetaInv * costheta;
955 const double offdiag = cosphi * sinphi * wt * (1.0 - tt);
956 Hessian[0][0] += wt * (cosphiSq + tt * sinphiSq);
957 Hessian[0][1] += offdiag;
958 Hessian[1][0] += offdiag;
959 Hessian[1][1] += wt * (sinphiSq + tt * cosphiSq);
969 xDisplocal[0] * Clocalx + xDisplocal[1] * Clocaly;
972 const Point<3> candidateOld = candidate;
977 if ((candidate - candidateOld).norm_square() < tolerance * tolerance)
987template <
int dim,
int spacedim>
1007 const Point<3> & candidate_point)
const
1009 return do_get_new_point(directions, distances, weights, candidate_point);
1020 const Point<3> & candidate_point)
const
1022 return do_get_new_point(directions, distances, weights, candidate_point);
1033 const Point<3> & candidate_point)
const
1035 return do_get_new_point(directions, distances, weights, candidate_point);
1043template <
int dim,
int spacedim>
1045 const double tolerance)
1053 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1058template <
int dim,
int spacedim>
1062 const double tolerance)
1065 , direction(direction / direction.norm())
1066 , point_on_axis(point_on_axis)
1067 , tolerance(tolerance)
1068 , dxn(cross_product_3d(this->direction, normal_direction))
1073 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1078template <
int dim,
int spacedim>
1079std::unique_ptr<Manifold<dim, spacedim>>
1082 return std::make_unique<CylindricalManifold<dim, spacedim>>(*this);
1087template <
int dim,
int spacedim>
1094 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1098 double average_length = 0.;
1099 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
1101 middle += surrounding_points[i] * weights[i];
1102 average_length += surrounding_points[i].square() * weights[i];
1104 middle -= point_on_axis;
1105 const double lambda = middle * direction;
1107 if ((middle - direction * lambda).square() < tolerance * average_length)
1108 return point_on_axis + direction * lambda;
1116template <
int dim,
int spacedim>
1122 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1126 const double lambda = normalized_point * direction;
1127 const Point<spacedim> projection = point_on_axis + direction * lambda;
1142template <
int dim,
int spacedim>
1148 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1151 const double sine_r =
std::sin(chart_point(1)) * chart_point(0);
1152 const double cosine_r =
std::cos(chart_point(1)) * chart_point(0);
1154 normal_direction * cosine_r + dxn * sine_r;
1157 return point_on_axis + direction * chart_point(2) + intermediate;
1162template <
int dim,
int spacedim>
1168 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1173 const double sine =
std::sin(chart_point(1));
1174 const double cosine =
std::cos(chart_point(1));
1176 normal_direction * cosine + dxn * sine;
1179 constexpr int s0 = 0 % spacedim;
1180 constexpr int s1 = 1 % spacedim;
1181 constexpr int s2 = 2 % spacedim;
1184 derivatives[s0][s0] = intermediate[s0];
1185 derivatives[s1][s0] = intermediate[s1];
1186 derivatives[s2][s0] = intermediate[s2];
1189 derivatives[s0][s1] = -normal_direction[s0] * sine + dxn[s0] * cosine;
1190 derivatives[s1][s1] = -normal_direction[s1] * sine + dxn[s1] * cosine;
1191 derivatives[s2][s1] = -normal_direction[s2] * sine + dxn[s2] * cosine;
1194 derivatives[s0][s2] = direction[s0];
1195 derivatives[s1][s2] = direction[s1];
1196 derivatives[s2][s2] = direction[s2];
1206template <
int dim,
int spacedim>
1210 const double eccentricity)
1213 , direction(major_axis_direction)
1215 , cosh_u(1.0 / eccentricity)
1216 , sinh_u(
std::sqrt(cosh_u * cosh_u - 1.0))
1223 "Invalid eccentricity: It must satisfy 0 < eccentricity < 1."));
1225 Assert(direction_norm != 0.0,
1227 "Invalid major axis direction vector: Null vector not allowed."));
1233template <
int dim,
int spacedim>
1234std::unique_ptr<Manifold<dim, spacedim>>
1237 return std::make_unique<EllipticalManifold<dim, spacedim>>(*this);
1242template <
int dim,
int spacedim>
1255template <
int dim,
int spacedim>
1269 const double cs =
std::cos(chart_point[1]);
1270 const double sn =
std::sin(chart_point[1]);
1273 const double x = chart_point[0] * cosh_u * cs;
1274 const double y = chart_point[0] * sinh_u * sn;
1276 const Point<2> p(direction[0] * x - direction[1] * y,
1277 direction[1] * x + direction[0] * y);
1283template <
int dim,
int spacedim>
1298 const double x0 = space_point[0] -
center[0];
1299 const double y0 = space_point[1] -
center[1];
1300 const double x = direction[0] * x0 + direction[1] * y0;
1301 const double y = -direction[1] * x0 + direction[0] * y0;
1303 std::sqrt((x * x) / (cosh_u * cosh_u) + (y * y) / (sinh_u * sinh_u));
1309 double cos_eta = x / (pt0 * cosh_u);
1318 const double eta = std::acos(cos_eta);
1319 const double pt1 = (std::signbit(y) ? 2.0 *
numbers::PI - eta : eta);
1325template <
int dim,
int spacedim>
1341 const double cs =
std::cos(chart_point[1]);
1342 const double sn =
std::sin(chart_point[1]);
1344 dX[0][0] = cosh_u * cs;
1345 dX[0][1] = -chart_point[0] * cosh_u * sn;
1346 dX[1][0] = sinh_u * sn;
1347 dX[1][1] = chart_point[0] * sinh_u * cs;
1351 {{+direction[0], -direction[1]}, {direction[1], direction[0]}}};
1361template <
int dim,
int spacedim,
int chartdim>
1366 const double tolerance)
1369 , push_forward_function(&push_forward_function)
1370 , pull_back_function(&pull_back_function)
1371 , tolerance(tolerance)
1372 , owns_pointers(false)
1373 , finite_difference_step(0)
1381template <
int dim,
int spacedim,
int chartdim>
1386 const double tolerance)
1389 , push_forward_function(push_forward.release())
1390 , pull_back_function(pull_back.release())
1391 , tolerance(tolerance)
1392 , owns_pointers(true)
1393 , finite_difference_step(0)
1401template <
int dim,
int spacedim,
int chartdim>
1403 const std::string push_forward_expression,
1404 const std::string pull_back_expression,
1407 const std::string chart_vars,
1408 const std::string space_vars,
1409 const double tolerance,
1412 , const_map(const_map)
1413 , tolerance(tolerance)
1414 , owns_pointers(true)
1415 , push_forward_expression(push_forward_expression)
1416 , pull_back_expression(pull_back_expression)
1417 , chart_vars(chart_vars)
1418 , space_vars(space_vars)
1419 , finite_difference_step(h)
1431template <
int dim,
int spacedim,
int chartdim>
1434 if (owns_pointers ==
true)
1437 push_forward_function =
nullptr;
1441 pull_back_function =
nullptr;
1448template <
int dim,
int spacedim,
int chartdim>
1449std::unique_ptr<Manifold<dim, spacedim>>
1464 if (!(push_forward_expression.empty() && pull_back_expression.empty()))
1466 return std::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1467 push_forward_expression,
1468 pull_back_expression,
1469 this->get_periodicity(),
1474 finite_difference_step);
1478 return std::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1479 *push_forward_function,
1480 *pull_back_function,
1481 this->get_periodicity(),
1488template <
int dim,
int spacedim,
int chartdim>
1495 push_forward_function->vector_value(chart_point, pf);
1496 for (
unsigned int i = 0; i < spacedim; ++i)
1501 pull_back_function->vector_value(result, pb);
1502 for (
unsigned int i = 0; i < chartdim; ++i)
1504 (chart_point.
norm() > tolerance &&
1505 (
std::abs(pb[i] - chart_point[i]) < tolerance * chart_point.
norm())) ||
1506 (
std::abs(pb[i] - chart_point[i]) < tolerance),
1508 "The push forward is not the inverse of the pull back! Bailing out."));
1516template <
int dim,
int spacedim,
int chartdim>
1522 for (
unsigned int i = 0; i < spacedim; ++i)
1524 const auto gradient = push_forward_function->gradient(chart_point, i);
1525 for (
unsigned int j = 0; j < chartdim; ++j)
1526 DF[i][j] = gradient[j];
1533template <
int dim,
int spacedim,
int chartdim>
1540 pull_back_function->vector_value(space_point, pb);
1541 for (
unsigned int i = 0; i < chartdim; ++i)
1558 double phi = std::atan2(y, x);
1559 double theta = std::atan2(z,
std::sqrt(x * x + y * y) - R);
1562 Utilities::fixed_power<2>(x -
std::cos(phi) * R) + z * z) /
1564 return {phi, theta, w};
1573 double phi = chart_point(0);
1574 double theta = chart_point(1);
1575 double w = chart_point(2);
1591 ExcMessage(
"Outer radius R must be greater than the inner "
1599std::unique_ptr<Manifold<dim, 3>>
1602 return std::make_unique<TorusManifold<dim>>(R, r);
1613 double phi = chart_point(0);
1614 double theta = chart_point(1);
1615 double w = chart_point(2);
1622 DX[1][1] = r * w *
std::cos(theta);
1637template <
int dim,
int spacedim>
1648template <
int dim,
int spacedim>
1650 spacedim>::~TransfiniteInterpolationManifold()
1652 if (clear_signal.connected())
1653 clear_signal.disconnect();
1658template <
int dim,
int spacedim>
1659std::unique_ptr<Manifold<dim, spacedim>>
1665 return std::unique_ptr<Manifold<dim, spacedim>>(ptr);
1670template <
int dim,
int spacedim>
1677 clear_signal.disconnect();
1678 clear_signal =
triangulation.signals.clear.connect([&]() ->
void {
1679 this->triangulation =
nullptr;
1680 this->level_coarse = -1;
1683 coarse_cell_is_flat.resize(
triangulation.n_cells(level_coarse),
false);
1684 quadratic_approximation.clear();
1694 std::vector<Point<dim>> unit_points =
1696 std::vector<Point<spacedim>> real_points(unit_points.size());
1698 for (
const auto &cell :
triangulation.active_cell_iterators())
1700 bool cell_is_flat =
true;
1701 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
1702 if (cell->line(l)->manifold_id() != cell->manifold_id() &&
1704 cell_is_flat =
false;
1706 for (
unsigned int q = 0; q < GeometryInfo<dim>::quads_per_cell; ++q)
1707 if (cell->quad(q)->manifold_id() != cell->manifold_id() &&
1709 cell_is_flat =
false;
1711 coarse_cell_is_flat.size());
1712 coarse_cell_is_flat[cell->index()] = cell_is_flat;
1715 for (
unsigned int i = 0; i < unit_points.size(); ++i)
1716 real_points[i] = push_forward(cell, unit_points[i]);
1717 quadratic_approximation.emplace_back(real_points, unit_points);
1726 template <
typename AccessorType>
1728 compute_transfinite_interpolation(
const AccessorType &cell,
1732 return cell.vertex(0) * (1. - chart_point[0]) +
1733 cell.vertex(1) * chart_point[0];
1737 template <
typename AccessorType>
1739 compute_transfinite_interpolation(
const AccessorType &cell,
1741 const bool cell_is_flat)
1743 const unsigned int dim = AccessorType::dimension;
1744 const unsigned int spacedim = AccessorType::space_dimension;
1752 const std::array<Point<spacedim>, 4>
vertices{
1753 {cell.vertex(0), cell.vertex(1), cell.vertex(2), cell.vertex(3)}};
1758 std::array<double, 4> weights_vertices{
1759 {(1. - chart_point[0]) * (1. - chart_point[1]),
1760 chart_point[0] * (1. - chart_point[1]),
1761 (1. - chart_point[0]) * chart_point[1],
1762 chart_point[0] * chart_point[1]}};
1767 new_point += weights_vertices[v] *
vertices[v];
1779 std::array<double, GeometryInfo<2>::vertices_per_face> weights;
1782 const auto weights_view =
1784 const auto points_view =
make_array_view(points.begin(), points.end());
1786 for (
unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1789 const double my_weight =
1790 (line % 2) ? chart_point[line / 2] : 1 - chart_point[line / 2];
1791 const double line_point = chart_point[1 - line / 2];
1798 cell.line(line)->manifold_id();
1799 if (line_manifold_id == my_manifold_id ||
1804 my_weight * (1. - line_point);
1807 my_weight * line_point;
1815 weights[0] = 1. - line_point;
1816 weights[1] = line_point;
1825 new_point -= weights_vertices[v] *
vertices[v];
1834 static constexpr unsigned int face_to_cell_vertices_3d[6][4] = {{0, 2, 4, 6},
1844 static constexpr unsigned int face_to_cell_lines_3d[6][4] = {{8, 10, 0, 4},
1852 template <
typename AccessorType>
1854 compute_transfinite_interpolation(
const AccessorType &cell,
1856 const bool cell_is_flat)
1858 const unsigned int dim = AccessorType::dimension;
1859 const unsigned int spacedim = AccessorType::space_dimension;
1865 const std::array<Point<spacedim>, 8>
vertices{{cell.vertex(0),
1877 double linear_shapes[10];
1878 for (
unsigned int d = 0;
d < 3; ++
d)
1880 linear_shapes[2 *
d] = 1. - chart_point[
d];
1881 linear_shapes[2 *
d + 1] = chart_point[
d];
1885 for (
unsigned int d = 6;
d < 10; ++
d)
1886 linear_shapes[d] = linear_shapes[d - 6];
1888 std::array<double, 8> weights_vertices;
1889 for (
unsigned int i2 = 0, v = 0; i2 < 2; ++i2)
1890 for (
unsigned int i1 = 0; i1 < 2; ++i1)
1891 for (
unsigned int i0 = 0; i0 < 2; ++i0, ++v)
1892 weights_vertices[v] =
1893 (linear_shapes[4 + i2] * linear_shapes[2 + i1]) * linear_shapes[i0];
1897 for (
unsigned int v = 0; v < 8; ++v)
1898 new_point += weights_vertices[v] *
vertices[v];
1904 std::array<double, GeometryInfo<3>::lines_per_cell> weights_lines;
1905 std::fill(weights_lines.begin(), weights_lines.end(), 0.0);
1908 std::array<double, GeometryInfo<2>::vertices_per_cell> weights;
1911 const auto weights_view =
1913 const auto points_view =
make_array_view(points.begin(), points.end());
1915 for (
const unsigned int face :
GeometryInfo<3>::face_indices())
1917 const double my_weight = linear_shapes[face];
1918 const unsigned int face_even = face - face % 2;
1926 cell.face(face)->manifold_id();
1927 if (face_manifold_id == my_manifold_id ||
1930 for (
unsigned int line = 0;
1931 line < GeometryInfo<2>::lines_per_cell;
1934 const double line_weight =
1935 linear_shapes[face_even + 2 + line];
1936 weights_lines[face_to_cell_lines_3d[face][line]] +=
1937 my_weight * line_weight;
1943 weights_vertices[face_to_cell_vertices_3d[face][0]] -=
1944 linear_shapes[face_even + 2] *
1945 (linear_shapes[face_even + 4] * my_weight);
1946 weights_vertices[face_to_cell_vertices_3d[face][1]] -=
1947 linear_shapes[face_even + 3] *
1948 (linear_shapes[face_even + 4] * my_weight);
1949 weights_vertices[face_to_cell_vertices_3d[face][2]] -=
1950 linear_shapes[face_even + 2] *
1951 (linear_shapes[face_even + 5] * my_weight);
1952 weights_vertices[face_to_cell_vertices_3d[face][3]] -=
1953 linear_shapes[face_even + 3] *
1954 (linear_shapes[face_even + 5] * my_weight);
1959 points[v] =
vertices[face_to_cell_vertices_3d[face][v]];
1961 linear_shapes[face_even + 2] * linear_shapes[face_even + 4];
1963 linear_shapes[face_even + 3] * linear_shapes[face_even + 4];
1965 linear_shapes[face_even + 2] * linear_shapes[face_even + 5];
1967 linear_shapes[face_even + 3] * linear_shapes[face_even + 5];
1975 const auto weights_view_line =
1977 const auto points_view_line =
1979 for (
unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
1982 const double line_point =
1983 (line < 8 ? chart_point[1 - (line % 4) / 2] : chart_point[2]);
1984 double my_weight = 0.;
1986 my_weight = linear_shapes[line % 4] * linear_shapes[4 + line / 4];
1989 const unsigned int subline = line - 8;
1991 linear_shapes[subline % 2] * linear_shapes[2 + subline / 2];
1993 my_weight -= weights_lines[line];
1999 cell.line(line)->manifold_id();
2000 if (line_manifold_id == my_manifold_id ||
2005 my_weight * (1. - line_point);
2008 my_weight * (line_point);
2016 weights[0] = 1. - line_point;
2017 weights[1] = line_point;
2026 new_point += weights_vertices[v] *
vertices[v];
2034template <
int dim,
int spacedim>
2046 ExcMessage(
"chart_point is not in unit interval"));
2048 return compute_transfinite_interpolation(*cell,
2050 coarse_cell_is_flat[cell->index()]);
2055template <
int dim,
int spacedim>
2064 for (
unsigned int d = 0; d < dim; ++d)
2067 const double step = chart_point[d] > 0.5 ? -1e-8 : 1e-8;
2070 modified[d] += step;
2072 compute_transfinite_interpolation(*cell,
2074 coarse_cell_is_flat[cell->index()]) -
2075 pushed_forward_chart_point;
2076 for (
unsigned int e = 0; e < spacedim; ++e)
2077 grad[e][d] = difference[e] / step;
2084template <
int dim,
int spacedim>
2092 for (
unsigned int d = 0; d < dim; ++d)
2106 compute_transfinite_interpolation(*cell,
2108 coarse_cell_is_flat[cell->index()]);
2109 const double tolerance = 1e-21 * Utilities::fixed_power<2>(cell->diameter());
2110 double residual_norm_square = residual.
norm_square();
2112 bool must_recompute_jacobian =
true;
2113 for (
unsigned int i = 0; i < 100; ++i)
2115 if (residual_norm_square < tolerance)
2122 for (
unsigned int d = 0; d < spacedim; ++d)
2123 for (
unsigned int e = 0; e < dim; ++e)
2124 update[e] += inv_grad[d][e] * residual[d];
2125 return chart_point + update;
2137 if (must_recompute_jacobian || i % 9 == 0)
2145 push_forward_gradient(cell,
2151 must_recompute_jacobian =
false;
2154 for (
unsigned int d = 0; d < spacedim; ++d)
2155 for (
unsigned int e = 0; e < dim; ++e)
2156 update[e] += inv_grad[d][e] * residual[d];
2170 while (alpha > 1e-4)
2172 Point<dim> guess = chart_point + alpha * update;
2174 point - compute_transfinite_interpolation(
2175 *cell, guess, coarse_cell_is_flat[cell->index()]);
2176 const double residual_norm_new = residual_guess.
norm_square();
2177 if (residual_norm_new < residual_norm_square)
2179 residual = residual_guess;
2180 residual_norm_square = residual_norm_new;
2181 chart_point += alpha * update;
2194 must_recompute_jacobian =
true;
2208 for (
unsigned int d = 0; d < spacedim; ++d)
2209 for (
unsigned int e = 0; e < dim; ++e)
2210 Jinv_deltaf[e] += inv_grad[d][e] * delta_f[d];
2217 if (
std::abs(delta_x * Jinv_deltaf) > 1e-12 && !must_recompute_jacobian)
2220 (delta_x - Jinv_deltaf) / (delta_x * Jinv_deltaf);
2222 for (
unsigned int d = 0; d < spacedim; ++d)
2223 for (
unsigned int e = 0; e < dim; ++e)
2224 jac_update[d] += delta_x[e] * inv_grad[d][e];
2225 for (
unsigned int d = 0; d < spacedim; ++d)
2226 for (
unsigned int e = 0; e < dim; ++e)
2227 inv_grad[d][e] += factor[e] * jac_update[d];
2235template <
int dim,
int spacedim>
2236std::array<unsigned int, 20>
2246 ExcMessage(
"The manifold was initialized with level " +
2247 std::to_string(level_coarse) +
" but there are now" +
2248 "active cells on a lower level. Coarsening the mesh is " +
2249 "currently not supported"));
2259 boost::container::small_vector<std::pair<double, unsigned int>, 200>
2260 distances_and_cells;
2261 for (; cell != endc; ++cell)
2264 if (&cell->get_manifold() !=
this)
2271 vertices[vertex_n] = cell->vertex(vertex_n);
2281 double radius_square = 0.;
2285 bool inside_circle =
true;
2286 for (
unsigned int i = 0; i < points.size(); ++i)
2287 if ((
center - points[i]).norm_square() > radius_square * 1.5)
2289 inside_circle =
false;
2292 if (inside_circle ==
false)
2296 double current_distance = 0;
2297 for (
unsigned int i = 0; i < points.size(); ++i)
2300 quadratic_approximation[cell->index()].compute(points[i]);
2303 distances_and_cells.push_back(
2304 std::make_pair(current_distance, cell->index()));
2309 std::sort(distances_and_cells.begin(), distances_and_cells.end());
2310 std::array<unsigned int, 20> cells;
2312 for (
unsigned int i = 0; i < distances_and_cells.size() && i < cells.size();
2314 cells[i] = distances_and_cells[i].
second;
2321template <
int dim,
int spacedim>
2327 Assert(surrounding_points.size() == chart_points.size(),
2328 ExcMessage(
"The chart points array view must be as large as the "
2329 "surrounding points array view."));
2331 std::array<unsigned int, 20> nearby_cells =
2332 get_possible_cells_around_points(surrounding_points);
2356 auto guess_chart_point_structdim_2 = [&](
const unsigned int i) ->
Point<dim> {
2357 Assert(surrounding_points.size() == 8 && 2 < i && i < 8,
2358 ExcMessage(
"This function assumes that there are eight surrounding "
2359 "points around a two-dimensional object. It also assumes "
2360 "that the first three chart points have already been "
2370 return chart_points[1] + (chart_points[2] - chart_points[0]);
2372 return 0.5 * (chart_points[0] + chart_points[2]);
2374 return 0.5 * (chart_points[1] + chart_points[3]);
2376 return 0.5 * (chart_points[0] + chart_points[1]);
2378 return 0.5 * (chart_points[2] + chart_points[3]);
2407 auto guess_chart_point_structdim_3 = [&](
const unsigned int i) ->
Point<dim> {
2408 Assert(surrounding_points.size() == 8 && 4 < i && i < 8,
2409 ExcMessage(
"This function assumes that there are eight surrounding "
2410 "points around a three-dimensional object. It also "
2411 "assumes that the first five chart points have already "
2413 return chart_points[i - 4] + (chart_points[4] - chart_points[0]);
2417 bool use_structdim_2_guesses =
false;
2418 bool use_structdim_3_guesses =
false;
2423 if (surrounding_points.size() == 8)
2426 surrounding_points[6] - surrounding_points[0];
2428 surrounding_points[7] - surrounding_points[2];
2431 const double cosine = scalar_product(v06, v27) /
2436 use_structdim_2_guesses =
true;
2437 else if (spacedim == 3)
2440 use_structdim_3_guesses =
true;
2443 Assert((!use_structdim_2_guesses && !use_structdim_3_guesses) ||
2444 (use_structdim_2_guesses ^ use_structdim_3_guesses),
2449 auto compute_chart_point =
2451 const unsigned int point_index) {
2457 bool used_quadratic_approximation =
false;
2460 if (point_index == 3 && surrounding_points.size() >= 8)
2461 guess = chart_points[1] + (chart_points[2] - chart_points[0]);
2462 else if (use_structdim_2_guesses && 3 < point_index)
2463 guess = guess_chart_point_structdim_2(point_index);
2464 else if (use_structdim_3_guesses && 4 < point_index)
2465 guess = guess_chart_point_structdim_3(point_index);
2466 else if (dim == 3 && point_index > 7 && surrounding_points.size() == 26)
2468 if (point_index < 20)
2471 point_index - 8, 0)] +
2473 point_index - 8, 1)]);
2477 point_index - 20, 0)] +
2479 point_index - 20, 1)] +
2481 point_index - 20, 2)] +
2483 point_index - 20, 3)]);
2487 guess = quadratic_approximation[cell->index()].compute(
2488 surrounding_points[point_index]);
2489 used_quadratic_approximation =
true;
2491 chart_points[point_index] =
2492 pull_back(cell, surrounding_points[point_index], guess);
2497 if (chart_points[point_index][0] ==
2499 !used_quadratic_approximation)
2501 guess = quadratic_approximation[cell->index()].compute(
2502 surrounding_points[point_index]);
2503 chart_points[point_index] =
2504 pull_back(cell, surrounding_points[point_index], guess);
2507 if (chart_points[point_index][0] ==
2510 for (
unsigned int d = 0; d < dim; ++d)
2512 chart_points[point_index] =
2513 pull_back(cell, surrounding_points[point_index], guess);
2518 for (
unsigned int c = 0; c < nearby_cells.size(); ++c)
2522 bool inside_unit_cell =
true;
2523 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2525 compute_chart_point(cell, i);
2532 inside_unit_cell =
false;
2536 if (inside_unit_cell ==
true)
2543 if (c == nearby_cells.size() - 1 ||
2548 std::ostringstream message;
2549 for (
unsigned int b = 0; b <= c; ++b)
2553 message <<
"Looking at cell " << cell->id()
2554 <<
" with vertices: " << std::endl;
2556 message << cell->vertex(v) <<
" ";
2557 message << std::endl;
2558 message <<
"Transformation to chart coordinates: " << std::endl;
2559 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2561 compute_chart_point(cell, i);
2562 message << surrounding_points[i] <<
" -> " << chart_points[i]
2582template <
int dim,
int spacedim>
2588 boost::container::small_vector<Point<dim>, 100> chart_points(
2589 surrounding_points.size());
2592 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2595 chart_manifold.get_new_point(chart_points_view, weights);
2597 return push_forward(cell, p_chart);
2602template <
int dim,
int spacedim>
2612 boost::container::small_vector<Point<dim>, 100> chart_points(
2613 surrounding_points.size());
2616 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2618 boost::container::small_vector<Point<dim>, 100> new_points_on_chart(
2620 chart_manifold.get_new_points(chart_points_view,
2623 new_points_on_chart.end()));
2625 for (
unsigned int row = 0; row < weights.size(0); ++row)
2626 new_points[row] = push_forward(cell, new_points_on_chart[row]);
2632#include "manifold_lib.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
virtual Point< spacedim > push_forward(const Point< 3 > &chart_point) const override
virtual DerivativeForm< 1, 3, spacedim > push_forward_gradient(const Point< 3 > &chart_point) const override
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
CylindricalManifold(const unsigned int axis=0, const double tolerance=1e-10)
virtual Point< 3 > pull_back(const Point< spacedim > &space_point) const override
Tensor< 1, spacedim > direction
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< spacedim > push_forward(const Point< spacedim > &chart_point) const override
EllipticalManifold(const Point< spacedim > ¢er, const Tensor< 1, spacedim > &major_axis_direction, const double eccentricity)
virtual DerivativeForm< 1, spacedim, spacedim > push_forward_gradient(const Point< spacedim > &chart_point) const override
virtual Point< spacedim > pull_back(const Point< spacedim > &space_point) const override
static Tensor< 1, spacedim > get_periodicity()
SmartPointer< const Function< spacedim >, FunctionManifold< dim, spacedim, chartdim > > pull_back_function
const FunctionParser< spacedim >::ConstMap const_map
const std::string chart_vars
const std::string pull_back_expression
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const override
SmartPointer< const Function< chartdim >, FunctionManifold< dim, spacedim, chartdim > > push_forward_function
const std::string space_vars
virtual Point< spacedim > push_forward(const Point< chartdim > &chart_point) const override
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)
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
const std::string push_forward_expression
virtual Point< chartdim > pull_back(const Point< spacedim > &space_point) const override
virtual ~FunctionManifold() override
virtual void initialize(const std::string &vars, const std::vector< std::string > &expressions, const ConstMap &constants, const bool time_dependent=false) override
std::map< std::string, double > ConstMap
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face > FaceVertexNormals
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const
Abstract base class for mapping classes.
PolarManifold(const Point< spacedim > center=Point< spacedim >())
static Tensor< 1, spacedim > get_periodicity()
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
virtual Point< spacedim > push_forward(const Point< spacedim > &chart_point) const override
virtual Point< spacedim > pull_back(const Point< spacedim > &space_point) const override
const std::vector< Point< dim > > & get_points() const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &vertices, const ArrayView< const double > &weights) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, typename Manifold< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override
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
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const override
SphericalManifold(const Point< spacedim > center=Point< spacedim >())
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
numbers::NumberTraits< Number >::real_type norm() const
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
virtual Point< 3 > pull_back(const Point< 3 > &p) const override
TorusManifold(const double R, const double r)
virtual DerivativeForm< 1, 3, 3 > push_forward_gradient(const Point< 3 > &chart_point) const override
virtual Point< 3 > push_forward(const Point< 3 > &chart_point) const override
virtual std::unique_ptr< Manifold< dim, 3 > > clone() const override
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
void initialize(const Triangulation< dim, spacedim > &triangulation)
Triangulation< dim, spacedim >::cell_iterator compute_chart_points(const ArrayView< const Point< spacedim > > &surrounding_points, ArrayView< Point< dim > > chart_points) const
std::array< unsigned int, 20 > get_possible_cells_around_points(const ArrayView< const Point< spacedim > > &surrounding_points) const
Point< dim > pull_back(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_guess) 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
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
Point< spacedim > push_forward(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &chart_point) const
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Triangulation< dim, spacedim > & get_triangulation()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Number signed_angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b, const Tensor< 1, spacedim, Number > &axis)
Tensor< 1, 3 > projected_direction(const Tensor< 1, 3 > &u, const Tensor< 1, 3 > &v)
Tensor< 1, 3 > apply_exponential_map(const Tensor< 1, 3 > &u, const Tensor< 1, 3 > &dir)
static constexpr double invalid_pull_back_coordinate
Point< spacedim > compute_normal(const Tensor< 1, spacedim > &, bool=false)
static constexpr double PI
static const unsigned int invalid_unsigned_int
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static double distance_to_unit_cell(const Point< dim > &p)
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)
static Point< dim, Number > project_to_unit_cell(const Point< dim, Number > &p)
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
const ::Triangulation< dim, spacedim > & tria