Reference documentation for deal.II version 9.0.0
|
#include <deal.II/grid/manifold.h>
Public Member Functions | |
FlatManifold (const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10) | |
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 |
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 > | project_to_manifold (const ArrayView< const Point< spacedim >> &points, const Point< spacedim > &candidate) const override |
virtual Tensor< 1, spacedim > | get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) 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 |
const Tensor< 1, spacedim > & | get_periodicity () const |
Public Member Functions inherited from Manifold< dim, spacedim > | |
virtual | ~Manifold ()=default |
virtual Point< spacedim > | get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const |
virtual Point< spacedim > | get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const |
virtual Point< spacedim > | get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const |
virtual Point< spacedim > | get_new_point_on_hex (const typename Triangulation< dim, spacedim >::hex_iterator &hex) const |
Point< spacedim > | get_new_point_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face) const |
Point< spacedim > | get_new_point_on_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const |
virtual void | get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Private Member Functions | |
static ::ExceptionBase & | ExcPeriodicBox (int arg1, Point< spacedim > arg2, double arg3) |
Private Attributes | |
const Tensor< 1, spacedim > | periodicity |
const double | tolerance |
Additional Inherited Members | |
Public Types inherited from Manifold< dim, spacedim > | |
typedef Tensor< 1, spacedim > | FaceVertexNormals[GeometryInfo< dim >::vertices_per_face] |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Specialization of Manifold<dim,spacedim>, which represent a possibly periodic Euclidean space of dimension dim
embedded in the Euclidean space of spacedim
dimensions. The main characteristic of this Manifold is the fact that the function FlatManifold<dim,spacedim>::project_to_manifold() is the identity function.
Definition at line 669 of file manifold.h.
FlatManifold< dim, spacedim >::FlatManifold | ( | const Tensor< 1, spacedim > & | periodicity = Tensor<1,spacedim>() , |
const double | tolerance = 1e-10 |
||
) |
Default constructor. The optional argument can be used to specify the periodicity of the spacedim-dimensional manifold (one period per direction). A periodicity value of zero means that along that direction there is no periodicity. By default no periodicity is assumed.
Periodicity affects the way a middle point is computed. It is assumed that if two points are more than half period distant, then the distance should be computed by crossing the periodicity boundary, i.e., the average is computed by adding a full period to the sum of the two. For example, if along direction 0 we have 2*pi periodicity, then the average of (2*pi-eps) and (eps) is not pi, but 2*pi (or zero), since, on a periodic manifold, these two points are at distance 2*eps and not (2*pi- eps). Special cases are taken into account, to ensure that the behavior is always as expected. The third argument is used as a relative tolerance when computing distances.
Periodicity will be intended in the following way: the domain is considered to be the box contained in [Point<spacedim>(), periodicity) where the right extreme is excluded. If any of the components of this box has zero length, then no periodicity is assumed in that direction. Whenever a function that tries to compute averages is called, an exception will be thrown if one of the points which you are using for the average lies outside the periodicity box. The return points are guaranteed to lie in the periodicity box plus or minus tolerance*periodicity.norm().
Definition at line 507 of file manifold.cc.
|
overridevirtual |
Return a copy of this manifold.
Implements Manifold< dim, spacedim >.
Reimplemented in TorusBoundary< dim, spacedim >, HalfHyperShellBoundary< dim >, HyperShellBoundary< dim >, HalfHyperBallBoundary< dim >, HyperBallBoundary< dim, spacedim >, HyperBallBoundary< dim >, StraightBoundary< dim, spacedim >, StraightBoundary< dim, dim >, StraightBoundary< dim >, OpenCASCADE::NormalToMeshProjectionBoundary< dim, spacedim >, ConeBoundary< dim >, OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >, Boundary< dim, spacedim >, Boundary< dim, dim >, OpenCASCADE::NormalProjectionBoundary< dim, spacedim >, and CylinderBoundary< dim, spacedim >.
Definition at line 518 of file manifold.cc.
|
overridevirtual |
Let the new point be the average sum of surrounding vertices.
This particular implementation constructs the weighted average of the surrounding points, and then calls internally the function project_to_manifold(). The reason why we do it this way, is to allow lazy programmers to implement only the project_to_manifold() function for their own Manifold classes which are small (or trivial) perturbations of a flat manifold. This is the case whenever the coarse mesh is a decent approximation of the manifold geometry. In this case, the middle point of a cell is close to true middle point of the manifold, and a projection may suffice.
For most simple geometries, it is possible to get reasonable results by deriving your own Manifold class from FlatManifold, and write a new interface only for the project_to_manifold function. You will have good approximations also with large deformations, as long as in the coarsest mesh size you are trying to refine, the middle point is not too far from the manifold mid point, i.e., as long as the coarse mesh size is small enough.
Reimplemented from Manifold< dim, spacedim >.
Definition at line 528 of file manifold.cc.
|
overridevirtual |
Compute a new set of points that interpolate between the given points surrounding_points
. weights
is a table with as many columns as surrounding_points.size()
. The number of rows in weights
must match the length of new_points
.
For this particular implementation, the interpolation of the surrounding_points
according to the weights
is simply performed in Cartesian space.
Reimplemented from Manifold< dim, spacedim >.
Definition at line 583 of file manifold.cc.
|
overridevirtual |
Project to FlatManifold. This is the identity function for flat, Euclidean spaces. Note however that this function can be overloaded by derived classes, which will then benefit from the logic behind the get_new_point() function which are often very similar (if not identical) to the one implemented in this class.
Reimplemented from Manifold< dim, spacedim >.
Reimplemented in OpenCASCADE::NormalToMeshProjectionBoundary< dim, spacedim >, OpenCASCADE::DirectionalProjectionBoundary< dim, spacedim >, and OpenCASCADE::NormalProjectionBoundary< dim, spacedim >.
Definition at line 659 of file manifold.cc.
|
overridevirtual |
Return a vector that, at \(\mathbf x_1\), is tangential to the geodesic that connects two points \(\mathbf x_1,\mathbf x_2\). For the current class, we assume that the manifold is flat, so the geodesic is the straight line between the two points, and we return \(\mathbf x_2-\mathbf x_1\). The normalization of the vector is chosen so that it fits the convention described in Manifold::get_tangent_vector().
x1 | The first point that describes the geodesic, and the one at which the "direction" is to be evaluated. |
x2 | The second point that describes the geodesic. |
Reimplemented from Manifold< dim, spacedim >.
Definition at line 678 of file manifold.cc.
|
overridevirtual |
Return the normal vector to the given face at point p taking into account that quadrilateral faces of hexahedral cells in 3d may not be planar. In those cases, the face is assumed to have a geometry described by a bilinear function, and the normal vector is computed by embedding this bilinear form into a Cartesian space with a flat metric.
Reimplemented from Manifold< dim, spacedim >.
Reimplemented in HyperBallBoundary< dim, spacedim >, StraightBoundary< dim, spacedim >, and StraightBoundary< dim, dim >.
Definition at line 841 of file manifold.cc.
|
overridevirtual |
Compute the normal vectors to the boundary at each vertex of the given face taking into account that quadrilateral faces of hexahedral cells in 3d may not be planar. In those cases, the face is assumed to have a geometry described by a bilinear function, and the normal vector is computed by embedding this bilinear form into a Cartesian space with a flat metric.
const Tensor< 1, spacedim > & FlatManifold< dim, spacedim >::get_periodicity | ( | ) | const |
Return the periodicity of this Manifold.
Definition at line 669 of file manifold.cc.
|
private |
The periodicity of this Manifold. Periodicity affects the way a middle point is computed. It is assumed that if two points are more than half period distant, then the distance should be computed by crossing the periodicity boundary, i.e., the average is computed by adding a full period to the sum of the two. For example, if along direction 0 we have 2*pi periodicity, then the average of (2*pi-eps) and (eps) is not pi, but 2*pi (or zero), since, on a periodic manifold, these two points are at distance 2*eps and not (2*pi-eps).
A periodicity 0 along one direction means no periodicity. This is the default value for all directions.
Definition at line 832 of file manifold.h.
|
private |
Relative tolerance. This tolerance is used to compute distances in double precision.
Definition at line 842 of file manifold.h.