Reference documentation for deal.II version 9.0.0
|
#include <deal.II/grid/manifold_lib.h>
Public Member Functions | |
TransfiniteInterpolationManifold () | |
virtual | ~TransfiniteInterpolationManifold () override |
virtual std::unique_ptr< Manifold< dim, spacedim > > | clone () const override |
void | initialize (const Triangulation< dim, spacedim > &triangulation) |
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 |
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 > | project_to_manifold (const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) 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 Tensor< 1, spacedim > | get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const |
virtual Tensor< 1, spacedim > | normal_vector (const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) 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) |
Private Member Functions | |
std::array< unsigned int, 20 > | get_possible_cells_around_points (const ArrayView< const Point< spacedim >> &surrounding_points) const |
Triangulation< dim, spacedim >::cell_iterator | compute_chart_points (const ArrayView< const Point< spacedim >> &surrounding_points, ArrayView< Point< dim >> chart_points) const |
Point< dim > | pull_back (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_guess) const |
Point< spacedim > | push_forward (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &chart_point) const |
DerivativeForm< 1, dim, spacedim > | push_forward_gradient (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &chart_point, const Point< spacedim > &pushed_forward_chart_point) const |
Private Attributes | |
const Triangulation< dim, spacedim > * | triangulation |
int | level_coarse |
std::vector< bool > | coarse_cell_is_flat |
FlatManifold< dim > | chart_manifold |
boost::signals2::connection | clear_signal |
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) |
A mapping class that extends curved boundary descriptions into the interior of the computational domain. The outer curved boundary description is assumed to be given by another manifold (e.g. a polar manifold on a circle). The mechanism to extend the boundary information is a so-called transfinite interpolation.
The formula for extending such a description in 2D is, for example, described on Wikipedia. Given a point \((u,v)\) on the chart, the image of this point in real space is given by
\begin{align*} \mathbf S(u,v) &= (1-v)\mathbf c_0(u)+v \mathbf c_1(u) + (1-u)\mathbf c_2(v) + u \mathbf c_3(v) \\ &\quad - \left[(1-u)(1-v) \mathbf x_0 + u(1-v) \mathbf x_1 + (1-u)v \mathbf x_2 + uv \mathbf x_3 \right] \end{align*}
where \(\bf x_0, \bf x_1, \bf x_2, \bf x_3\) denote the four bounding vertices bounding the image space and \(\bf c_0, \bf c_1, \bf c_2, \bf c_3\) are the four curves describing the lines of the cell. If a curved manifold is attached to any of these lines, the evaluation is done according to Manifold::get_new_point() with the two end points of the line and appropriate weight. In 3D, the generalization of this formula is implemented, creating a weighted sum of the vertices (positive contribution), the lines (negative), and the faces (positive contribution).
This manifold is usually attached to a coarse mesh and then places new points as a combination of the descriptions on the boundaries, weighted appropriately according to the position of the point in the original chart coordinates \((u,v)\). This manifold should be preferred over setting only a curved manifold on the boundary of a mesh in most situations as it yields more uniform mesh distributions as the mesh is refined because it switches from a curved description to a straight description over all children of the initial coarse cell this manifold was attached to. This way, the curved nature of the manifold that is originally contained in one coarse mesh layer will be applied to more than one fine mesh layer once the mesh gets refined. Note that the mechanisms of TransfiniteInterpolationManifold are also built into the MappingQGeneric class when only a surface of a cell is subject to a curved description, ensuring that even the default case without this manifold gets optimal convergence rates when applying curved boundary descriptions.
If no curved boundaries surround a coarse cell, this class reduces to a flat manifold description.
To give an example of using this class, the following code attaches a transfinite manifold to a circle:
In this code, we first set all manifold ids to the id of the transfinite interpolation, and then re-set the manifold ids on the boundary to identify the curved boundary described by the polar manifold. With this code, one gets a really nice mesh:
which is obviously much nicer than the polar manifold applied to just the boundary:
In the implementation of this class, the manifolds surrounding a coarse cell are queried repeatedly to compute points on their interior. For optimal mesh quality, those manifolds should be compatible with a chart notion. For example, computing a point that is 0.25 along the line between two vertices using the weights 0.25 and 0.75 for the two vertices should give the same result as first computing the mid point at 0.5 and then again compute the midpoint between the first vertex and coarse mid point. This is the case for most of the manifold classes provided by deal.II, such as SphericalManifold or PolarManifold, but it might be violated by naive implementations. In case the quality of the manifold is not good enough, upon mesh refinement it may happen that the transformation to a chart inside the get_new_point() or get_new_points() methods produces points that are outside the unit cell. Then this class throws an exception of type Mapping::ExcTransformationFailed. In that case, the mesh should be refined before attaching this class, as done in the following example:
Definition at line 826 of file manifold_lib.h.
TransfiniteInterpolationManifold< dim, spacedim >::TransfiniteInterpolationManifold | ( | ) |
Constructor.
Definition at line 1340 of file manifold_lib.cc.
|
overridevirtual |
Destructor.
Definition at line 1351 of file manifold_lib.cc.
|
overridevirtual |
Make a clone of this Manifold object.
Implements Manifold< dim, spacedim >.
Definition at line 1361 of file manifold_lib.cc.
void TransfiniteInterpolationManifold< dim, spacedim >::initialize | ( | const Triangulation< dim, spacedim > & | triangulation | ) |
Initializes the manifold with a coarse mesh. The prerequisite for using this class is that the input triangulation is uniformly refined and the manifold is later attached to the same triangulation.
Whenever the assignment of manifold ids changes on the level of the triangulation which this class was initialized with, initialize() must be called again to update the manifold ids connected to the coarse cells.
Definition at line 1374 of file manifold_lib.cc.
|
overridevirtual |
Return the point which shall become the new vertex surrounded by the given points surrounding_points
. weights
contains appropriate weights for the surrounding points according to which the manifold determines the new point's position.
The implementation in this class overrides the method in the base class and computes the new point by a transfinite interpolation. The first step in the implementation is to identify the coarse cell on which the surrounding points are located. Then, the coordinates are transformed to the unit coordinates on the coarse cell by a Newton iteration, where the new point is then computed according to the weights. Finally, it is pushed forward to the real space according to the transfinite interpolation.
Reimplemented from Manifold< dim, spacedim >.
Definition at line 2169 of file manifold_lib.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 columns in weights
must match the length of new_points
.
The implementation in this class overrides the method in the base class and computes the new point by a transfinite interpolation. The first step in the implementation is to identify the coarse cell on which the surrounding points are located. Then, the coordinates are transformed to the unit coordinates on the coarse cell by a Newton iteration, where the new points are then computed according to the weights. Finally, the is pushed forward to the real space according to the transfinite interpolation.
The implementation does not allow for surrounding_points
and new_points
to point to the same vector, so make sure to pass different objects into the function.
Reimplemented from Manifold< dim, spacedim >.
Definition at line 2187 of file manifold_lib.cc.
|
private |
Internal function to identify the most suitable cells (=charts) where the given surrounding points are located. We use a cheap algorithm to identify the cells and rank the cells by probability before we actually do the search inside the relevant cells. The cells are sorted by the distance of a Q1 approximation of the inverse mapping to the unit cell of the surrounding points. We expect at most 20 cells (it should be up to 8 candidates on a 3D structured mesh and a bit more on unstructured ones, typically we only get two or three), so get an array with 20 entries of a the indices cell->index()
.
Definition at line 1871 of file manifold_lib.cc.
|
private |
Finalizes the identification of the correct chart and populates chart_points
with the pullbacks of the surrounding points. This method internally calls get_possible_cells_around_points()
.
Return an iterator to the cell on which the chart is defined.
Definition at line 1947 of file manifold_lib.cc.
|
private |
Pull back operation into the unit coordinates on the given coarse cell.
This method is currently based on a Newton-like iteration to find the point in the origin. One may speed up the iteration by providing a good initial guess as the third argument. If no better point is known, use cell->real_to_unit_cell_affine_approximation(p)
Definition at line 1747 of file manifold_lib.cc.
|
private |
Push forward operation.
Definition at line 1699 of file manifold_lib.cc.
|
private |
Gradient of the push_forward method.
Definition at line 1719 of file manifold_lib.cc.
|
private |
The underlying triangulation.
Definition at line 984 of file manifold_lib.h.
|
private |
The level of the mesh cells where the transfinite approximation is applied, usually level 0.
Definition at line 990 of file manifold_lib.h.
|
private |
In case there all surrounding manifolds are the transfinite manifold or have default (invalid) manifold id, the manifold degenerates to a flat manifold and we can choose cheaper algorithms for the push_forward method.
Definition at line 997 of file manifold_lib.h.
|
private |
A flat manifold used to compute new points in the chart space where we use a FlatManifold description.
Definition at line 1003 of file manifold_lib.h.
|
private |
The connection to Triangulation::signals::clear that must be reset once this class goes out of scope.
Definition at line 1009 of file manifold_lib.h.