27#include <boost/container/small_vector.hpp>
38template <
int dim,
int spacedim>
50template <
int dim,
int spacedim>
56 const std::array<Point<spacedim>, 2>
vertices{{p1, p2}};
58 w * p2 + (1 -
w) * p1);
63template <
int dim,
int spacedim>
69 const double tol = 1
e-10;
70 const unsigned int n_points = surrounding_points.size();
76 "There should be as many surrounding points as weights given."));
80 ExcMessage(
"The weights for the individual points should sum to 1!"));
85 boost::container::small_vector<unsigned int, 100> permutation(n_points);
86 std::iota(permutation.begin(), permutation.end(), 0u);
87 std::sort(permutation.begin(),
89 [&weights](
const std::size_t a,
const std::size_t
b) {
90 return weights[a] < weights[b];
95 double w = weights[permutation[0]];
97 for (
unsigned int i = 1; i < n_points; ++i)
100 if (
std::abs(weights[permutation[i]] +
w) < tol)
103 weight =
w / (weights[permutation[i]] +
w);
108 surrounding_points[permutation[i]],
113 p = surrounding_points[permutation[i]];
115 w += weights[permutation[i]];
123template <
int dim,
int spacedim>
132 for (
unsigned int row = 0; row < weights.size(0); ++row)
136 surrounding_points.end()),
148 const int spacedim = 2;
153 ((p - face->vertex(0)).norm_square() > (p - face->vertex(1)).norm_square() ?
154 -get_tangent_vector(p, face->vertex(0)) :
155 get_tangent_vector(p, face->vertex(1)));
159 return normal / normal.
norm();
169 const int spacedim = 3;
171 const std::array<Point<spacedim>, 4>
vertices{
172 {face->vertex(0), face->vertex(1), face->vertex(2), face->vertex(3)}};
178 std::max(distances[2], distances[3]));
189 for (
unsigned int i = 0; i < 3; ++i)
190 if (distances[i] > 1
e-8 * max_distance)
191 for (
unsigned int j = i + 1; j < 4; ++j)
192 if (distances[j] > 1
e-8 * max_distance)
195 (distances[i] * distances[j]);
198 if (
std::abs(new_angle) < 0.999 * abs_cos_angle)
200 abs_cos_angle =
std::abs(new_angle);
206 ExcMessage(
"The search for possible directions did not succeed."));
217 "Manifold::normal_vector was unable to find a suitable combination "
218 "of vertices to compute a normal on this face. We chose the secants "
219 "that are as orthogonal as possible, but tangents appear to be "
220 "linearly dependent. Check for distorted faces in your triangulation."));
229 if (reference_normal * normal < 0.0)
232 return normal / normal.
norm();
237template <
int dim,
int spacedim>
255 n[0] =
cross_product_2d(get_tangent_vector(face->vertex(0), face->vertex(1)));
259 for (
unsigned int i = 0; i < 2; ++i)
263 "computed normals have "
277 n[0] =
cross_product_3d(get_tangent_vector(face->vertex(0), face->vertex(1)),
278 get_tangent_vector(face->vertex(0), face->vertex(2)));
280 n[1] =
cross_product_3d(get_tangent_vector(face->vertex(1), face->vertex(3)),
281 get_tangent_vector(face->vertex(1), face->vertex(0)));
283 n[2] =
cross_product_3d(get_tangent_vector(face->vertex(2), face->vertex(0)),
284 get_tangent_vector(face->vertex(2), face->vertex(3)));
286 n[3] =
cross_product_3d(get_tangent_vector(face->vertex(3), face->vertex(2)),
287 get_tangent_vector(face->vertex(3), face->vertex(1)));
289 for (
unsigned int i = 0; i < 4; ++i)
293 "computed normals have "
301template <
int dim,
int spacedim>
305 FaceVertexNormals & n)
const
307 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
316template <
int dim,
int spacedim>
323 points_weights.first.end()),
325 points_weights.second.end()));
330template <
int dim,
int spacedim>
337 points_weights.first.end()),
339 points_weights.second.end()));
344template <
int dim,
int spacedim>
364template <
int dim,
int spacedim>
450template <
int dim,
int spacedim>
468 points_weights.first.end()),
470 points_weights.second.end()));
475template <
int dim,
int spacedim>
482 const std::array<Point<spacedim>, 2> points{{x1, x2}};
487 return (neighbor_point - x1) /
epsilon;
497 normalized_alternating_product(
const Tensor<1, 3> (&)[1])
509 normalized_alternating_product(
const Tensor<1, 3> (&basis_vectors)[2])
512 return tmp / tmp.
norm();
518template <
int dim,
int spacedim>
521 const double tolerance)
522 : periodicity(periodicity)
523 , tolerance(tolerance)
528template <
int dim,
int spacedim>
529std::unique_ptr<Manifold<dim, spacedim>>
537template <
int dim,
int spacedim>
545 ExcMessage(
"The weights for the individual points should sum to 1!"));
552 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
553 p += surrounding_points[i] * weights[i];
559 for (
unsigned int d = 0;
d < spacedim; ++
d)
561 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
563 minP[
d] =
std::min(minP[
d], surrounding_points[i][
d]);
564 Assert((surrounding_points[i][
d] <
566 (surrounding_points[i][
d] >=
573 for (
unsigned int i = 0; i < surrounding_points.size(); ++i)
576 for (
unsigned int d = 0;
d < spacedim; ++
d)
583 p += (surrounding_points[i] + dp) * weights[i];
587 for (
unsigned int d = 0;
d < spacedim; ++
d)
598template <
int dim,
int spacedim>
606 if (weights.size(0) == 0)
609 const std::size_t n_points = surrounding_points.size();
612 for (
unsigned int d = 0;
d < spacedim; ++
d)
614 for (
unsigned int i = 0; i < n_points; ++i)
616 minP[
d] =
std::min(minP[
d], surrounding_points[i][
d]);
617 Assert((surrounding_points[i][
d] <
625 const Point<spacedim> *surrounding_points_start = surrounding_points.data();
627 boost::container::small_vector<Point<spacedim>, 200> modified_points;
628 bool adjust_periodicity =
false;
629 for (
unsigned int d = 0;
d < spacedim; ++
d)
631 for (
unsigned int i = 0; i < n_points; ++i)
634 adjust_periodicity =
true;
637 if (adjust_periodicity ==
true)
639 modified_points.resize(surrounding_points.size());
641 surrounding_points.end(),
642 modified_points.begin());
643 for (
unsigned int d = 0;
d < spacedim; ++
d)
645 for (
unsigned int i = 0; i < n_points; ++i)
648 surrounding_points_start = modified_points.data();
652 for (
unsigned int row = 0; row < weights.size(0); ++row)
656 std::accumulate(&weights(row, 0), &weights(row, 0) + n_points, 0.0) -
658 ExcMessage(
"The weights for the individual points should sum to 1!"));
660 for (
unsigned int p = 0; p < n_points; ++p)
661 new_point += surrounding_points_start[p] * weights(row, p);
664 for (
unsigned int d = 0;
d < spacedim; ++
d)
666 if (new_point[
d] < 0)
673 surrounding_points.end()),
680template <
int dim,
int spacedim>
691template <
int dim,
int spacedim>
700template <
int dim,
int spacedim>
711 for (
unsigned int d = 0;
d < spacedim; ++
d)
764 const Tensor<1, 2> tangent = face->vertex(1) - face->vertex(0);
765 for (
unsigned int vertex = 0; vertex < GeometryInfo<2>::vertices_per_face;
768 face_vertex_normals[vertex] =
Point<2>(tangent[1], -tangent[0]);
792 static const unsigned int neighboring_vertices[4][2] = {{1, 2},
796 for (
unsigned int vertex = 0; vertex < vertices_per_face; ++vertex)
801 face->vertex(neighboring_vertices[vertex][0]) - face->vertex(vertex),
802 face->vertex(neighboring_vertices[vertex][1]) - face->vertex(vertex)};
858template <
int dim,
int spacedim>
894 const unsigned int facedim = dim - 1;
898 const auto face_reference_cell = face->reference_cell();
900 if (face_reference_cell == ReferenceCells::get_hypercube<facedim>())
902 for (
unsigned int i = 0; i < facedim; ++i)
907 for (
unsigned int i = 0; i < facedim; ++i)
911 const double eps = 1
e-12;
913 unsigned int iteration = 0;
917 for (
const unsigned int v : face->vertex_indices())
919 face->vertex(v) * face_reference_cell.d_linear_shape_function(xi, v);
921 for (
unsigned int i = 0; i < facedim; ++i)
924 for (
const unsigned int v : face->vertex_indices())
927 face_reference_cell.d_linear_shape_function_gradient(xi, v)[i];
931 for (
unsigned int i = 0; i < facedim; ++i)
932 for (
unsigned int j = 0; j < spacedim; ++j)
933 J[i] += grad_F[i][j] * (
F - p)[j];
936 for (
unsigned int i = 0; i < facedim; ++i)
937 for (
unsigned int j = 0; j < facedim; ++j)
938 for (
unsigned int k = 0; k < spacedim; ++k)
939 H[i][j] += grad_F[i][k] * grad_F[j][k];
946 ExcMessage(
"The Newton iteration to find the reference point "
947 "did not converge in 10 iterations. Do you have a "
948 "deformed cell? (See the glossary for a definition "
949 "of what a deformed cell is. You may want to output "
950 "the vertices of your cell."));
959 const double normalized_delta_world = (
F - p).
norm() / face->diameter();
961 if (delta_xi.
norm() <
eps || normalized_delta_world <
eps)
969 return internal::normalized_alternating_product(grad_F);
974template <
int dim,
int spacedim,
int chartdim>
977 : sub_manifold(periodicity)
982template <
int dim,
int spacedim,
int chartdim>
987 const double w)
const
989 const std::array<Point<spacedim>, 2> points{{p1, p2}};
990 const std::array<double, 2> weights{{1. -
w,
w}};
997template <
int dim,
int spacedim,
int chartdim>
1003 const std::size_t n_points = surrounding_points.size();
1005 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1007 for (
unsigned int i = 0; i < n_points; ++i)
1008 chart_points[i] =
pull_back(surrounding_points[i]);
1018template <
int dim,
int spacedim,
int chartdim>
1028 const std::size_t n_points = surrounding_points.size();
1030 boost::container::small_vector<Point<chartdim>, 200> chart_points(n_points);
1031 for (std::size_t i = 0; i < n_points; ++i)
1032 chart_points[i] =
pull_back(surrounding_points[i]);
1034 boost::container::small_vector<Point<chartdim>, 200> new_points_on_chart(
1036 sub_manifold.get_new_points(
1039 make_array_view(new_points_on_chart.begin(), new_points_on_chart.end()));
1041 for (std::size_t row = 0; row < weights.size(0); ++row)
1042 new_points[row] =
push_forward(new_points_on_chart[row]);
1047template <
int dim,
int spacedim,
int chartdim>
1060template <
int dim,
int spacedim,
int chartdim>
1076 1
e-12 * F_prime.
norm(),
1078 "The derivative of a chart function must not be singular."));
1084 for (
unsigned int i = 0; i < spacedim; ++i)
1085 result[i] += F_prime[i] * delta;
1092template <
int dim,
int spacedim,
int chartdim>
1096 return sub_manifold.get_periodicity();
1100#include "manifold.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
ChartManifold(const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
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 DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const
const Tensor< 1, chartdim > & get_periodicity() const
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) 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
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim > > &points, const Point< spacedim > &candidate) const override
const Tensor< 1, spacedim > & get_periodicity() const
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
const Tensor< 1, spacedim > periodicity
FlatManifold(const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
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_on_hex(const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const
virtual void get_normals_at_vertices(const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const
std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face > FaceVertexNormals
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const
Point< spacedim > get_new_point_on_face(const typename Triangulation< dim, spacedim >::face_iterator &face) const
virtual Point< spacedim > get_new_point_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
virtual Tensor< 1, spacedim > normal_vector(const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
Point< spacedim > get_new_point_on_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
numbers::NumberTraits< Number >::real_type norm() const
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
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)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcPeriodicBox(int arg1, Point< spacedim > arg2, double arg3)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcPureFunctionCalled()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
typename IteratorSelector::quad_iterator quad_iterator
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
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)
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
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)
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)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
void copy(const T *begin, const T *end, U *dest)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)