27 #include <boost/container/small_vector.hpp>
36 template <
int dim,
int spacedim>
48 template <
int dim,
int spacedim>
54 const std::array<Point<spacedim>, 2>
vertices{{p1, p2}};
56 w * p2 + (1 -
w) * p1);
61 template <
int dim,
int spacedim>
67 const double tol = 1
e-10;
68 const unsigned int n_points = surrounding_points.size();
74 "There should be as many surrounding points as weights given."));
76 Assert(std::abs(std::accumulate(weights.
begin(), weights.
end(), 0.0) - 1.0) <
78 ExcMessage(
"The weights for the individual points should sum to 1!"));
83 boost::container::small_vector<unsigned int, 100> permutation(n_points);
84 std::iota(permutation.begin(), permutation.end(), 0u);
85 std::sort(permutation.begin(),
87 [&weights](
const std::size_t a,
const std::size_t
b) {
88 return weights[a] < weights[b];
93 double w = weights[permutation[0]];
95 for (
unsigned int i = 1; i < n_points; ++i)
98 if (std::abs(weights[permutation[i]] +
w) < tol)
101 weight =
w / (weights[permutation[i]] +
w);
103 if (std::abs(weight) > 1
e-14)
105 p = get_intermediate_point(p,
106 surrounding_points[permutation[i]],
111 p = surrounding_points[permutation[i]];
113 w += weights[permutation[i]];
121 template <
int dim,
int spacedim>
130 for (
unsigned int row = 0; row < weights.
size(0); ++row)
134 surrounding_points.end()),
146 const int spacedim = 2;
151 ((p - face->vertex(0)).norm_square() > (p - face->vertex(1)).norm_square() ?
152 -get_tangent_vector(p, face->vertex(0)) :
153 get_tangent_vector(p, face->vertex(1)));
157 return normal / normal.
norm();
167 const int spacedim = 3;
169 const std::array<Point<spacedim>, 4>
vertices{
170 {face->vertex(0), face->vertex(1), face->vertex(2), face->vertex(3)}};
176 std::max(distances[2], distances[3]));
187 for (
unsigned int i = 0; i < 3; ++i)
188 if (distances[i] > 1
e-8 * max_distance)
189 for (
unsigned int j = i + 1; j < 4; ++j)
190 if (distances[j] > 1
e-8 * max_distance)
193 (distances[i] * distances[j]);
196 if (std::abs(new_angle) < 0.999 * abs_cos_angle)
198 abs_cos_angle = std::abs(new_angle);
204 ExcMessage(
"The search for possible directions did not succeed."));
215 "Manifold::normal_vector was unable to find a suitable combination "
216 "of vertices to compute a normal on this face. We chose the secants "
217 "that are as orthogonal as possible, but tangents appear to be "
218 "linearly dependent. Check for distorted faces in your triangulation."));
227 if (reference_normal * normal < 0.0)
230 return normal / normal.
norm();
235 template <
int dim,
int spacedim>
253 n[0] =
cross_product_2d(get_tangent_vector(face->vertex(0), face->vertex(1)));
257 for (
unsigned int i = 0; i < 2; ++i)
261 "computed normals have "
275 n[0] =
cross_product_3d(get_tangent_vector(face->vertex(0), face->vertex(1)),
276 get_tangent_vector(face->vertex(0), face->vertex(2)));
278 n[1] =
cross_product_3d(get_tangent_vector(face->vertex(1), face->vertex(3)),
279 get_tangent_vector(face->vertex(1), face->vertex(0)));
281 n[2] =
cross_product_3d(get_tangent_vector(face->vertex(2), face->vertex(0)),
282 get_tangent_vector(face->vertex(2), face->vertex(3)));
284 n[3] =
cross_product_3d(get_tangent_vector(face->vertex(3), face->vertex(2)),
285 get_tangent_vector(face->vertex(3), face->vertex(1)));
287 for (
unsigned int i = 0; i < 4; ++i)
291 "computed normals have "
299 template <
int dim,
int spacedim>
305 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
307 n[v] = normal_vector(face, face->vertex(v));
314 template <
int dim,
int spacedim>
321 points_weights.first.end()),
323 points_weights.second.end()));
328 template <
int dim,
int spacedim>
335 points_weights.first.end()),
337 points_weights.second.end()));
342 template <
int dim,
int spacedim>
352 return get_new_point_on_line(face);
354 return get_new_point_on_quad(face);
362 template <
int dim,
int spacedim>
370 return get_new_point_on_line(cell);
372 return get_new_point_on_quad(cell);
374 return get_new_point_on_hex(cell);
448 template <
int dim,
int spacedim>
466 points_weights.first.end()),
468 points_weights.second.end()));
473 template <
int dim,
int spacedim>
480 const std::array<Point<spacedim>, 2> points{{x1, x2}};
485 return (neighbor_point - x1) /
epsilon;
495 normalized_alternating_product(
const Tensor<1, 3> (&)[1])
507 normalized_alternating_product(
const Tensor<1, 3> (&basis_vectors)[2])
510 return tmp / tmp.
norm();
516 template <
int dim,
int spacedim>
519 const double tolerance)
520 : periodicity(periodicity)
521 , tolerance(tolerance)
526 template <
int dim,
int spacedim>
527 std::unique_ptr<Manifold<dim, spacedim>>
530 return std_cxx14::make_unique<FlatManifold<dim, spacedim>>(periodicity,
536 template <
int dim,
int spacedim>
542 Assert(std::abs(std::accumulate(weights.
begin(), weights.
end(), 0.0) - 1.0) <
544 ExcMessage(
"The weights for the individual points should sum to 1!"));
551 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
552 p += surrounding_points[i] * weights[i];
558 for (
unsigned int d = 0;
d < spacedim; ++
d)
559 if (periodicity[
d] > 0)
560 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
562 minP[
d] =
std::min(minP[
d], surrounding_points[i][
d]);
563 Assert((surrounding_points[i][
d] <
564 periodicity[
d] + tolerance * periodicity[
d]) ||
565 (surrounding_points[i][
d] >=
566 -tolerance * periodicity[
d]),
567 ExcPeriodicBox(
d, surrounding_points[i], periodicity[
d]));
572 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
575 for (
unsigned int d = 0;
d < spacedim; ++
d)
576 if (periodicity[
d] > 0)
578 ((surrounding_points[i][
d] - minP[
d]) > periodicity[
d] / 2.0 ?
582 p += (surrounding_points[i] + dp) * weights[i];
586 for (
unsigned int d = 0;
d < spacedim; ++
d)
587 if (periodicity[
d] > 0)
589 p[
d] += periodicity[
d];
592 return project_to_manifold(surrounding_points, p);
597 template <
int dim,
int spacedim>
605 if (weights.
size(0) == 0)
608 const std::size_t n_points = surrounding_points.size();
611 for (
unsigned int d = 0;
d < spacedim; ++
d)
612 if (periodicity[
d] > 0)
613 for (
unsigned int i = 0; i < n_points; ++i)
615 minP[
d] =
std::min(minP[
d], surrounding_points[i][
d]);
616 Assert((surrounding_points[i][
d] <
617 periodicity[
d] + tolerance * periodicity[
d]) ||
618 (surrounding_points[i][
d] >= -tolerance * periodicity[
d]),
619 ExcPeriodicBox(
d, surrounding_points[i], periodicity[i]));
624 const Point<spacedim> *surrounding_points_start = surrounding_points.data();
626 boost::container::small_vector<Point<spacedim>, 200> modified_points;
627 bool adjust_periodicity =
false;
628 for (
unsigned int d = 0;
d < spacedim; ++
d)
629 if (periodicity[
d] > 0)
630 for (
unsigned int i = 0; i < n_points; ++i)
631 if ((surrounding_points[i][
d] - minP[
d]) > periodicity[
d] / 2.0)
633 adjust_periodicity =
true;
636 if (adjust_periodicity ==
true)
638 modified_points.resize(surrounding_points.size());
640 surrounding_points.end(),
641 modified_points.begin());
642 for (
unsigned int d = 0;
d < spacedim; ++
d)
643 if (periodicity[
d] > 0)
644 for (
unsigned int i = 0; i < n_points; ++i)
645 if ((surrounding_points[i][
d] - minP[
d]) > periodicity[
d] / 2.0)
646 modified_points[i][
d] -= periodicity[
d];
647 surrounding_points_start = modified_points.data();
651 for (
unsigned int row = 0; row < weights.
size(0); ++row)
655 std::accumulate(&weights(row, 0), &weights(row, 0) + n_points, 0.0) -
657 ExcMessage(
"The weights for the individual points should sum to 1!"));
659 for (
unsigned int p = 0; p < n_points; ++p)
660 new_point += surrounding_points_start[p] * weights(row, p);
663 for (
unsigned int d = 0;
d < spacedim; ++
d)
664 if (periodicity[
d] > 0)
665 if (new_point[
d] < 0)
666 new_point[
d] += periodicity[
d];
672 surrounding_points.end()),
679 template <
int dim,
int spacedim>
690 template <
int dim,
int spacedim>
699 template <
int dim,
int spacedim>
710 for (
unsigned int d = 0;
d < spacedim; ++
d)
711 if (periodicity[
d] > tolerance)
713 if (direction[
d] < -periodicity[
d] / 2)
714 direction[
d] += periodicity[
d];
715 else if (direction[
d] > periodicity[
d] / 2)
716 direction[
d] -= periodicity[
d];
763 const Tensor<1, 2> tangent = face->vertex(1) - face->vertex(0);
764 for (
unsigned int vertex = 0; vertex < GeometryInfo<2>::vertices_per_face;
767 face_vertex_normals[vertex] =
Point<2>(tangent[1], -tangent[0]);
791 static const unsigned int neighboring_vertices[4][2] = {{1, 2},
795 for (
unsigned int vertex = 0; vertex < vertices_per_face; ++vertex)
800 face->vertex(neighboring_vertices[vertex][0]) - face->vertex(vertex),
801 face->vertex(neighboring_vertices[vertex][1]) - face->vertex(vertex)};
857 template <
int dim,
int spacedim>
893 const unsigned int facedim = dim - 1;
896 for (
unsigned int i = 0; i < facedim; ++i)
899 const double eps = 1
e-12;
901 unsigned int iteration = 0;
906 F += face->vertex(v) *
909 for (
unsigned int i = 0; i < facedim; ++i)
919 for (
unsigned int i = 0; i < facedim; ++i)
920 for (
unsigned int j = 0; j < spacedim; ++j)
921 J[i] += grad_F[i][j] * (
F - p)[j];
924 for (
unsigned int i = 0; i < facedim; ++i)
925 for (
unsigned int j = 0; j < facedim; ++j)
926 for (
unsigned int k = 0; k < spacedim; ++k)
927 H[i][j] += grad_F[i][k] * grad_F[j][k];
934 ExcMessage(
"The Newton iteration to find the reference point "
935 "did not converge in 10 iterations. Do you have a "
936 "deformed cell? (See the glossary for a definition "
937 "of what a deformed cell is. You may want to output "
938 "the vertices of your cell."));
947 const double normalized_delta_world = (
F - p).
norm() / face->diameter();
949 if (delta_xi.
norm() <
eps || normalized_delta_world <
eps)
957 return internal::normalized_alternating_product(grad_F);
962 template <
int dim,
int spacedim,
int chartdim>
965 : sub_manifold(periodicity)
970 template <
int dim,
int spacedim,
int chartdim>
975 const double w)
const
977 const std::array<Point<spacedim>, 2> points{{p1, p2}};
978 const std::array<double, 2> weights{{1. -
w,
w}};
985 template <
int dim,
int spacedim,
int chartdim>
991 const std::size_t n_points = surrounding_points.size();
993 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
995 for (
unsigned int i = 0; i < n_points; ++i)
996 chart_points[i] =
pull_back(surrounding_points[i]);
1006 template <
int dim,
int spacedim,
int chartdim>
1016 const std::size_t n_points = surrounding_points.size();
1018 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1019 for (std::size_t i = 0; i < n_points; ++i)
1020 chart_points[i] =
pull_back(surrounding_points[i]);
1022 boost::container::small_vector<Point<chartdim>, 200> new_points_on_chart(
1024 sub_manifold.get_new_points(
1027 make_array_view(new_points_on_chart.begin(), new_points_on_chart.end()));
1029 for (std::size_t row = 0; row < weights.
size(0); ++row)
1030 new_points[row] =
push_forward(new_points_on_chart[row]);
1035 template <
int dim,
int spacedim,
int chartdim>
1048 template <
int dim,
int spacedim,
int chartdim>
1064 1
e-12 * F_prime.
norm(),
1066 "The derivative of a chart function must not be singular."));
1072 for (
unsigned int i = 0; i < spacedim; ++i)
1073 result[i] += F_prime[i] * delta;
1080 template <
int dim,
int spacedim,
int chartdim>
1084 return sub_manifold.get_periodicity();
1088 #include "manifold.inst"