16 #include <deal.II/base/std_cxx14/memory.h> 17 #include <deal.II/base/table.h> 18 #include <deal.II/base/tensor.h> 20 #include <deal.II/fe/fe_q.h> 22 #include <deal.II/grid/manifold.h> 23 #include <deal.II/grid/tria.h> 24 #include <deal.II/grid/tria_accessor.h> 25 #include <deal.II/grid/tria_iterator.h> 27 #include <boost/container/small_vector.hpp> 31 DEAL_II_NAMESPACE_OPEN
36 template <
int dim,
int spacedim>
48 template <
int dim,
int spacedim>
54 const std::array<Point<spacedim>, 2> vertices{{p1, p2}};
55 return project_to_manifold(make_array_view(vertices.begin(), vertices.end()),
56 w * p2 + (1 - w) * p1);
61 template <
int dim,
int spacedim>
67 const double tol = 1e-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 ((weights[permutation[i]] + w) < tol)
101 weight = w / (weights[permutation[i]] + w);
103 if (std::abs(weight) > 1e-14)
104 p = get_intermediate_point(p,
105 surrounding_points[permutation[i]],
107 w += weights[permutation[i]];
115 template <
int dim,
int spacedim>
124 for (
unsigned int row = 0; row < weights.
size(0); ++row)
127 get_new_point(make_array_view(surrounding_points.begin(),
128 surrounding_points.end()),
129 make_array_view(weights, row));
140 const int spacedim = 2;
145 ((p - face->vertex(0)).norm_square() > (p - face->vertex(1)).norm_square() ?
146 -get_tangent_vector(p, face->vertex(0)) :
147 get_tangent_vector(p, face->vertex(1)));
151 return normal / normal.
norm();
161 const int spacedim = 3;
163 const std::array<Point<spacedim>, 4> vertices{
164 {face->vertex(0), face->vertex(1), face->vertex(2), face->vertex(3)}};
165 const std::array<double, 4> distances{{vertices[0].
distance(p),
169 const double max_distance = std::max(std::max(distances[0], distances[1]),
170 std::max(distances[2], distances[3]));
178 double abs_cos_angle = std::numeric_limits<double>::max();
181 for (
unsigned int i = 0; i < 3; ++i)
182 if (distances[i] > 1e-8 * max_distance)
183 for (
unsigned int j = i + 1; j < 4; ++j)
184 if (distances[j] > 1e-8 * max_distance)
186 const double new_angle = (p - vertices[i]) * (p - vertices[j]) /
187 (distances[i] * distances[j]);
190 if (std::abs(new_angle) < 0.999 * abs_cos_angle)
192 abs_cos_angle = std::abs(new_angle);
198 ExcMessage(
"The search for possible directions did not succeed."));
206 normal.
norm_square() > 1e4 * std::numeric_limits<double>::epsilon() *
209 "Manifold::normal_vector was unable to find a suitable combination " 210 "of vertices to compute a normal on this face. We chose the secants " 211 "that are as orthogonal as possible, but tangents appear to be " 212 "linearly dependent. Check for distorted faces in your triangulation."));
221 if (reference_normal * normal < 0.0)
224 return normal / normal.
norm();
229 template <
int dim,
int spacedim>
245 FaceVertexNormals & n)
const 247 n[0] = cross_product_2d(get_tangent_vector(face->vertex(0), face->vertex(1)));
249 -cross_product_2d(get_tangent_vector(face->vertex(1), face->vertex(0)));
251 for (
unsigned int i = 0; i < 2; ++i)
255 "computed normals have " 267 FaceVertexNormals & n)
const 269 n[0] = cross_product_3d(get_tangent_vector(face->vertex(0), face->vertex(1)),
270 get_tangent_vector(face->vertex(0), face->vertex(2)));
272 n[1] = cross_product_3d(get_tangent_vector(face->vertex(1), face->vertex(3)),
273 get_tangent_vector(face->vertex(1), face->vertex(0)));
275 n[2] = cross_product_3d(get_tangent_vector(face->vertex(2), face->vertex(0)),
276 get_tangent_vector(face->vertex(2), face->vertex(3)));
278 n[3] = cross_product_3d(get_tangent_vector(face->vertex(3), face->vertex(2)),
279 get_tangent_vector(face->vertex(3), face->vertex(1)));
281 for (
unsigned int i = 0; i < 4; ++i)
285 "computed normals have " 293 template <
int dim,
int spacedim>
299 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
301 n[v] = normal_vector(face, face->vertex(v));
308 template <
int dim,
int spacedim>
314 return get_new_point(make_array_view(points_weights.first.begin(),
315 points_weights.first.end()),
316 make_array_view(points_weights.second.begin(),
317 points_weights.second.end()));
322 template <
int dim,
int spacedim>
328 return get_new_point(make_array_view(points_weights.first.begin(),
329 points_weights.first.end()),
330 make_array_view(points_weights.second.begin(),
331 points_weights.second.end()));
336 template <
int dim,
int spacedim>
346 return get_new_point_on_line(face);
348 return get_new_point_on_quad(face);
356 template <
int dim,
int spacedim>
364 return get_new_point_on_line(cell);
366 return get_new_point_on_quad(cell);
368 return get_new_point_on_hex(cell);
442 template <
int dim,
int spacedim>
459 return get_new_point(make_array_view(points_weights.first.begin(),
460 points_weights.first.end()),
461 make_array_view(points_weights.second.begin(),
462 points_weights.second.end()));
467 template <
int dim,
int spacedim>
472 const double epsilon = 1e-8;
474 const std::array<Point<spacedim>, 2> points{{x1, x2}};
475 const std::array<double, 2> weights{{epsilon, 1.0 - epsilon}};
477 get_new_point(make_array_view(points.begin(), points.end()),
478 make_array_view(weights.begin(), weights.end()));
479 return (neighbor_point - x1) / epsilon;
489 normalized_alternating_product(
const Tensor<1, 3> (&)[1])
501 normalized_alternating_product(
const Tensor<1, 3> (&basis_vectors)[2])
503 Tensor<1, 3> tmp = cross_product_3d(basis_vectors[0], basis_vectors[1]);
504 return tmp / tmp.
norm();
510 template <
int dim,
int spacedim>
513 const double tolerance)
514 : periodicity(periodicity)
515 , tolerance(tolerance)
520 template <
int dim,
int spacedim>
521 std::unique_ptr<Manifold<dim, spacedim>>
524 return std_cxx14::make_unique<FlatManifold<dim, spacedim>>(periodicity,
530 template <
int dim,
int spacedim>
536 Assert(std::abs(std::accumulate(weights.
begin(), weights.
end(), 0.0) - 1.0) <
538 ExcMessage(
"The weights for the individual points should sum to 1!"));
545 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
546 p += surrounding_points[i] * weights[i];
552 for (
unsigned int d = 0; d < spacedim; ++d)
553 if (periodicity[d] > 0)
554 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
556 minP[d] = std::min(minP[d], surrounding_points[i][d]);
557 Assert((surrounding_points[i][d] <
558 periodicity[d] + tolerance * periodicity[d]) ||
559 (surrounding_points[i][d] >=
560 -tolerance * periodicity[d]),
561 ExcPeriodicBox(d, surrounding_points[i], periodicity[d]));
566 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
569 for (
unsigned int d = 0; d < spacedim; ++d)
570 if (periodicity[d] > 0)
572 ((surrounding_points[i][d] - minP[d]) > periodicity[d] / 2.0 ?
576 p += (surrounding_points[i] + dp) * weights[i];
580 for (
unsigned int d = 0; d < spacedim; ++d)
581 if (periodicity[d] > 0)
583 p[d] += periodicity[d];
586 return project_to_manifold(surrounding_points, p);
591 template <
int dim,
int spacedim>
599 if (weights.
size(0) == 0)
602 const std::size_t n_points = surrounding_points.size();
605 for (
unsigned int d = 0; d < spacedim; ++d)
606 if (periodicity[d] > 0)
607 for (
unsigned int i = 0; i < n_points; ++i)
609 minP[d] = std::min(minP[d], surrounding_points[i][d]);
610 Assert((surrounding_points[i][d] <
611 periodicity[d] + tolerance * periodicity[d]) ||
612 (surrounding_points[i][d] >= -tolerance * periodicity[d]),
613 ExcPeriodicBox(d, surrounding_points[i], periodicity[i]));
618 const Point<spacedim> *surrounding_points_start = surrounding_points.data();
620 boost::container::small_vector<Point<spacedim>, 200> modified_points;
621 bool adjust_periodicity =
false;
622 for (
unsigned int d = 0; d < spacedim; ++d)
623 if (periodicity[d] > 0)
624 for (
unsigned int i = 0; i < n_points; ++i)
625 if ((surrounding_points[i][d] - minP[d]) > periodicity[d] / 2.0)
627 adjust_periodicity =
true;
630 if (adjust_periodicity ==
true)
632 modified_points.resize(surrounding_points.size());
633 std::copy(surrounding_points.begin(),
634 surrounding_points.end(),
635 modified_points.begin());
636 for (
unsigned int d = 0; d < spacedim; ++d)
637 if (periodicity[d] > 0)
638 for (
unsigned int i = 0; i < n_points; ++i)
639 if ((surrounding_points[i][d] - minP[d]) > periodicity[d] / 2.0)
640 modified_points[i][d] -= periodicity[d];
641 surrounding_points_start = modified_points.data();
645 for (
unsigned int row = 0; row < weights.
size(0); ++row)
649 std::accumulate(&weights(row, 0), &weights(row, 0) + n_points, 0.0) -
651 ExcMessage(
"The weights for the individual points should sum to 1!"));
653 for (
unsigned int p = 0; p < n_points; ++p)
654 new_point += surrounding_points_start[p] * weights(row, p);
657 for (
unsigned int d = 0; d < spacedim; ++d)
658 if (periodicity[d] > 0)
659 if (new_point[d] < 0)
660 new_point[d] += periodicity[d];
666 surrounding_points.end()),
673 template <
int dim,
int spacedim>
684 template <
int dim,
int spacedim>
693 template <
int dim,
int spacedim>
704 for (
unsigned int d = 0; d < spacedim; ++d)
705 if (periodicity[d] > tolerance)
707 if (direction[d] < -periodicity[d] / 2)
708 direction[d] += periodicity[d];
709 else if (direction[d] > periodicity[d] / 2)
710 direction[d] -= periodicity[d];
757 const Tensor<1, 2> tangent = face->vertex(1) - face->vertex(0);
758 for (
unsigned int vertex = 0; vertex < GeometryInfo<2>::vertices_per_face;
761 face_vertex_normals[vertex] =
Point<2>(tangent[1], -tangent[0]);
785 static const unsigned int neighboring_vertices[4][2] = {{1, 2},
789 for (
unsigned int vertex = 0; vertex < vertices_per_face; ++vertex)
794 face->vertex(neighboring_vertices[vertex][0]) - face->vertex(vertex),
795 face->vertex(neighboring_vertices[vertex][1]) - face->vertex(vertex)};
851 template <
int dim,
int spacedim>
887 const unsigned int facedim = dim - 1;
890 for (
unsigned int i = 0; i < facedim; ++i)
893 const double eps = 1e-12;
895 unsigned int iteration = 0;
899 for (
unsigned int v = 0; v < GeometryInfo<facedim>::vertices_per_cell;
901 F += face->vertex(v) *
904 for (
unsigned int i = 0; i < facedim; ++i)
907 for (
unsigned int v = 0; v < GeometryInfo<facedim>::vertices_per_cell;
915 for (
unsigned int i = 0; i < facedim; ++i)
916 for (
unsigned int j = 0; j < spacedim; ++j)
917 J[i] += grad_F[i][j] * (F - p)[j];
920 for (
unsigned int i = 0; i < facedim; ++i)
921 for (
unsigned int j = 0; j < facedim; ++j)
922 for (
unsigned int k = 0; k < spacedim; ++k)
923 H[i][j] += grad_F[i][k] * grad_F[j][k];
930 ExcMessage(
"The Newton iteration to find the reference point " 931 "did not converge in 10 iterations. Do you have a " 932 "deformed cell? (See the glossary for a definition " 933 "of what a deformed cell is. You may want to output " 934 "the vertices of your cell."));
943 const double normalized_delta_world = (F - p).norm() / face->diameter();
945 if (delta_xi.
norm() < eps || normalized_delta_world < eps)
953 return internal::normalized_alternating_product(grad_F);
958 template <
int dim,
int spacedim,
int chartdim>
961 : sub_manifold(periodicity)
966 template <
int dim,
int spacedim,
int chartdim>
971 const double w)
const 973 const std::array<Point<spacedim>, 2> points{{p1, p2}};
974 const std::array<double, 2> weights{{1. - w, w}};
981 template <
int dim,
int spacedim,
int chartdim>
987 const std::size_t n_points = surrounding_points.size();
989 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
991 for (
unsigned int i = 0; i < n_points; ++i)
992 chart_points[i] = pull_back(surrounding_points[i]);
997 return push_forward(p_chart);
1002 template <
int dim,
int spacedim,
int chartdim>
1012 const std::size_t n_points = surrounding_points.size();
1014 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1015 for (std::size_t i = 0; i < n_points; ++i)
1016 chart_points[i] = pull_back(surrounding_points[i]);
1018 boost::container::small_vector<Point<chartdim>, 200> new_points_on_chart(
1020 sub_manifold.get_new_points(
1023 make_array_view(new_points_on_chart.begin(), new_points_on_chart.end()));
1025 for (std::size_t row = 0; row < weights.
size(0); ++row)
1026 new_points[row] = push_forward(new_points_on_chart[row]);
1031 template <
int dim,
int spacedim,
int chartdim>
1044 template <
int dim,
int spacedim,
int chartdim>
1051 push_forward_gradient(pull_back(x1));
1060 1e-12 * F_prime.
norm(),
1062 "The derivative of a chart function must not be singular."));
1065 sub_manifold.get_tangent_vector(pull_back(x1), pull_back(x2));
1068 for (
unsigned int i = 0; i < spacedim; ++i)
1069 result[i] += F_prime[i] * delta;
1076 template <
int dim,
int spacedim,
int chartdim>
1080 return sub_manifold.get_periodicity();
1084 #include "manifold.inst" 1086 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcPureFunctionCalled()
static const unsigned int invalid_unsigned_int
FlatManifold(const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
std::pair< std::array< Point< MeshIteratorType::AccessorType::space_dimension >, n_default_points_per_cell< MeshIteratorType >)>, std::array< double, n_default_points_per_cell< MeshIteratorType >)> > get_default_points_and_weights(const MeshIteratorType &iterator, const bool with_interpolation=false)
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
#define AssertDimension(dim1, dim2)
typename IteratorSelector::line_iterator line_iterator
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &t)
const Tensor< 1, spacedim > & get_periodicity() const
virtual Point< spacedim > get_new_point_on_hex(const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
typename IteratorSelector::hex_iterator hex_iterator
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim >>> &Du)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
numbers::NumberTraits< Number >::real_type norm() const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const
#define Assert(cond, exc)
numbers::NumberTraits< Number >::real_type norm_square() const
const Tensor< 1, chartdim > & get_periodicity() const
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim >> &points, const Point< spacedim > &candidate) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
typename IteratorSelector::quad_iterator quad_iterator
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
Tensor< 1, dim, Number > cross_product_3d(const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2)
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 Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
virtual void get_new_points(const ArrayView< const Point< spacedim >> &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim >> new_points) const
size_type size(const unsigned int i) const
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
static ::ExceptionBase & ExcEmptyObject()
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::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face > FaceVertexNormals
static ::ExceptionBase & ExcNotImplemented()
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) 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
ChartManifold(const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
Point< spacedim > get_new_point_on_face(const typename Triangulation< dim, spacedim >::face_iterator &face) const
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
static ::ExceptionBase & ExcInternalError()
Point< spacedim > get_new_point_on_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const