Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
manifold.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_tria_manifold_h
17#define dealii_tria_manifold_h
18
19
20/*---------------------------- manifold.h ---------------------------*/
21
22#include <deal.II/base/config.h>
23
26#include <deal.II/base/point.h>
29
30#include <deal.II/grid/tria.h>
31
33
34// forward declaration
35#ifndef DOXYGEN
36template <int, typename>
37class Table;
38#endif
39
44namespace Manifolds
45{
52 template <typename MeshIteratorType>
53 inline constexpr std::size_t
55 {
56 // Note that in C++11 a constexpr function can only have a return
57 // statement, so we cannot alias the structure dimension
59 vertices_per_cell +
61 lines_per_cell +
63 quads_per_cell +
65 hexes_per_cell -
66 1; // don't count the cell itself, just the bounding objects
67 }
68
111 template <typename MeshIteratorType>
112 std::pair<std::array<Point<MeshIteratorType::AccessorType::space_dimension>,
113 n_default_points_per_cell<MeshIteratorType>()>,
114 std::array<double, n_default_points_per_cell<MeshIteratorType>()>>
115 get_default_points_and_weights(const MeshIteratorType &iterator,
116 const bool with_interpolation = false);
117} // namespace Manifolds
118
119
120
285template <int dim, int spacedim = dim>
286class Manifold : public Subscriptor
287{
288public:
289 // explicitly check for sensible template arguments
290 static_assert(dim <= spacedim,
291 "The dimension <dim> of a Manifold must be less than or "
292 "equal to the space dimension <spacedim> in which it lives.");
293
294
307 std::array<Tensor<1, spacedim>, GeometryInfo<dim>::vertices_per_face>;
308
309
314 virtual ~Manifold() override = default;
315
321 virtual std::unique_ptr<Manifold<dim, spacedim>>
322 clone() const = 0;
323
328
346 virtual Point<spacedim>
348 const Point<spacedim> &p2,
349 const double w) const;
350
367 virtual Point<spacedim>
368 get_new_point(const ArrayView<const Point<spacedim>> &surrounding_points,
369 const ArrayView<const double> & weights) const;
370
371
392 virtual void
393 get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
394 const Table<2, double> & weights,
395 ArrayView<Point<spacedim>> new_points) const;
396
408 virtual Point<spacedim>
410 const ArrayView<const Point<spacedim>> &surrounding_points,
411 const Point<spacedim> & candidate) const;
412
426 virtual Point<spacedim>
428 const typename Triangulation<dim, spacedim>::line_iterator &line) const;
429
447 virtual Point<spacedim>
449 const typename Triangulation<dim, spacedim>::quad_iterator &quad) const;
450
469 virtual Point<spacedim>
471 const typename Triangulation<dim, spacedim>::hex_iterator &hex) const;
472
473
482 const typename Triangulation<dim, spacedim>::face_iterator &face) const;
483
484
493 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
494
496
501
538 virtual Tensor<1, spacedim>
540 const Point<spacedim> &x2) const;
541
543
548
594 virtual Tensor<1, spacedim>
597 const Point<spacedim> & p) const;
598
613 virtual void
616 FaceVertexNormals &face_vertex_normals) const;
617
619};
620
621
631template <int dim, int spacedim = dim>
632class FlatManifold : public Manifold<dim, spacedim>
633{
634public:
663 const double tolerance = 1e-10);
664
668 virtual std::unique_ptr<Manifold<dim, spacedim>>
669 clone() const override;
670
692 virtual Point<spacedim>
693 get_new_point(const ArrayView<const Point<spacedim>> &surrounding_points,
694 const ArrayView<const double> & weights) const override;
695
696
707 virtual void
708 get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
709 const Table<2, double> & weights,
710 ArrayView<Point<spacedim>> new_points) const override;
711
719 virtual Point<spacedim>
721 const Point<spacedim> &candidate) const override;
722
744 virtual Tensor<1, spacedim>
746 const Point<spacedim> &x2) const override;
747
755 virtual Tensor<1, spacedim>
758 const Point<spacedim> &p) const override;
759
768 virtual void
771 typename Manifold<dim, spacedim>::FaceVertexNormals &face_vertex_normals)
772 const override;
773
777 const Tensor<1, spacedim> &
779
780private:
795
797 int,
799 double,
800 << "The component number " << arg1 << " of the point [ "
801 << arg2 << " ] is not in the interval [ 0, " << arg3
802 << "), bailing out.");
803
808 const double tolerance;
809};
810
811
898template <int dim, int spacedim = dim, int chartdim = dim>
899class ChartManifold : public Manifold<dim, spacedim>
900{
901public:
902 // explicitly check for sensible template arguments
903 static_assert(dim <= spacedim,
904 "The dimension <dim> of a ChartManifold must be less than or "
905 "equal to the space dimension <spacedim> in which it lives.");
906
922
927 virtual ~ChartManifold() override = default;
928
933 virtual Point<spacedim>
935 const Point<spacedim> &p2,
936 const double w) const override;
937
942 virtual Point<spacedim>
943 get_new_point(const ArrayView<const Point<spacedim>> &surrounding_points,
944 const ArrayView<const double> & weights) const override;
945
967 virtual void
968 get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
969 const Table<2, double> & weights,
970 ArrayView<Point<spacedim>> new_points) const override;
977 virtual Point<chartdim>
978 pull_back(const Point<spacedim> &space_point) const = 0;
979
986 virtual Point<spacedim>
987 push_forward(const Point<chartdim> &chart_point) const = 0;
988
1006 push_forward_gradient(const Point<chartdim> &chart_point) const;
1007
1063 virtual Tensor<1, spacedim>
1065 const Point<spacedim> &x2) const override;
1066
1070 const Tensor<1, chartdim> &
1071 get_periodicity() const;
1072
1073private:
1086};
1087
1088
1089/* -------------- declaration of explicit specializations ------------- */
1090
1091#ifndef DOXYGEN
1092
1093template <>
1097
1098template <>
1102
1103
1104template <>
1108
1109
1110template <>
1114
1115template <>
1119
1120
1121template <>
1125
1126
1127template <>
1130 const Triangulation<3, 3>::hex_iterator &) const;
1131
1132/*---Templated functions---*/
1133
1134namespace Manifolds
1135{
1136 template <typename MeshIteratorType>
1137 std::pair<std::array<Point<MeshIteratorType::AccessorType::space_dimension>,
1138 n_default_points_per_cell<MeshIteratorType>()>,
1139 std::array<double, n_default_points_per_cell<MeshIteratorType>()>>
1140 get_default_points_and_weights(const MeshIteratorType &iterator,
1141 const bool with_interpolation)
1142 {
1143 const int dim = MeshIteratorType::AccessorType::structure_dimension;
1144 const int spacedim = MeshIteratorType::AccessorType::space_dimension;
1145 constexpr std::size_t points_per_cell =
1146 n_default_points_per_cell<MeshIteratorType>();
1147
1148 std::pair<std::array<Point<spacedim>, points_per_cell>,
1149 std::array<double, points_per_cell>>
1150 points_weights;
1151
1152
1153 // note that the exact weights are chosen such as to minimize the
1154 // distortion of the four new quads from the optimal shape; their
1155 // derivation and values is copied over from the
1156 // interpolation function in the mapping
1157 switch (dim)
1158 {
1159 case 1:
1160 Assert(points_weights.first.size() == 2, ExcInternalError());
1161 Assert(points_weights.second.size() == 2, ExcInternalError());
1162 points_weights.first[0] = iterator->vertex(0);
1163 points_weights.second[0] = .5;
1164 points_weights.first[1] = iterator->vertex(1);
1165 points_weights.second[1] = .5;
1166 break;
1167 case 2:
1168 Assert(points_weights.first.size() == 8, ExcInternalError());
1169 Assert(points_weights.second.size() == 8, ExcInternalError());
1170
1171 for (unsigned int i = 0; i < 4; ++i)
1172 {
1173 points_weights.first[i] = iterator->vertex(i);
1174 points_weights.first[4 + i] =
1175 (iterator->line(i)->has_children() ?
1176 iterator->line(i)->child(0)->vertex(1) :
1177 iterator->line(i)->get_manifold().get_new_point_on_line(
1178 iterator->line(i)));
1179 }
1180
1181 if (with_interpolation)
1182 {
1183 std::fill(points_weights.second.begin(),
1184 points_weights.second.begin() + 4,
1185 -0.25);
1186 std::fill(points_weights.second.begin() + 4,
1187 points_weights.second.end(),
1188 0.5);
1189 }
1190 else
1191 std::fill(points_weights.second.begin(),
1192 points_weights.second.end(),
1193 1.0 / 8.0);
1194 break;
1195 case 3:
1196 {
1198 static_cast<TriaIterator<TriaAccessor<3, 3, 3>>>(iterator);
1199 const unsigned int np = GeometryInfo<dim>::vertices_per_cell +
1202 Assert(points_weights.first.size() == np, ExcInternalError());
1203 Assert(points_weights.second.size() == np, ExcInternalError());
1204 auto *sp3 = reinterpret_cast<
1205 std::array<Point<3>, n_default_points_per_cell<decltype(hex)>()>
1206 *>(&points_weights.first);
1207
1208 unsigned int j = 0;
1209
1210 // note that the exact weights are chosen such as to minimize the
1211 // distortion of the eight new hexes from the optimal shape through
1212 // transfinite interpolation from the faces and vertices, see
1213 // TransfiniteInterpolationManifold for a deeper explanation of the
1214 // mechanisms
1215 if (with_interpolation)
1216 {
1217 for (unsigned int i = 0;
1218 i < GeometryInfo<dim>::vertices_per_cell;
1219 ++i, ++j)
1220 {
1221 (*sp3)[j] = hex->vertex(i);
1222 points_weights.second[j] = 1.0 / 8.0;
1223 }
1224 for (unsigned int i = 0; i < GeometryInfo<dim>::lines_per_cell;
1225 ++i, ++j)
1226 {
1227 (*sp3)[j] =
1228 (hex->line(i)->has_children() ?
1229 hex->line(i)->child(0)->vertex(1) :
1230 hex->line(i)->get_manifold().get_new_point_on_line(
1231 hex->line(i)));
1232 points_weights.second[j] = -1.0 / 4.0;
1233 }
1234 for (unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
1235 ++i, ++j)
1236 {
1237 (*sp3)[j] =
1238 (hex->quad(i)->has_children() ?
1239 hex->quad(i)->isotropic_child(0)->vertex(3) :
1240 hex->quad(i)->get_manifold().get_new_point_on_quad(
1241 hex->quad(i)));
1242 points_weights.second[j] = 1.0 / 2.0;
1243 }
1244 }
1245 else
1246 // Overwrite the weights with 1/np if we don't want to use
1247 // interpolation.
1248 std::fill(points_weights.second.begin(),
1249 points_weights.second.end(),
1250 1.0 / np);
1251 }
1252 break;
1253 default:
1254 Assert(false, ExcInternalError());
1255 break;
1256 }
1257 return points_weights;
1258 }
1259} // namespace Manifolds
1260
1261#endif // DOXYGEN
1262
1264
1265#endif
virtual Point< spacedim > push_forward(const Point< chartdim > &chart_point) const =0
const FlatManifold< chartdim, chartdim > sub_manifold
Definition: manifold.h:1085
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const override
Definition: manifold.cc:1020
virtual Point< chartdim > pull_back(const Point< spacedim > &space_point) const =0
virtual Tensor< 1, spacedim > get_tangent_vector(const Point< spacedim > &x1, const Point< spacedim > &x2) const override
Definition: manifold.cc:1062
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient(const Point< chartdim > &chart_point) const
Definition: manifold.cc:1049
const Tensor< 1, chartdim > & get_periodicity() const
Definition: manifold.cc:1094
virtual ~ChartManifold() override=default
virtual Point< spacedim > get_intermediate_point(const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const override
Definition: manifold.cc:984
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const override
Definition: manifold.cc:999
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
Definition: manifold.h:794
FlatManifold(const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
const double tolerance
Definition: manifold.h:808
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 ~Manifold() override=default
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
Definition: manifold.h:307
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 std::unique_ptr< Manifold< dim, spacedim > > clone() const =0
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
Definition: point.h:111
Definition: tensor.h:503
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcPeriodicBox(int arg1, Point< spacedim > arg2, double arg3)
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcInternalError()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:555
typename IteratorSelector::hex_iterator hex_iterator
Definition: tria.h:1490
typename IteratorSelector::quad_iterator quad_iterator
Definition: tria.h:1466
typename IteratorSelector::line_iterator line_iterator
Definition: tria.h:1442
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)
constexpr std::size_t n_default_points_per_cell()
Definition: manifold.h:54