48 const double theta = dir.
norm();
57 std::cos(theta) * u + std::sin(theta) * dirUnit;
58 return tmp / tmp.
norm();
75 template <
int spacedim>
87 ExcMessage(
"The direction parameter must not be zero!"));
89 if (std::abs(vector[0]) >= std::abs(vector[1]) &&
90 std::abs(vector[0]) >= std::abs(vector[2]))
94 normal[0] = (vector[1] + vector[2]) / vector[0];
96 else if (std::abs(vector[1]) >= std::abs(vector[0]) &&
97 std::abs(vector[1]) >= std::abs(vector[2]))
101 normal[1] = (vector[0] + vector[2]) / vector[1];
107 normal[2] = (vector[0] + vector[1]) / vector[2];
110 normal /= normal.
norm();
121 template <
int dim,
int spacedim>
130 template <
int dim,
int spacedim>
131 std::unique_ptr<Manifold<dim, spacedim>>
134 return std_cxx14::make_unique<PolarManifold<dim, spacedim>>(
center);
139 template <
int dim,
int spacedim>
156 template <
int dim,
int spacedim>
161 Assert(spherical_point[0] >= 0.0,
162 ExcMessage(
"Negative radius for given point."));
163 const double rho = spherical_point[0];
164 const double theta = spherical_point[1];
171 p[0] = rho * std::cos(theta);
172 p[1] = rho * std::sin(theta);
176 const double phi = spherical_point[2];
177 p[0] = rho * std::sin(theta) * std::cos(phi);
178 p[1] = rho * std::sin(theta) * std::sin(phi);
179 p[2] = rho * std::cos(theta);
190 template <
int dim,
int spacedim>
196 const double rho = R.
norm();
213 const double z = R[2];
217 p[1] =
std::atan2(std::sqrt(R[0] * R[0] + R[1] * R[1]), z);
229 template <
int dim,
int spacedim>
234 Assert(spherical_point[0] >= 0.0,
235 ExcMessage(
"Negative radius for given point."));
236 const double rho = spherical_point[0];
237 const double theta = spherical_point[1];
245 DX[0][0] = std::cos(theta);
246 DX[0][1] = -rho * std::sin(theta);
247 DX[1][0] = std::sin(theta);
248 DX[1][1] = rho * std::cos(theta);
254 const double phi = spherical_point[2];
255 DX[0][0] = std::sin(theta) * std::cos(phi);
256 DX[0][1] = rho * std::cos(theta) * std::cos(phi);
257 DX[0][2] = -rho * std::sin(theta) * std::sin(phi);
259 DX[1][0] = std::sin(theta) * std::sin(phi);
260 DX[1][1] = rho * std::cos(theta) * std::sin(phi);
261 DX[1][2] = rho * std::sin(theta) * std::cos(phi);
263 DX[2][0] = std::cos(theta);
264 DX[2][1] = -rho * std::sin(theta);
279 template <
int dim,
int spacedim>
281 spherical_face_is_horizontal(
291 constexpr
unsigned int n_vertices =
293 std::array<double, n_vertices> sqr_distances_to_center;
294 std::array<
double, n_vertices - 1> sqr_distances_to_first_vertex;
295 sqr_distances_to_center[0] =
296 (face->vertex(0) - manifold_center).norm_square();
297 for (
unsigned int i = 1; i < n_vertices; ++i)
299 sqr_distances_to_center[i] =
300 (face->vertex(i) - manifold_center).norm_square();
301 sqr_distances_to_first_vertex[i - 1] =
302 (face->vertex(i) - face->vertex(0)).norm_square();
304 const auto minmax_sqr_distance =
305 std::minmax_element(sqr_distances_to_center.begin(),
306 sqr_distances_to_center.end());
307 const auto min_sqr_distance_to_first_vertex =
308 std::min_element(sqr_distances_to_first_vertex.begin(),
309 sqr_distances_to_first_vertex.end());
311 return (*minmax_sqr_distance.second - *minmax_sqr_distance.first <
312 1.e-10 * *min_sqr_distance_to_first_vertex);
318 template <
int dim,
int spacedim>
328 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
335 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
336 return normalized_spherical_normal;
352 template <
int dim,
int spacedim>
361 template <
int dim,
int spacedim>
362 std::unique_ptr<Manifold<dim, spacedim>>
365 return std_cxx14::make_unique<SphericalManifold<dim, spacedim>>(
center);
370 template <
int dim,
int spacedim>
375 const double w)
const
377 const double tol = 1
e-10;
379 if ((p1 - p2).norm_square() < tol * tol || std::abs(
w) < tol)
381 else if (std::abs(
w - 1.0) < tol)
391 const double r1 = v1.
norm();
392 const double r2 = v2.
norm();
394 Assert(r1 > tol && r2 > tol,
395 ExcMessage(
"p1 and p2 cannot coincide with the center."));
401 const double cosgamma = e1 * e2;
419 const double n_norm = n.
norm();
422 "Probably, this means v1==v2 or v2==0."));
436 template <
int dim,
int spacedim>
442 const double tol = 1
e-10;
449 const double r1 = v1.
norm();
450 const double r2 = v2.
norm();
460 const double cosgamma = e1 * e2;
463 ExcMessage(
"p1 and p2 cannot lie on the same diameter and be opposite "
464 "respect to the center."));
472 const double n_norm = n.
norm();
475 "Probably, this means v1==v2 or v2==0."));
482 return (r2 - r1) * e1 + r1 *
gamma * n;
487 template <
int dim,
int spacedim>
497 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
504 unnormalized_spherical_normal / unnormalized_spherical_normal.
norm();
505 return normalized_spherical_normal;
539 template <
int dim,
int spacedim>
550 if (spherical_face_is_horizontal<dim, spacedim>(face,
center))
555 for (
unsigned int vertex = 0;
556 vertex < GeometryInfo<spacedim>::vertices_per_face;
558 face_vertex_normals[vertex] = face->vertex(vertex) -
center;
566 template <
int dim,
int spacedim>
576 get_new_points(surrounding_points,
make_array_view(weights), new_points);
583 template <
int dim,
int spacedim>
601 template <
int dim,
int spacedim>
609 new_points.size() * surrounding_points.size());
610 const unsigned int weight_rows = new_points.size();
611 const unsigned int weight_columns = surrounding_points.size();
613 if (surrounding_points.size() == 2)
615 for (
unsigned int row = 0; row < weight_rows; ++row)
617 get_intermediate_point(surrounding_points[0],
618 surrounding_points[1],
619 weights[row * weight_columns + 1]);
623 boost::container::small_vector<std::pair<double, Tensor<1, spacedim>>, 100>
624 new_candidates(new_points.size());
625 boost::container::small_vector<Tensor<1, spacedim>, 100> directions(
627 boost::container::small_vector<double, 100> distances(
628 surrounding_points.size(), 0.0);
629 double max_distance = 0.;
630 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
632 directions[i] = surrounding_points[i] -
center;
633 distances[i] = directions[i].
norm();
635 if (distances[i] != 0.)
636 directions[i] /= distances[i];
639 ExcMessage(
"One of the vertices coincides with the center. "
640 "This is not allowed!"));
644 for (
unsigned int k = 0; k < i; ++k)
646 const double squared_distance =
647 (directions[i] - directions[k]).norm_square();
648 max_distance =
std::max(max_distance, squared_distance);
654 const double tolerance = 1
e-10;
655 boost::container::small_vector<bool, 100> accurate_point_was_found(
656 new_points.size(),
false);
661 for (
unsigned int row = 0; row < weight_rows; ++row)
663 new_candidates[row] =
664 guess_new_point(array_directions,
671 if (new_candidates[row].
first == 0.0)
674 accurate_point_was_found[row] =
true;
681 new_points[row] = polar_manifold.get_new_point(
696 if (max_distance < 2
e-2)
698 for (
unsigned int row = 0; row < weight_rows; ++row)
710 boost::container::small_vector<double, 1000> merged_weights(weights.
size());
711 boost::container::small_vector<Tensor<1, spacedim>, 100> merged_directions(
713 boost::container::small_vector<double, 100> merged_distances(
714 surrounding_points.size(), 0.0);
716 unsigned int n_unique_directions = 0;
717 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
719 bool found_duplicate =
false;
723 for (
unsigned int j = 0; j < n_unique_directions; ++j)
725 const double squared_distance =
726 (directions[i] - directions[j]).norm_square();
727 if (!found_duplicate && squared_distance < 1
e-28)
729 found_duplicate =
true;
730 for (
unsigned int row = 0; row < weight_rows; ++row)
731 merged_weights[row * weight_columns + j] +=
732 weights[row * weight_columns + i];
736 if (found_duplicate ==
false)
738 merged_directions[n_unique_directions] = directions[i];
739 merged_distances[n_unique_directions] = distances[i];
740 for (
unsigned int row = 0; row < weight_rows; ++row)
741 merged_weights[row * weight_columns + n_unique_directions] =
742 weights[row * weight_columns + i];
744 ++n_unique_directions;
750 boost::container::small_vector<unsigned int, 100> merged_weights_index(
752 for (
unsigned int row = 0; row < weight_rows; ++row)
754 for (
unsigned int existing_row = 0; existing_row < row; ++existing_row)
756 bool identical_weights =
true;
758 for (
unsigned int weight_index = 0;
759 weight_index < n_unique_directions;
761 if (std::abs(merged_weights[row * weight_columns + weight_index] -
762 merged_weights[existing_row * weight_columns +
763 weight_index]) > tolerance)
765 identical_weights =
false;
769 if (identical_weights)
771 merged_weights_index[row] = existing_row;
781 merged_directions.begin() + n_unique_directions);
784 merged_distances.begin() + n_unique_directions);
786 for (
unsigned int row = 0; row < weight_rows; ++row)
787 if (!accurate_point_was_found[row])
792 &merged_weights[row * weight_columns], n_unique_directions);
793 new_candidates[row].second =
794 get_new_point(array_merged_directions,
795 array_merged_distances,
796 array_merged_weights,
800 new_candidates[row].second =
801 new_candidates[merged_weights_index[row]].second;
804 center + new_candidates[row].first * new_candidates[row].second;
810 template <
int dim,
int spacedim>
811 std::pair<double, Tensor<1, spacedim>>
817 const double tolerance = 1
e-10;
822 double total_weights = 0.;
823 for (
unsigned int i = 0; i < directions.size(); i++)
826 if (std::abs(1 - weights[i]) < tolerance)
827 return std::make_pair(distances[i], directions[i]);
829 rho += distances[i] * weights[i];
830 candidate += directions[i] * weights[i];
831 total_weights += weights[i];
835 const double norm = candidate.
norm();
839 rho /= total_weights;
841 return std::make_pair(rho, candidate);
847 template <
int spacedim>
870 Point<3> candidate = candidate_point;
871 const unsigned int n_merged_points = directions.size();
872 const double tolerance = 1
e-10;
873 const int max_iterations = 10;
878 for (
unsigned int i = 0; i < n_merged_points; ++i)
880 const double squared_distance =
881 (candidate - directions[i]).norm_square();
882 if (squared_distance < tolerance * tolerance)
888 if (n_merged_points == 2)
891 Assert(std::abs(weights[0] + weights[1] - 1.0) < 1
e-13,
907 for (
unsigned int i = 0; i < max_iterations; ++i)
921 for (
unsigned int i = 0; i < n_merged_points; ++i)
922 if (std::abs(weights[i]) > 1.
e-15)
926 const double sintheta =
std::sqrt(sinthetaSq);
927 if (sintheta < tolerance)
929 Hessian[0][0] += weights[i];
930 Hessian[1][1] += weights[i];
934 const double costheta = (directions[i]) * candidate;
935 const double theta =
std::atan2(sintheta, costheta);
936 const double sincthetaInv = theta / sintheta;
938 const double cosphi = vPerp * Clocalx;
939 const double sinphi = vPerp * Clocaly;
941 gradlocal[0] = cosphi;
942 gradlocal[1] = sinphi;
943 gradient += (weights[i] * sincthetaInv) * gradlocal;
945 const double wt = weights[i] / sinthetaSq;
946 const double sinphiSq = sinphi * sinphi;
947 const double cosphiSq = cosphi * cosphi;
948 const double tt = sincthetaInv * costheta;
949 const double offdiag = cosphi * sinphi * wt * (1.0 - tt);
950 Hessian[0][0] += wt * (cosphiSq + tt * sinphiSq);
951 Hessian[0][1] += offdiag;
952 Hessian[1][0] += offdiag;
953 Hessian[1][1] += wt * (sinphiSq + tt * cosphiSq);
961 const Tensor<1, 2> xDisplocal = inverse_Hessian * gradient;
963 xDisplocal[0] * Clocalx + xDisplocal[1] * Clocaly;
966 const Point<3> candidateOld = candidate;
971 if ((candidate - candidateOld).norm_square() < tolerance * tolerance)
981 template <
int dim,
int spacedim>
1001 const Point<3> & candidate_point)
const
1003 return do_get_new_point(directions, distances, weights, candidate_point);
1014 const Point<3> & candidate_point)
const
1016 return do_get_new_point(directions, distances, weights, candidate_point);
1027 const Point<3> & candidate_point)
const
1029 return do_get_new_point(directions, distances, weights, candidate_point);
1037 template <
int dim,
int spacedim>
1039 const double tolerance)
1047 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1052 template <
int dim,
int spacedim>
1056 const double tolerance)
1059 , direction(direction / direction.norm())
1060 , point_on_axis(point_on_axis)
1061 , tolerance(tolerance)
1066 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1071 template <
int dim,
int spacedim>
1072 std::unique_ptr<Manifold<dim, spacedim>>
1075 return std_cxx14::make_unique<CylindricalManifold<dim, spacedim>>(
1076 direction, point_on_axis, tolerance);
1081 template <
int dim,
int spacedim>
1088 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1092 double average_length = 0.;
1093 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
1095 middle += surrounding_points[i] * weights[i];
1096 average_length += surrounding_points[i].square() * weights[i];
1098 middle -= point_on_axis;
1099 const double lambda = middle * direction;
1101 if ((middle - direction * lambda).square() < tolerance * average_length)
1102 return point_on_axis + direction * lambda;
1110 template <
int dim,
int spacedim>
1116 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1120 const double lambda = normalized_point * direction;
1121 const Point<spacedim> projection = point_on_axis + direction * lambda;
1126 const double dot = normal_direction * p_diff;
1131 return Point<3>(p_diff.norm(), phi, lambda);
1136 template <
int dim,
int spacedim>
1142 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1145 const double sine_r = std::sin(chart_point(1)) * chart_point(0);
1146 const double cosine_r = std::cos(chart_point(1)) * chart_point(0);
1149 normal_direction * cosine_r + dxn * sine_r;
1152 return point_on_axis + direction * chart_point(2) + intermediate;
1157 template <
int dim,
int spacedim>
1163 ExcMessage(
"CylindricalManifold can only be used for spacedim==3!"));
1168 const double sine = std::sin(chart_point(1));
1169 const double cosine = std::cos(chart_point(1));
1172 normal_direction * cosine + dxn * sine;
1175 constexpr
int s0 = 0 % spacedim;
1176 constexpr
int s1 = 1 % spacedim;
1177 constexpr
int s2 = 2 % spacedim;
1180 derivatives[s0][s0] = intermediate[s0];
1181 derivatives[s1][s0] = intermediate[s1];
1182 derivatives[s2][s0] = intermediate[s2];
1185 derivatives[s0][s1] = -normal_direction[s0] * sine + dxn[s0] * cosine;
1186 derivatives[s1][s1] = -normal_direction[s1] * sine + dxn[s1] * cosine;
1187 derivatives[s2][s1] = -normal_direction[s2] * sine + dxn[s2] * cosine;
1190 derivatives[s0][s2] = direction[s0];
1191 derivatives[s1][s2] = direction[s1];
1192 derivatives[s2][s2] = direction[s2];
1202 template <
int dim,
int spacedim>
1206 const double eccentricity)
1209 , direction(major_axis_direction)
1211 , cosh_u(1.0 / eccentricity)
1212 , sinh_u(std::
sqrt(cosh_u * cosh_u - 1.0))
1219 "Invalid eccentricity: It must satisfy 0 < eccentricity < 1."));
1221 Assert(direction_norm != 0.0,
1223 "Invalid major axis direction vector: Null vector not allowed."));
1229 template <
int dim,
int spacedim>
1230 std::unique_ptr<Manifold<dim, spacedim>>
1233 const double eccentricity = 1.0 / cosh_u;
1234 return std_cxx14::make_unique<EllipticalManifold<dim, spacedim>>(
1235 center, direction, eccentricity);
1240 template <
int dim,
int spacedim>
1253 template <
int dim,
int spacedim>
1267 const double cs = std::cos(chart_point[1]);
1268 const double sn = std::sin(chart_point[1]);
1271 const double x = chart_point[0] * cosh_u * cs;
1272 const double y = chart_point[0] * sinh_u * sn;
1274 const Point<2> p(direction[0] * x - direction[1] * y,
1275 direction[1] * x + direction[0] * y);
1281 template <
int dim,
int spacedim>
1296 const double x0 = space_point[0] -
center[0];
1297 const double y0 = space_point[1] -
center[1];
1298 const double x = direction[0] * x0 + direction[1] * y0;
1299 const double y = -direction[1] * x0 + direction[0] * y0;
1301 std::sqrt((x * x) / (cosh_u * cosh_u) + (y * y) / (sinh_u * sinh_u));
1307 double cos_eta = x / (pt0 * cosh_u);
1317 const double pt1 = (std::signbit(y) ? 2.0 *
numbers::PI - eta : eta);
1323 template <
int dim,
int spacedim>
1339 const double cs = std::cos(chart_point[1]);
1340 const double sn = std::sin(chart_point[1]);
1342 dX[0][0] = cosh_u * cs;
1343 dX[0][1] = -chart_point[0] * cosh_u * sn;
1344 dX[1][0] = sinh_u * sn;
1345 dX[1][1] = chart_point[0] * sinh_u * cs;
1349 {{+direction[0], -direction[1]}, {direction[1], direction[0]}}};
1359 template <
int dim,
int spacedim,
int chartdim>
1364 const double tolerance)
1367 , push_forward_function(&push_forward_function)
1368 , pull_back_function(&pull_back_function)
1369 , tolerance(tolerance)
1370 , owns_pointers(false)
1371 , finite_difference_step(0)
1379 template <
int dim,
int spacedim,
int chartdim>
1384 const double tolerance)
1388 , pull_back_function(
pull_back.release())
1389 , tolerance(tolerance)
1390 , owns_pointers(true)
1391 , finite_difference_step(0)
1399 template <
int dim,
int spacedim,
int chartdim>
1401 const std::string push_forward_expression,
1402 const std::string pull_back_expression,
1405 const std::string chart_vars,
1406 const std::string space_vars,
1407 const double tolerance,
1410 , const_map(const_map)
1411 , tolerance(tolerance)
1412 , owns_pointers(true)
1413 , push_forward_expression(push_forward_expression)
1414 , pull_back_expression(pull_back_expression)
1415 , chart_vars(chart_vars)
1416 , space_vars(space_vars)
1417 , finite_difference_step(h)
1429 template <
int dim,
int spacedim,
int chartdim>
1432 if (owns_pointers ==
true)
1435 push_forward_function =
nullptr;
1439 pull_back_function =
nullptr;
1446 template <
int dim,
int spacedim,
int chartdim>
1447 std::unique_ptr<Manifold<dim, spacedim>>
1462 if (!(push_forward_expression.empty() && pull_back_expression.empty()))
1464 return std_cxx14::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1465 push_forward_expression,
1466 pull_back_expression,
1467 this->get_periodicity(),
1472 finite_difference_step);
1476 return std_cxx14::make_unique<FunctionManifold<dim, spacedim, chartdim>>(
1477 *push_forward_function,
1478 *pull_back_function,
1479 this->get_periodicity(),
1486 template <
int dim,
int spacedim,
int chartdim>
1493 push_forward_function->vector_value(chart_point, pf);
1494 for (
unsigned int i = 0; i < spacedim; ++i)
1499 pull_back_function->vector_value(result, pb);
1500 for (
unsigned int i = 0; i < chartdim; ++i)
1502 (chart_point.
norm() > tolerance &&
1503 (std::abs(pb[i] - chart_point[i]) < tolerance * chart_point.
norm())) ||
1504 (std::abs(pb[i] - chart_point[i]) < tolerance),
1506 "The push forward is not the inverse of the pull back! Bailing out."));
1514 template <
int dim,
int spacedim,
int chartdim>
1520 for (
unsigned int i = 0; i < spacedim; ++i)
1522 const auto gradient = push_forward_function->gradient(chart_point, i);
1523 for (
unsigned int j = 0; j < chartdim; ++j)
1524 DF[i][j] = gradient[j];
1531 template <
int dim,
int spacedim,
int chartdim>
1538 pull_back_function->vector_value(space_point, pb);
1539 for (
unsigned int i = 0; i < chartdim; ++i)
1557 double theta =
std::atan2(z, std::sqrt(x * x + y * y) - R);
1558 double w = std::sqrt(std::pow(y - std::sin(phi) * R, 2.0) +
1559 std::pow(x - std::cos(phi) * R, 2.0) + z * z) /
1561 return {phi, theta,
w};
1570 double phi = chart_point(0);
1571 double theta = chart_point(1);
1572 double w = chart_point(2);
1574 return {std::cos(phi) * R + r *
w * std::cos(theta) * std::cos(phi),
1575 r *
w * std::sin(theta),
1576 std::sin(phi) * R + r *
w * std::cos(theta) * std::sin(phi)};
1588 ExcMessage(
"Outer radius R must be greater than the inner "
1596 std::unique_ptr<Manifold<dim, 3>>
1599 return std_cxx14::make_unique<TorusManifold<dim>>(R, r);
1610 double phi = chart_point(0);
1611 double theta = chart_point(1);
1612 double w = chart_point(2);
1614 DX[0][0] = -std::sin(phi) * R - r *
w * std::cos(theta) * std::sin(phi);
1615 DX[0][1] = -r *
w * std::sin(theta) * std::cos(phi);
1616 DX[0][2] = r * std::cos(theta) * std::cos(phi);
1619 DX[1][1] = r *
w * std::cos(theta);
1620 DX[1][2] = r * std::sin(theta);
1622 DX[2][0] = std::cos(phi) * R + r *
w * std::cos(theta) * std::cos(phi);
1623 DX[2][1] = -r *
w * std::sin(theta) * std::sin(phi);
1624 DX[2][2] = r * std::cos(theta) * std::sin(phi);
1634 template <
int dim,
int spacedim>
1645 template <
int dim,
int spacedim>
1647 spacedim>::~TransfiniteInterpolationManifold()
1649 if (clear_signal.connected())
1650 clear_signal.disconnect();
1655 template <
int dim,
int spacedim>
1656 std::unique_ptr<Manifold<dim, spacedim>>
1662 return std::unique_ptr<Manifold<dim, spacedim>>(ptr);
1667 template <
int dim,
int spacedim>
1674 clear_signal =
triangulation.signals.clear.connect([&]() ->
void {
1675 this->triangulation =
nullptr;
1676 this->level_coarse = -1;
1679 coarse_cell_is_flat.resize(
triangulation.n_cells(level_coarse),
false);
1683 for (; cell != endc; ++cell)
1685 bool cell_is_flat =
true;
1686 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1687 if (cell->line(
l)->manifold_id() != cell->manifold_id() &&
1689 cell_is_flat =
false;
1691 for (
unsigned int q = 0; q < GeometryInfo<dim>::quads_per_cell; ++q)
1692 if (cell->quad(q)->manifold_id() != cell->manifold_id() &&
1694 cell_is_flat =
false;
1696 coarse_cell_is_flat.size());
1697 coarse_cell_is_flat[cell->index()] = cell_is_flat;
1706 template <
typename AccessorType>
1708 compute_transfinite_interpolation(
const AccessorType &cell,
1712 return cell.vertex(0) * (1. - chart_point[0]) +
1713 cell.vertex(1) * chart_point[0];
1717 template <
typename AccessorType>
1719 compute_transfinite_interpolation(
const AccessorType &cell,
1721 const bool cell_is_flat)
1723 const unsigned int dim = AccessorType::dimension;
1724 const unsigned int spacedim = AccessorType::space_dimension;
1732 const std::array<Point<spacedim>, 4>
vertices{
1733 {cell.vertex(0), cell.vertex(1), cell.vertex(2), cell.vertex(3)}};
1738 std::array<double, 4> weights_vertices{
1739 {(1. - chart_point[0]) * (1. - chart_point[1]),
1740 chart_point[0] * (1. - chart_point[1]),
1741 (1. - chart_point[0]) * chart_point[1],
1742 chart_point[0] * chart_point[1]}};
1747 new_point += weights_vertices[v] *
vertices[v];
1759 std::array<double, GeometryInfo<2>::vertices_per_face> weights;
1762 const auto weights_view =
1764 const auto points_view =
make_array_view(points.begin(), points.end());
1766 for (
unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1769 const double my_weight =
1770 (line % 2) ? chart_point[line / 2] : 1 - chart_point[line / 2];
1771 const double line_point = chart_point[1 - line / 2];
1778 cell.line(line)->manifold_id();
1779 if (line_manifold_id == my_manifold_id ||
1784 my_weight * (1. - line_point);
1787 my_weight * line_point;
1795 weights[0] = 1. - line_point;
1796 weights[1] = line_point;
1799 .get_new_point(points_view, weights_view);
1805 new_point -= weights_vertices[v] *
vertices[v];
1814 static constexpr
unsigned int face_to_cell_vertices_3d[6][4] = {{0, 2, 4, 6},
1824 static constexpr
unsigned int face_to_cell_lines_3d[6][4] = {{8, 10, 0, 4},
1832 template <
typename AccessorType>
1834 compute_transfinite_interpolation(
const AccessorType &cell,
1836 const bool cell_is_flat)
1838 const unsigned int dim = AccessorType::dimension;
1839 const unsigned int spacedim = AccessorType::space_dimension;
1845 const std::array<Point<spacedim>, 8>
vertices{{cell.vertex(0),
1857 double linear_shapes[10];
1858 for (
unsigned int d = 0;
d < 3; ++
d)
1860 linear_shapes[2 *
d] = 1. - chart_point[
d];
1861 linear_shapes[2 *
d + 1] = chart_point[
d];
1865 for (
unsigned int d = 6;
d < 10; ++
d)
1866 linear_shapes[
d] = linear_shapes[
d - 6];
1868 std::array<double, 8> weights_vertices;
1869 for (
unsigned int i2 = 0, v = 0; i2 < 2; ++i2)
1870 for (
unsigned int i1 = 0; i1 < 2; ++i1)
1871 for (
unsigned int i0 = 0; i0 < 2; ++i0, ++v)
1872 weights_vertices[v] =
1873 (linear_shapes[4 + i2] * linear_shapes[2 + i1]) * linear_shapes[i0];
1877 for (
unsigned int v = 0; v < 8; ++v)
1878 new_point += weights_vertices[v] *
vertices[v];
1884 std::array<double, GeometryInfo<3>::lines_per_cell> weights_lines;
1885 std::fill(weights_lines.begin(), weights_lines.end(), 0.0);
1888 std::array<double, GeometryInfo<2>::vertices_per_cell> weights;
1891 const auto weights_view =
1893 const auto points_view =
make_array_view(points.begin(), points.end());
1897 const double my_weight = linear_shapes[face];
1898 const unsigned int face_even = face - face % 2;
1900 if (std::abs(my_weight) < 1
e-13)
1906 cell.face(face)->manifold_id();
1907 if (face_manifold_id == my_manifold_id ||
1910 for (
unsigned int line = 0;
1911 line < GeometryInfo<2>::lines_per_cell;
1914 const double line_weight =
1915 linear_shapes[face_even + 2 + line];
1916 weights_lines[face_to_cell_lines_3d[face][line]] +=
1917 my_weight * line_weight;
1923 weights_vertices[face_to_cell_vertices_3d[face][0]] -=
1924 linear_shapes[face_even + 2] *
1925 (linear_shapes[face_even + 4] * my_weight);
1926 weights_vertices[face_to_cell_vertices_3d[face][1]] -=
1927 linear_shapes[face_even + 3] *
1928 (linear_shapes[face_even + 4] * my_weight);
1929 weights_vertices[face_to_cell_vertices_3d[face][2]] -=
1930 linear_shapes[face_even + 2] *
1931 (linear_shapes[face_even + 5] * my_weight);
1932 weights_vertices[face_to_cell_vertices_3d[face][3]] -=
1933 linear_shapes[face_even + 3] *
1934 (linear_shapes[face_even + 5] * my_weight);
1939 points[v] =
vertices[face_to_cell_vertices_3d[face][v]];
1941 linear_shapes[face_even + 2] * linear_shapes[face_even + 4];
1943 linear_shapes[face_even + 3] * linear_shapes[face_even + 4];
1945 linear_shapes[face_even + 2] * linear_shapes[face_even + 5];
1947 linear_shapes[face_even + 3] * linear_shapes[face_even + 5];
1950 .get_new_point(points_view, weights_view);
1955 const auto weights_view_line =
1957 const auto points_view_line =
1959 for (
unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
1962 const double line_point =
1963 (line < 8 ? chart_point[1 - (line % 4) / 2] : chart_point[2]);
1964 double my_weight = 0.;
1966 my_weight = linear_shapes[line % 4] * linear_shapes[4 + line / 4];
1969 const unsigned int subline = line - 8;
1971 linear_shapes[subline % 2] * linear_shapes[2 + subline / 2];
1973 my_weight -= weights_lines[line];
1975 if (std::abs(my_weight) < 1
e-13)
1979 cell.line(line)->manifold_id();
1980 if (line_manifold_id == my_manifold_id ||
1985 my_weight * (1. - line_point);
1988 my_weight * (line_point);
1996 weights[0] = 1. - line_point;
1997 weights[1] = line_point;
1998 new_point -= my_weight * tria.
get_manifold(line_manifold_id)
1999 .get_new_point(points_view_line,
2006 new_point += weights_vertices[v] *
vertices[v];
2014 template <
int dim,
int spacedim>
2026 ExcMessage(
"chart_point is not in unit interval"));
2028 return compute_transfinite_interpolation(*cell,
2030 coarse_cell_is_flat[cell->index()]);
2035 template <
int dim,
int spacedim>
2044 for (
unsigned int d = 0;
d < dim; ++
d)
2047 const double step = chart_point[
d] > 0.5 ? -1
e-8 : 1
e-8;
2050 modified[
d] += step;
2052 compute_transfinite_interpolation(*cell,
2054 coarse_cell_is_flat[cell->index()]) -
2055 pushed_forward_chart_point;
2056 for (
unsigned int e = 0;
e < spacedim; ++
e)
2057 grad[
e][
d] = difference[
e] / step;
2064 template <
int dim,
int spacedim>
2072 for (
unsigned int d = 0;
d < dim; ++
d)
2086 compute_transfinite_interpolation(*cell,
2088 coarse_cell_is_flat[cell->index()]);
2089 const double tolerance = 1
e-21 * Utilities::fixed_power<2>(cell->diameter());
2090 double residual_norm_square = residual.
norm_square();
2092 bool must_recompute_jacobian =
true;
2093 for (
unsigned int i = 0; i < 100; ++i)
2095 if (residual_norm_square < tolerance)
2102 for (
unsigned int d = 0;
d < spacedim; ++
d)
2103 for (
unsigned int e = 0;
e < dim; ++
e)
2104 update[
e] += inv_grad[
d][
e] * residual[
d];
2105 return chart_point + update;
2117 if (must_recompute_jacobian || i % 9 == 0)
2125 push_forward_gradient(cell,
2131 must_recompute_jacobian =
false;
2134 for (
unsigned int d = 0;
d < spacedim; ++
d)
2135 for (
unsigned int e = 0;
e < dim; ++
e)
2136 update[
e] += inv_grad[
d][
e] * residual[
d];
2150 while (alpha > 1
e-4)
2152 Point<dim> guess = chart_point + alpha * update;
2154 point - compute_transfinite_interpolation(
2155 *cell, guess, coarse_cell_is_flat[cell->index()]);
2156 const double residual_norm_new = residual.
norm_square();
2157 if (residual_norm_new < residual_norm_square)
2159 residual_norm_square = residual_norm_new;
2160 chart_point += alpha * update;
2173 must_recompute_jacobian =
true;
2187 for (
unsigned int d = 0;
d < spacedim; ++
d)
2188 for (
unsigned int e = 0;
e < dim; ++
e)
2189 Jinv_deltaf[
e] += inv_grad[
d][
e] * delta_f[
d];
2196 if (std::abs(delta_x * Jinv_deltaf) > 1
e-12)
2199 (delta_x - Jinv_deltaf) / (delta_x * Jinv_deltaf);
2201 for (
unsigned int d = 0;
d < spacedim; ++
d)
2202 for (
unsigned int e = 0;
e < dim; ++
e)
2203 jac_update[
d] += delta_x[
e] * inv_grad[
d][
e];
2204 for (
unsigned int d = 0;
d < spacedim; ++
d)
2205 for (
unsigned int e = 0;
e < dim; ++
e)
2206 inv_grad[
d][
e] += factor[
e] * jac_update[
d];
2214 template <
int dim,
int spacedim>
2215 std::array<unsigned int, 20>
2225 ExcMessage(
"The manifold was initialized with level " +
2227 "active cells on a lower level. Coarsening the mesh is " +
2228 "currently not supported"));
2238 boost::container::small_vector<std::pair<double, unsigned int>, 200>
2239 distances_and_cells;
2240 for (; cell != endc; ++cell)
2243 if (&cell->get_manifold() !=
this)
2250 vertices[vertex_n] = cell->vertex(vertex_n);
2260 double radius_square = 0.;
2264 bool inside_circle =
true;
2265 for (
unsigned int i = 0; i < points.size(); ++i)
2266 if ((
center - points[i]).norm_square() > radius_square * 1.5)
2268 inside_circle =
false;
2271 if (inside_circle ==
false)
2275 double current_distance = 0;
2276 for (
unsigned int i = 0; i < points.size(); ++i)
2279 cell->real_to_unit_cell_affine_approximation(points[i]);
2282 distances_and_cells.push_back(
2283 std::make_pair(current_distance, cell->index()));
2288 std::sort(distances_and_cells.begin(), distances_and_cells.end());
2289 std::array<unsigned int, 20> cells;
2291 for (
unsigned int i = 0; i < distances_and_cells.size() && i < cells.size();
2293 cells[i] = distances_and_cells[i].second;
2300 template <
int dim,
int spacedim>
2306 Assert(surrounding_points.size() == chart_points.size(),
2307 ExcMessage(
"The chart points array view must be as large as the "
2308 "surrounding points array view."));
2310 std::array<unsigned int, 20> nearby_cells =
2311 get_possible_cells_around_points(surrounding_points);
2335 auto guess_chart_point_structdim_2 = [&](
const unsigned int i) ->
Point<dim> {
2336 Assert(surrounding_points.size() == 8 && 2 < i && i < 8,
2337 ExcMessage(
"This function assumes that there are eight surrounding "
2338 "points around a two-dimensional object. It also assumes "
2339 "that the first three chart points have already been "
2349 return chart_points[1] + (chart_points[2] - chart_points[0]);
2351 return 0.5 * (chart_points[0] + chart_points[2]);
2353 return 0.5 * (chart_points[1] + chart_points[3]);
2355 return 0.5 * (chart_points[0] + chart_points[1]);
2357 return 0.5 * (chart_points[2] + chart_points[3]);
2386 auto guess_chart_point_structdim_3 = [&](
const unsigned int i) ->
Point<dim> {
2387 Assert(surrounding_points.size() == 8 && 4 < i && i < 8,
2388 ExcMessage(
"This function assumes that there are eight surrounding "
2389 "points around a three-dimensional object. It also "
2390 "assumes that the first five chart points have already "
2392 return chart_points[i - 4] + (chart_points[4] - chart_points[0]);
2396 bool use_structdim_2_guesses =
false;
2397 bool use_structdim_3_guesses =
false;
2402 if (surrounding_points.size() == 8)
2405 surrounding_points[6] - surrounding_points[0];
2407 surrounding_points[7] - surrounding_points[2];
2415 use_structdim_2_guesses =
true;
2416 else if (spacedim == 3)
2419 use_structdim_3_guesses =
true;
2422 Assert((!use_structdim_2_guesses && !use_structdim_3_guesses) ||
2423 (use_structdim_2_guesses ^ use_structdim_3_guesses),
2428 auto compute_chart_point =
2430 const unsigned int point_index) {
2436 bool used_affine_approximation =
false;
2439 if (point_index == 3 && surrounding_points.size() >= 8)
2440 guess = chart_points[1] + (chart_points[2] - chart_points[0]);
2441 else if (use_structdim_2_guesses && 3 < point_index)
2442 guess = guess_chart_point_structdim_2(point_index);
2443 else if (use_structdim_3_guesses && 4 < point_index)
2444 guess = guess_chart_point_structdim_3(point_index);
2445 else if (dim == 3 && point_index > 7 && surrounding_points.size() == 26)
2447 if (point_index < 20)
2450 point_index - 8, 0)] +
2452 point_index - 8, 1)]);
2456 point_index - 20, 0)] +
2458 point_index - 20, 1)] +
2460 point_index - 20, 2)] +
2462 point_index - 20, 3)]);
2466 guess = cell->real_to_unit_cell_affine_approximation(
2467 surrounding_points[point_index]);
2468 used_affine_approximation =
true;
2470 chart_points[point_index] =
2471 pull_back(cell, surrounding_points[point_index], guess);
2476 if (chart_points[point_index][0] ==
2478 !used_affine_approximation)
2480 guess = cell->real_to_unit_cell_affine_approximation(
2481 surrounding_points[point_index]);
2482 chart_points[point_index] =
2483 pull_back(cell, surrounding_points[point_index], guess);
2486 if (chart_points[point_index][0] ==
2489 for (
unsigned int d = 0;
d < dim; ++
d)
2491 chart_points[point_index] =
2492 pull_back(cell, surrounding_points[point_index], guess);
2497 for (
unsigned int c = 0; c < nearby_cells.size(); ++c)
2501 bool inside_unit_cell =
true;
2502 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2504 compute_chart_point(cell, i);
2511 inside_unit_cell =
false;
2515 if (inside_unit_cell ==
true)
2522 if (c == nearby_cells.size() - 1 ||
2527 std::ostringstream message;
2528 for (
unsigned int b = 0;
b <= c; ++
b)
2532 message <<
"Looking at cell " << cell->id()
2533 <<
" with vertices: " << std::endl;
2535 message << cell->vertex(v) <<
" ";
2536 message << std::endl;
2537 message <<
"Transformation to chart coordinates: " << std::endl;
2538 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
2540 compute_chart_point(cell, i);
2541 message << surrounding_points[i] <<
" -> " << chart_points[i]
2561 template <
int dim,
int spacedim>
2567 boost::container::small_vector<Point<dim>, 100> chart_points(
2568 surrounding_points.size());
2571 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2574 chart_manifold.get_new_point(chart_points_view, weights);
2581 template <
int dim,
int spacedim>
2591 boost::container::small_vector<Point<dim>, 100> chart_points(
2592 surrounding_points.size());
2595 const auto cell = compute_chart_points(surrounding_points, chart_points_view);
2597 boost::container::small_vector<Point<dim>, 100> new_points_on_chart(
2599 chart_manifold.get_new_points(chart_points_view,
2602 new_points_on_chart.end()));
2604 for (
unsigned int row = 0; row < weights.
size(0); ++row)
2605 new_points[row] =
push_forward(cell, new_points_on_chart[row]);
2611 #include "manifold_lib.inst"