30#include <boost/container/small_vector.hpp>
53 const double theta = dir.
norm();
63 return tmp / tmp.
norm();
80 template <
int spacedim>
92 ExcMessage(
"The direction parameter must not be zero!"));
99 normal[0] = (vector[1] + vector[2]) / vector[0];
106 normal[1] = (vector[0] + vector[2]) / vector[1];
112 normal[2] = (vector[0] + vector[1]) / vector[2];
115 normal /= normal.
norm();
126template <
int dim,
int spacedim>
135template <
int dim,
int spacedim>
136std::unique_ptr<Manifold<dim, spacedim>>
139 return std::make_unique<PolarManifold<dim, spacedim>>(
center);
144template <
int dim,
int spacedim>
161template <
int dim,
int spacedim>
166 Assert(spherical_point[0] >= 0.0,
167 ExcMessage(
"Negative radius for given point."));
168 const double rho = spherical_point[0];
169 const double theta = spherical_point[1];
181 const double phi = spherical_point[2];
195template <
int dim,
int spacedim>
201 const double rho = R.
norm();
218 const double z = R[2];
234template <
int dim,
int spacedim>
239 Assert(spherical_point[0] >= 0.0,
240 ExcMessage(
"Negative radius for given point."));
241 const double rho = spherical_point[0];
242 const double theta = spherical_point[1];
259 const double phi = spherical_point[2];
284 template <
int dim,
int spacedim>
286 spherical_face_is_horizontal(
296 constexpr unsigned int n_vertices =
298 std::array<double, n_vertices> sqr_distances_to_center;
299 std::array<double, n_vertices - 1> sqr_distances_to_first_vertex;
300 sqr_distances_to_center[0] =
301 (face->vertex(0) - manifold_center).norm_square();
302 for (
unsigned int i = 1; i < n_vertices; ++i)
304 sqr_distances_to_center[i] =
305 (face->vertex(i) - manifold_center).norm_square();
306 sqr_distances_to_first_vertex[i - 1] =
307 (face->vertex(i) - face->vertex(0)).norm_square();
309 const auto minmax_sqr_distance =
310 std::minmax_element(sqr_distances_to_center.begin(),
311 sqr_distances_to_center.end());
312 const auto min_sqr_distance_to_first_vertex =
313 std::min_element(sqr_distances_to_first_vertex.begin(),
314 sqr_distances_to_first_vertex.end());
316 return (*minmax_sqr_distance.second - *minmax_sqr_distance.first <
317 1.e-10 * *min_sqr_distance_to_first_vertex);
323template <
int dim,
int spacedim>
333 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
340 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
341 return normalized_spherical_normal;
357template <
int dim,
int spacedim>
366template <
int dim,
int spacedim>
367std::unique_ptr<Manifold<dim, spacedim>>
370 return std::make_unique<SphericalManifold<dim, spacedim>>(
center);
375template <
int dim,
int spacedim>
380 const double w)
const
382 const double tol = 1
e-10;
384 if ((p1 - p2).norm_square() < tol * tol ||
std::abs(
w) < tol)
396 const double r1 =
v1.norm();
397 const double r2 = v2.
norm();
399 Assert(r1 > tol && r2 > tol,
400 ExcMessage(
"p1 and p2 cannot coincide with the center."));
406 const double cosgamma = e1 * e2;
424 const double n_norm = n.
norm();
427 "Probably, this means v1==v2 or v2==0."));
441template <
int dim,
int spacedim>
447 const double tol = 1
e-10;
454 const double r1 =
v1.norm();
455 const double r2 = v2.
norm();
465 const double cosgamma = e1 * e2;
468 ExcMessage(
"p1 and p2 cannot lie on the same diameter and be opposite "
469 "respect to the center."));
477 const double n_norm = n.
norm();
480 "Probably, this means v1==v2 or v2==0."));
487 return (r2 - r1) * e1 + r1 *
gamma * n;
492template <
int dim,
int spacedim>
502 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
509 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
510 return normalized_spherical_normal;
544template <
int dim,
int spacedim>
555 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
560 for (
unsigned int vertex = 0;
561 vertex < GeometryInfo<spacedim>::vertices_per_face;
563 face_vertex_normals[vertex] = face->vertex(vertex) -
center;
571template <
int dim,
int spacedim>
581 get_new_points(surrounding_points,
make_array_view(weights), new_points);
588template <
int dim,
int spacedim>
606template <
int dim,
int spacedim>
614 new_points.size() * surrounding_points.size());
615 const unsigned int weight_rows = new_points.size();
616 const unsigned int weight_columns = surrounding_points.size();
618 if (surrounding_points.size() == 2)
620 for (
unsigned int row = 0; row < weight_rows; ++row)
622 get_intermediate_point(surrounding_points[0],
623 surrounding_points[1],
624 weights[row * weight_columns + 1]);
628 boost::container::small_vector<std::pair<double, Tensor<1, spacedim>>, 100>
629 new_candidates(new_points.size());
630 boost::container::small_vector<Tensor<1, spacedim>, 100> directions(
632 boost::container::small_vector<double, 100> distances(
633 surrounding_points.size(), 0.0);
634 double max_distance = 0.;
635 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
637 directions[i] = surrounding_points[i] -
center;
638 distances[i] = directions[i].
norm();
640 if (distances[i] != 0.)
641 directions[i] /= distances[i];
644 ExcMessage(
"One of the vertices coincides with the center. "
645 "This is not allowed!"));
649 for (
unsigned int k = 0; k < i; ++k)
651 const double squared_distance =
652 (directions[i] - directions[k]).norm_square();
653 max_distance =
std::max(max_distance, squared_distance);
659 const double tolerance = 1
e-10;
660 boost::container::small_vector<bool, 100> accurate_point_was_found(
661 new_points.size(),
false);
666 for (
unsigned int row = 0; row < weight_rows; ++row)
668 new_candidates[row] =
669 guess_new_point(array_directions,
676 if (new_candidates[row].
first == 0.0)
679 accurate_point_was_found[row] =
true;
686 new_points[row] = polar_manifold.get_new_point(
701 if (max_distance < 2
e-2)
703 for (
unsigned int row = 0; row < weight_rows; ++row)
715 boost::container::small_vector<double, 1000> merged_weights(weights.
size());
716 boost::container::small_vector<Tensor<1, spacedim>, 100> merged_directions(
718 boost::container::small_vector<double, 100> merged_distances(
719 surrounding_points.size(), 0.0);
721 unsigned int n_unique_directions = 0;
722 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
724 bool found_duplicate =
false;
728 for (
unsigned int j = 0; j < n_unique_directions; ++j)
730 const double squared_distance =
731 (directions[i] - directions[j]).norm_square();
732 if (!found_duplicate && squared_distance < 1
e-28)
734 found_duplicate =
true;
735 for (
unsigned int row = 0; row < weight_rows; ++row)
736 merged_weights[row * weight_columns + j] +=
737 weights[row * weight_columns + i];
741 if (found_duplicate ==
false)
743 merged_directions[n_unique_directions] = directions[i];
744 merged_distances[n_unique_directions] = distances[i];
745 for (
unsigned int row = 0; row < weight_rows; ++row)
746 merged_weights[row * weight_columns + n_unique_directions] =
747 weights[row * weight_columns + i];
749 ++n_unique_directions;
755 boost::container::small_vector<unsigned int, 100> merged_weights_index(
757 for (
unsigned int row = 0; row < weight_rows; ++row)
759 for (
unsigned int existing_row = 0; existing_row < row; ++existing_row)
761 bool identical_weights =
true;
763 for (
unsigned int weight_index = 0;
764 weight_index < n_unique_directions;
766 if (
std::abs(merged_weights[row * weight_columns + weight_index] -
767 merged_weights[existing_row * weight_columns +
768 weight_index]) > tolerance)
770 identical_weights =
false;
774 if (identical_weights)
776 merged_weights_index[row] = existing_row;
786 merged_directions.begin() + n_unique_directions);
789 merged_distances.begin() + n_unique_directions);
791 for (
unsigned int row = 0; row < weight_rows; ++row)
792 if (!accurate_point_was_found[row])
797 &merged_weights[row * weight_columns], n_unique_directions);
798 new_candidates[row].second =
799 get_new_point(array_merged_directions,
800 array_merged_distances,
801 array_merged_weights,
805 new_candidates[row].second =
806 new_candidates[merged_weights_index[row]].second;
809 center + new_candidates[row].first * new_candidates[row].second;
815template <
int dim,
int spacedim>
816std::pair<double, Tensor<1, spacedim>>
822 const double tolerance = 1
e-10;
827 double total_weights = 0.;
828 for (
unsigned int i = 0; i < directions.size(); i++)
831 if (
std::abs(1 - weights[i]) < tolerance)
832 return std::make_pair(distances[i], directions[i]);
834 rho += distances[i] * weights[i];
835 candidate += directions[i] * weights[i];
836 total_weights += weights[i];
840 const double norm = candidate.
norm();
844 rho /= total_weights;
846 return std::make_pair(rho, candidate);
852 template <
int spacedim>
875 Point<3> candidate = candidate_point;
876 const unsigned int n_merged_points = directions.size();
877 const double tolerance = 1
e-10;
878 const int max_iterations = 10;
883 for (
unsigned int i = 0; i < n_merged_points; ++i)
885 const double squared_distance =
886 (candidate - directions[i]).norm_square();
887 if (squared_distance < tolerance * tolerance)
893 if (n_merged_points == 2)
912 for (
unsigned int i = 0; i < max_iterations; ++i)
926 for (
unsigned int i = 0; i < n_merged_points; ++i)
931 const double sintheta =
std::sqrt(sinthetaSq);
932 if (sintheta < tolerance)
934 Hessian[0][0] += weights[i];
935 Hessian[1][1] += weights[i];
939 const double costheta = (directions[i]) * candidate;
940 const double theta =
std::atan2(sintheta, costheta);
941 const double sincthetaInv = theta / sintheta;
943 const double cosphi = vPerp * Clocalx;
944 const double sinphi = vPerp * Clocaly;
946 gradlocal[0] = cosphi;
947 gradlocal[1] = sinphi;
948 gradient += (weights[i] * sincthetaInv) * gradlocal;
950 const double wt = weights[i] / sinthetaSq;
951 const double sinphiSq = sinphi * sinphi;
952 const double cosphiSq = cosphi * cosphi;
953 const double tt = sincthetaInv * costheta;
954 const double offdiag = cosphi * sinphi * wt * (1.0 - tt);
955 Hessian[0][0] += wt * (cosphiSq + tt * sinphiSq);
956 Hessian[0][1] += offdiag;
957 Hessian[1][0] += offdiag;
958 Hessian[1][1] += wt * (sinphiSq + tt * cosphiSq);
968 xDisplocal[0] * Clocalx + xDisplocal[1] * Clocaly;
971 const Point<3> candidateOld = candidate;
976 if ((candidate - candidateOld).norm_square() < tolerance * tolerance)
986template <
int dim,
int spacedim>
1006 const Point<3> & candidate_point)
const
1008 return do_get_new_point(directions, distances, weights, candidate_point);
1019 const Point<3> & candidate_point)
const
1021 return do_get_new_point(directions, distances, weights, candidate_point);
1032 const Point<3> & candidate_point)
const
1034 return do_get_new_point(directions, distances, weights, candidate_point);
1042template <
int dim,
int spacedim>
1044 const double tolerance)
1052 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1057template <
int dim,
int spacedim>
1061 const double tolerance)
1064 , direction(direction / direction.norm())
1065 , point_on_axis(point_on_axis)
1066 , tolerance(tolerance)
1071 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1076template <
int dim,
int spacedim>
1077std::unique_ptr<Manifold<dim, spacedim>>
1080 return std::make_unique<CylindricalManifold<dim, spacedim>>(direction,
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;
1132 const double dot = normal_direction * p_diff;
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);
1155 normal_direction * cosine_r + dxn * sine_r;
1158 return point_on_axis + direction * chart_point(2) + intermediate;
1163template <
int dim,
int spacedim>
1169 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1174 const double sine =
std::sin(chart_point(1));
1175 const double cosine =
std::cos(chart_point(1));
1178 normal_direction * cosine + dxn * sine;
1181 constexpr int s0 = 0 % spacedim;
1182 constexpr int s1 = 1 % spacedim;
1183 constexpr int s2 = 2 % spacedim;
1186 derivatives[s0][s0] = intermediate[s0];
1187 derivatives[s1][s0] = intermediate[s1];
1188 derivatives[s2][s0] = intermediate[s2];
1191 derivatives[s0][s1] = -normal_direction[s0] * sine + dxn[s0] * cosine;
1192 derivatives[s1][s1] = -normal_direction[s1] * sine + dxn[s1] * cosine;
1193 derivatives[s2][s1] = -normal_direction[s2] * sine + dxn[s2] * cosine;
1196 derivatives[s0][s2] = direction[s0];
1197 derivatives[s1][s2] = direction[s1];
1198 derivatives[s2][s2] = direction[s2];
1208template <
int dim,
int spacedim>
1212 const double eccentricity)
1215 , direction(major_axis_direction)
1217 , cosh_u(1.0 / eccentricity)
1218 , sinh_u(
std::
sqrt(cosh_u * cosh_u - 1.0))
1225 "Invalid eccentricity: It must satisfy 0 < eccentricity < 1."));
1227 Assert(direction_norm != 0.0,
1229 "Invalid major axis direction vector: Null vector not allowed."));
1235template <
int dim,
int spacedim>
1236std::unique_ptr<Manifold<dim, spacedim>>
1239 const double eccentricity = 1.0 / cosh_u;
1240 return std::make_unique<EllipticalManifold<dim, spacedim>>(
center,
1247template <
int dim,
int spacedim>
1260template <
int dim,
int spacedim>
1274 const double cs =
std::cos(chart_point[1]);
1275 const double sn =
std::sin(chart_point[1]);
1278 const double x = chart_point[0] * cosh_u * cs;
1279 const double y = chart_point[0] * sinh_u * sn;
1281 const Point<2> p(direction[0] * x - direction[1] * y,
1282 direction[1] * x + direction[0] * y);
1288template <
int dim,
int spacedim>
1303 const double x0 = space_point[0] -
center[0];
1304 const double y0 = space_point[1] -
center[1];
1305 const double x = direction[0] * x0 + direction[1] * y0;
1306 const double y = -direction[1] * x0 + direction[0] * y0;
1308 std::sqrt((x * x) / (cosh_u * cosh_u) + (y * y) / (sinh_u * sinh_u));
1314 double cos_eta = x / (pt0 * cosh_u);
1324 const double pt1 = (std::signbit(y) ? 2.0 *
numbers::PI - eta : eta);
1330template <
int dim,
int spacedim>
1346 const double cs =
std::cos(chart_point[1]);
1347 const double sn =
std::sin(chart_point[1]);
1349 dX[0][0] = cosh_u * cs;
1350 dX[0][1] = -chart_point[0] * cosh_u * sn;
1351 dX[1][0] = sinh_u * sn;
1352 dX[1][1] = chart_point[0] * sinh_u * cs;
1356 {{+direction[0], -direction[1]}, {direction[1], direction[0]}}};
1366template <
int dim,
int spacedim,
int chartdim>
1371 const double tolerance)
1374 , push_forward_function(&push_forward_function)
1375 , pull_back_function(&pull_back_function)
1376 , tolerance(tolerance)
1377 , owns_pointers(false)
1378 , finite_difference_step(0)
1386template <
int dim,
int spacedim,
int chartdim>
1391 const double tolerance)
1395 , pull_back_function(
pull_back.release())
1396 , tolerance(tolerance)
1397 , owns_pointers(true)
1398 , finite_difference_step(0)
1406template <
int dim,
int spacedim,
int chartdim>
1408 const std::string push_forward_expression,
1409 const std::string pull_back_expression,
1412 const std::string chart_vars,
1413 const std::string space_vars,
1414 const double tolerance,
1417 , const_map(const_map)
1418 , tolerance(tolerance)
1419 , owns_pointers(true)
1420 , push_forward_expression(push_forward_expression)
1421 , pull_back_expression(pull_back_expression)
1422 , chart_vars(chart_vars)
1423 , space_vars(space_vars)
1424 , finite_difference_step(h)
1436template <
int dim,
int spacedim,
int chartdim>
1439 if (owns_pointers ==
true)
1442 push_forward_function =
nullptr;
1446 pull_back_function =
nullptr;
1453template <
int dim,
int spacedim,
int chartdim>
1454std::unique_ptr<Manifold<dim, spacedim>>
1469 if (!(push_forward_expression.empty() && pull_back_expression.empty()))
1471 return std::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1472 push_forward_expression,
1473 pull_back_expression,
1474 this->get_periodicity(),
1479 finite_difference_step);
1483 return std::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1484 *push_forward_function,
1485 *pull_back_function,
1486 this->get_periodicity(),
1493template <
int dim,
int spacedim,
int chartdim>
1500 push_forward_function->vector_value(chart_point, pf);
1501 for (
unsigned int i = 0; i < spacedim; ++i)
1506 pull_back_function->vector_value(result, pb);
1507 for (
unsigned int i = 0; i < chartdim; ++i)
1509 (chart_point.
norm() > tolerance &&
1510 (
std::abs(pb[i] - chart_point[i]) < tolerance * chart_point.
norm())) ||
1511 (
std::abs(pb[i] - chart_point[i]) < tolerance),
1513 "The push forward is not the inverse of the pull back! Bailing out."));
1521template <
int dim,
int spacedim,
int chartdim>
1527 for (
unsigned int i = 0; i < spacedim; ++i)
1529 const auto gradient = push_forward_function->gradient(chart_point, i);
1530 for (
unsigned int j = 0; j < chartdim; ++j)
1531 DF[i][j] = gradient[j];
1538template <
int dim,
int spacedim,
int chartdim>
1545 pull_back_function->vector_value(space_point, pb);
1546 for (
unsigned int i = 0; i < chartdim; ++i)
1568 return {phi, theta,
w};
1577 double phi = chart_point(0);
1578 double theta = chart_point(1);
1579 double w = chart_point(2);
1595 ExcMessage(
"Outer radius R must be greater than the inner "
1603std::unique_ptr<Manifold<dim, 3>>
1606 return std::make_unique<TorusManifold<dim>>(R, r);
1617 double phi = chart_point(0);
1618 double theta = chart_point(1);
1619 double w = chart_point(2);
1641template <
int dim,
int spacedim>
1652template <
int dim,
int spacedim>
1654 spacedim>::~TransfiniteInterpolationManifold()
1656 if (clear_signal.connected())
1657 clear_signal.disconnect();
1662template <
int dim,
int spacedim>
1663std::unique_ptr<Manifold<dim, spacedim>>
1669 return std::unique_ptr<Manifold<dim, spacedim>>(ptr);
1674template <
int dim,
int spacedim>
1681 clear_signal.disconnect();
1682 clear_signal =
triangulation.signals.clear.connect([&]() ->
void {
1683 this->triangulation =
nullptr;
1684 this->level_coarse = -1;
1687 coarse_cell_is_flat.resize(
triangulation.n_cells(level_coarse),
false);
1688 quadratic_approximation.clear();
1690 std::vector<Point<dim>> unit_points =
1692 std::vector<Point<spacedim>> real_points(unit_points.size());
1694 for (
const auto &cell :
triangulation.active_cell_iterators())
1696 bool cell_is_flat =
true;
1697 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1698 if (cell->line(
l)->manifold_id() != cell->manifold_id() &&
1700 cell_is_flat =
false;
1702 for (
unsigned int q = 0; q < GeometryInfo<dim>::quads_per_cell; ++q)
1703 if (cell->quad(q)->manifold_id() != cell->manifold_id() &&
1705 cell_is_flat =
false;
1707 coarse_cell_is_flat.size());
1708 coarse_cell_is_flat[cell->index()] = cell_is_flat;
1711 for (
unsigned int i = 0; i < unit_points.size(); ++i)
1713 quadratic_approximation.emplace_back(real_points, unit_points);
1722 template <
typename AccessorType>
1724 compute_transfinite_interpolation(
const AccessorType &cell,
1728 return cell.vertex(0) * (1. - chart_point[0]) +
1729 cell.vertex(1) * chart_point[0];
1733 template <
typename AccessorType>
1735 compute_transfinite_interpolation(
const AccessorType &cell,
1737 const bool cell_is_flat)
1739 const unsigned int dim = AccessorType::dimension;
1740 const unsigned int spacedim = AccessorType::space_dimension;
1748 const std::array<Point<spacedim>, 4>
vertices{
1749 {cell.vertex(0), cell.vertex(1), cell.vertex(2), cell.vertex(3)}};
1754 std::array<double, 4> weights_vertices{
1755 {(1. - chart_point[0]) * (1. - chart_point[1]),
1756 chart_point[0] * (1. - chart_point[1]),
1757 (1. - chart_point[0]) * chart_point[1],
1758 chart_point[0] * chart_point[1]}};
1763 new_point += weights_vertices[v] *
vertices[v];
1775 std::array<double, GeometryInfo<2>::vertices_per_face> weights;
1778 const auto weights_view =
1780 const auto points_view =
make_array_view(points.begin(), points.end());
1782 for (
unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1785 const double my_weight =
1786 (line % 2) ? chart_point[line / 2] : 1 - chart_point[line / 2];
1787 const double line_point = chart_point[1 - line / 2];
1794 cell.line(line)->manifold_id();
1795 if (line_manifold_id == my_manifold_id ||
1800 my_weight * (1. - line_point);
1803 my_weight * line_point;
1811 weights[0] = 1. - line_point;
1812 weights[1] = line_point;
1821 new_point -= weights_vertices[v] *
vertices[v];
1830 static constexpr unsigned int face_to_cell_vertices_3d[6][4] = {{0, 2, 4, 6},
1840 static constexpr unsigned int face_to_cell_lines_3d[6][4] = {{8, 10, 0, 4},
1848 template <
typename AccessorType>
1850 compute_transfinite_interpolation(
const AccessorType &cell,
1852 const bool cell_is_flat)
1854 const unsigned int dim = AccessorType::dimension;
1855 const unsigned int spacedim = AccessorType::space_dimension;
1861 const std::array<Point<spacedim>, 8>
vertices{{cell.vertex(0),
1873 double linear_shapes[10];
1874 for (
unsigned int d = 0;
d < 3; ++
d)
1876 linear_shapes[2 *
d] = 1. - chart_point[
d];
1877 linear_shapes[2 *
d + 1] = chart_point[
d];
1881 for (
unsigned int d = 6;
d < 10; ++
d)
1882 linear_shapes[
d] = linear_shapes[
d - 6];
1884 std::array<double, 8> weights_vertices;
1885 for (
unsigned int i2 = 0, v = 0; i2 < 2; ++i2)
1886 for (
unsigned int i1 = 0; i1 < 2; ++i1)
1887 for (
unsigned int i0 = 0; i0 < 2; ++i0, ++v)
1888 weights_vertices[v] =
1889 (linear_shapes[4 + i2] * linear_shapes[2 + i1]) * linear_shapes[i0];
1893 for (
unsigned int v = 0; v < 8; ++v)
1894 new_point += weights_vertices[v] *
vertices[v];
1900 std::array<double, GeometryInfo<3>::lines_per_cell> weights_lines;
1901 std::fill(weights_lines.begin(), weights_lines.end(), 0.0);
1904 std::array<double, GeometryInfo<2>::vertices_per_cell> weights;
1907 const auto weights_view =
1909 const auto points_view =
make_array_view(points.begin(), points.end());
1913 const double my_weight = linear_shapes[face];
1914 const unsigned int face_even = face - face % 2;
1922 cell.face(face)->manifold_id();
1923 if (face_manifold_id == my_manifold_id ||
1926 for (
unsigned int line = 0;
1927 line < GeometryInfo<2>::lines_per_cell;
1930 const double line_weight =
1931 linear_shapes[face_even + 2 + line];
1932 weights_lines[face_to_cell_lines_3d[face][line]] +=
1933 my_weight * line_weight;
1939 weights_vertices[face_to_cell_vertices_3d[face][0]] -=
1940 linear_shapes[face_even + 2] *
1941 (linear_shapes[face_even + 4] * my_weight);
1942 weights_vertices[face_to_cell_vertices_3d[face][1]] -=
1943 linear_shapes[face_even + 3] *
1944 (linear_shapes[face_even + 4] * my_weight);
1945 weights_vertices[face_to_cell_vertices_3d[face][2]] -=
1946 linear_shapes[face_even + 2] *
1947 (linear_shapes[face_even + 5] * my_weight);
1948 weights_vertices[face_to_cell_vertices_3d[face][3]] -=
1949 linear_shapes[face_even + 3] *
1950 (linear_shapes[face_even + 5] * my_weight);
1955 points[v] =
vertices[face_to_cell_vertices_3d[face][v]];
1957 linear_shapes[face_even + 2] * linear_shapes[face_even + 4];
1959 linear_shapes[face_even + 3] * linear_shapes[face_even + 4];
1961 linear_shapes[face_even + 2] * linear_shapes[face_even + 5];
1963 linear_shapes[face_even + 3] * linear_shapes[face_even + 5];
1971 const auto weights_view_line =
1973 const auto points_view_line =
1975 for (
unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
1978 const double line_point =
1979 (line < 8 ? chart_point[1 - (line % 4) / 2] : chart_point[2]);
1980 double my_weight = 0.;
1982 my_weight = linear_shapes[line % 4] * linear_shapes[4 + line / 4];
1985 const unsigned int subline = line - 8;
1987 linear_shapes[subline % 2] * linear_shapes[2 + subline / 2];
1989 my_weight -= weights_lines[line];
1995 cell.line(line)->manifold_id();
1996 if (line_manifold_id == my_manifold_id ||
2001 my_weight * (1. - line_point);
2004 my_weight * (line_point);
2012 weights[0] = 1. - line_point;
2013 weights[1] = line_point;
2014 new_point -= my_weight * tria.
get_manifold(line_manifold_id)
2022 new_point += weights_vertices[v] *
vertices[v];
2030template <
int dim,
int spacedim>
2042 ExcMessage(
"chart_point is not in unit interval"));
2044 return compute_transfinite_interpolation(*cell,
2046 coarse_cell_is_flat[cell->index()]);
2051template <
int dim,
int spacedim>
2060 for (
unsigned int d = 0;
d < dim; ++
d)
2063 const double step = chart_point[
d] > 0.5 ? -1
e-8 : 1
e-8;
2066 modified[
d] += step;
2068 compute_transfinite_interpolation(*cell,
2070 coarse_cell_is_flat[cell->index()]) -
2071 pushed_forward_chart_point;
2072 for (
unsigned int e = 0;
e < spacedim; ++
e)
2073 grad[
e][
d] = difference[
e] / step;
2080template <
int dim,
int spacedim>
2088 for (
unsigned int d = 0;
d < dim; ++
d)
2102 compute_transfinite_interpolation(*cell,
2104 coarse_cell_is_flat[cell->index()]);
2105 const double tolerance = 1
e-21 * Utilities::fixed_power<2>(cell->diameter());
2106 double residual_norm_square = residual.
norm_square();
2108 bool must_recompute_jacobian =
true;
2109 for (
unsigned int i = 0; i < 100; ++i)
2111 if (residual_norm_square < tolerance)
2118 for (
unsigned int d = 0;
d < spacedim; ++
d)
2119 for (
unsigned int e = 0;
e < dim; ++
e)
2120 update[
e] += inv_grad[
d][
e] * residual[
d];
2121 return chart_point + update;
2133 if (must_recompute_jacobian || i % 9 == 0)
2141 push_forward_gradient(cell,
2147 must_recompute_jacobian =
false;
2150 for (
unsigned int d = 0;
d < spacedim; ++
d)
2151 for (
unsigned int e = 0;
e < dim; ++
e)
2152 update[
e] += inv_grad[
d][
e] * residual[
d];
2166 while (alpha > 1
e-4)
2168 Point<dim> guess = chart_point + alpha * update;
2170 point - compute_transfinite_interpolation(
2171 *cell, guess, coarse_cell_is_flat[cell->index()]);
2172 const double residual_norm_new = residual_guess.
norm_square();
2173 if (residual_norm_new < residual_norm_square)
2175 residual = residual_guess;
2176 residual_norm_square = residual_norm_new;
2177 chart_point += alpha * update;
2190 must_recompute_jacobian =
true;
2204 for (
unsigned int d = 0;
d < spacedim; ++
d)
2205 for (
unsigned int e = 0;
e < dim; ++
e)
2206 Jinv_deltaf[
e] += inv_grad[
d][
e] * delta_f[
d];
2213 if (
std::abs(delta_x * Jinv_deltaf) > 1
e-12 && !must_recompute_jacobian)
2216 (delta_x - Jinv_deltaf) / (delta_x * Jinv_deltaf);
2218 for (
unsigned int d = 0;
d < spacedim; ++
d)
2219 for (
unsigned int e = 0;
e < dim; ++
e)
2220 jac_update[
d] += delta_x[
e] * inv_grad[
d][
e];
2221 for (
unsigned int d = 0;
d < spacedim; ++
d)
2222 for (
unsigned int e = 0;
e < dim; ++
e)
2223 inv_grad[
d][
e] += factor[
e] * jac_update[
d];
2231template <
int dim,
int spacedim>
2232std::array<unsigned int, 20>
2242 ExcMessage(
"The manifold was initialized with level " +
2244 "active cells on a lower level. Coarsening the mesh is " +
2245 "currently not supported"));
2255 boost::container::small_vector<std::pair<double, unsigned int>, 200>
2256 distances_and_cells;
2257 for (; cell != endc; ++cell)
2267 vertices[vertex_n] = cell->vertex(vertex_n);
2277 double radius_square = 0.;
2281 bool inside_circle =
true;
2282 for (
unsigned int i = 0; i < points.size(); ++i)
2283 if ((
center - points[i]).norm_square() > radius_square * 1.5)
2285 inside_circle =
false;
2288 if (inside_circle ==
false)
2292 double current_distance = 0;
2293 for (
unsigned int i = 0; i < points.size(); ++i)
2296 quadratic_approximation[cell->index()].compute(points[i]);
2299 distances_and_cells.push_back(
2300 std::make_pair(current_distance, cell->index()));
2305 std::sort(distances_and_cells.begin(), distances_and_cells.end());
2306 std::array<unsigned int, 20> cells;
2308 for (
unsigned int i = 0; i < distances_and_cells.size() && i < cells.size();
2310 cells[i] = distances_and_cells[i].second;
2317template <
int dim,
int spacedim>
2323 Assert(surrounding_points.size() == chart_points.size(),
2324 ExcMessage(
"The chart points array view must be as large as the "
2325 "surrounding points array view."));
2327 std::array<unsigned int, 20> nearby_cells =
2328 get_possible_cells_around_points(surrounding_points);
2352 auto guess_chart_point_structdim_2 = [&](
const unsigned int i) ->
Point<dim> {
2353 Assert(surrounding_points.size() == 8 && 2 < i && i < 8,
2354 ExcMessage(
"This function assumes that there are eight surrounding "
2355 "points around a two-dimensional object. It also assumes "
2356 "that the first three chart points have already been "
2366 return chart_points[1] + (chart_points[2] - chart_points[0]);
2368 return 0.5 * (chart_points[0] + chart_points[2]);
2370 return 0.5 * (chart_points[1] + chart_points[3]);
2372 return 0.5 * (chart_points[0] + chart_points[1]);
2374 return 0.5 * (chart_points[2] + chart_points[3]);
2403 auto guess_chart_point_structdim_3 = [&](
const unsigned int i) ->
Point<dim> {
2404 Assert(surrounding_points.size() == 8 && 4 < i && i < 8,
2405 ExcMessage(
"This function assumes that there are eight surrounding "
2406 "points around a three-dimensional object. It also "
2407 "assumes that the first five chart points have already "
2409 return chart_points[i - 4] + (chart_points[4] - chart_points[0]);
2413 bool use_structdim_2_guesses =
false;
2414 bool use_structdim_3_guesses =
false;
2419 if (surrounding_points.size() == 8)
2422 surrounding_points[6] - surrounding_points[0];
2424 surrounding_points[7] - surrounding_points[2];
2432 use_structdim_2_guesses =
true;
2433 else if (spacedim == 3)
2436 use_structdim_3_guesses =
true;
2439 Assert((!use_structdim_2_guesses && !use_structdim_3_guesses) ||
2440 (use_structdim_2_guesses ^ use_structdim_3_guesses),
2445 auto compute_chart_point =
2447 const unsigned int point_index) {
2453 bool used_quadratic_approximation =
false;
2456 if (point_index == 3 && surrounding_points.size() >= 8)
2457 guess = chart_points[1] + (chart_points[2] - chart_points[0]);
2458 else if (use_structdim_2_guesses && 3 < point_index)
2459 guess = guess_chart_point_structdim_2(point_index);
2460 else if (use_structdim_3_guesses && 4 < point_index)
2461 guess = guess_chart_point_structdim_3(point_index);
2462 else if (dim == 3 && point_index > 7 && surrounding_points.size() == 26)
2464 if (point_index < 20)
2467 point_index - 8, 0)] +
2469 point_index - 8, 1)]);
2473 point_index - 20, 0)] +
2475 point_index - 20, 1)] +
2477 point_index - 20, 2)] +
2479 point_index - 20, 3)]);
2483 guess = quadratic_approximation[cell->index()].compute(
2484 surrounding_points[point_index]);
2485 used_quadratic_approximation =
true;
2487 chart_points[point_index] =
2488 pull_back(cell, surrounding_points[point_index], guess);
2493 if (chart_points[point_index][0] ==
2495 !used_quadratic_approximation)
2497 guess = quadratic_approximation[cell->index()].compute(
2498 surrounding_points[point_index]);
2499 chart_points[point_index] =
2500 pull_back(cell, surrounding_points[point_index], guess);
2503 if (chart_points[point_index][0] ==
2506 for (
unsigned int d = 0;
d < dim; ++
d)
2508 chart_points[point_index] =
2509 pull_back(cell, surrounding_points[point_index], guess);
2514 for (
unsigned int c = 0; c < nearby_cells.size(); ++c)
2518 bool inside_unit_cell =
true;
2519 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2521 compute_chart_point(cell, i);
2528 inside_unit_cell =
false;
2532 if (inside_unit_cell ==
true)
2539 if (c == nearby_cells.size() - 1 ||
2544 std::ostringstream message;
2545 for (
unsigned int b = 0;
b <= c; ++
b)
2549 message <<
"Looking at cell " << cell->id()
2550 <<
" with vertices: " << std::endl;
2552 message << cell->vertex(v) <<
" ";
2553 message << std::endl;
2554 message <<
"Transformation to chart coordinates: " << std::endl;
2555 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2557 compute_chart_point(cell, i);
2558 message << surrounding_points[i] <<
" -> " << chart_points[i]
2578template <
int dim,
int spacedim>
2584 boost::container::small_vector<Point<dim>, 100> chart_points(
2585 surrounding_points.size());
2588 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2591 chart_manifold.get_new_point(chart_points_view, weights);
2598template <
int dim,
int spacedim>
2608 boost::container::small_vector<Point<dim>, 100> chart_points(
2609 surrounding_points.size());
2612 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2614 boost::container::small_vector<Point<dim>, 100> new_points_on_chart(
2616 chart_manifold.get_new_points(chart_points_view,
2619 new_points_on_chart.end()));
2621 for (
unsigned int row = 0; row < weights.size(0); ++row)
2622 new_points[row] =
push_forward(cell, new_points_on_chart[row]);
2628#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
std::map< std::string, double > ConstMap
void initialize(const std::string &vars, const std::vector< std::string > &expressions, const ConstMap &constants, const bool time_dependent=false)
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
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &t)
constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &t)
constexpr ProductType< Number, OtherNumber >::type scalar_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, OtherNumber > &t2)
numbers::NumberTraits< Number >::real_type norm() const
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
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()
VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &x)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
std::string to_string(const T &t)
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
Expression atan2(const Expression &y, const Expression &x)
Expression acos(const Expression &x)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
long double gamma(const unsigned int n)
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 > pow(const ::VectorizedArray< Number, width > &, const Number p)
::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)