Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20:02+00:00
\(\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
Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
TransfiniteInterpolationManifold< dim, spacedim > Class Template Reference

#include <deal.II/grid/manifold_lib.h>

Inheritance diagram for TransfiniteInterpolationManifold< dim, spacedim >:
Inheritance graph
[legend]

Public Types

using FaceVertexNormals = std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face >
 

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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Computing the location of points.
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
 
Computing tangent vectors
virtual Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const
 
Computing normal vectors
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
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

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
 
void check_no_subscribers () const noexcept
 

Private Attributes

const Triangulation< dim, spacedim > * triangulation
 
int level_coarse
 
std::vector< boolcoarse_cell_is_flat
 
FlatManifold< dim > chart_manifold
 
std::vector< internal::MappingQImplementation::InverseQuadraticApproximation< dim, spacedim > > quadratic_approximation
 
boost::signals2::connection clear_signal
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<int dim, int spacedim = dim>
class TransfiniteInterpolationManifold< dim, spacedim >

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 use of this class is discussed extensively in step-65.

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 MappingQ 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:

PolarManifold<dim> polar_manifold;
triangulation.set_manifold (0, polar_manifold);
inner_manifold.initialize(triangulation);
triangulation.set_manifold (1, inner_manifold);
const Triangulation< dim, spacedim > * triangulation
void initialize(const Triangulation< dim, spacedim > &triangulation)
void refine_global(const unsigned int times=1)
void set_all_manifold_ids_on_boundary(const types::manifold_id number)
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
void set_all_manifold_ids(const types::manifold_id number)
void hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)

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:

This manifold is used in a few GridGenerator functions, including GridGenerator::channel_with_cylinder.

Implementation details

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:

Note
For performance and accuracy reasons, it is recommended to apply the transfinite manifold to as coarse a mesh as possible. Regarding accuracy, the curved description can only be applied to new points created from a given neighborhood, and the grid quality is typically higher when extending the curved description over as large a domain as possible. Regarding performance, the identification of the correct coarse cell in the get_new_point() method needs to pass all coarse cells, so expect a linear complexity in the number of coarse cells for each single mapping operation, i.e., at least quadratic in the number of coarse mesh cells for any global operation on the whole mesh. Thus, the current implementation is only economical when there are not more than a few hundreds of coarse cells. To make performance better for larger numbers of cells, one could extend the current implementation by a pre-identification of relevant cells with axis-aligned bounding boxes.

Definition at line 963 of file manifold_lib.h.

Member Typedef Documentation

◆ FaceVertexNormals

template<int dim, int spacedim = dim>
using Manifold< dim, spacedim >::FaceVertexNormals = std::array<Tensor<1, spacedim>, GeometryInfo<dim>::vertices_per_face>
inherited

Type keeping information about the normals at the vertices of a face of a cell. Thus, there are GeometryInfo<dim>::vertices_per_face normal vectors, that define the tangent spaces of the boundary at the vertices. Note that the vectors stored in this object are not required to be normalized, nor to actually point outward, as one often will only want to check for orthogonality to define the tangent plane; if a function requires the normals to be normalized, then it must do so itself.

For obvious reasons, this type is not useful in 1d.

Definition at line 305 of file manifold.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ TransfiniteInterpolationManifold()

template<int dim, int spacedim>
TransfiniteInterpolationManifold< dim, spacedim >::TransfiniteInterpolationManifold ( )

Constructor.

Definition at line 1604 of file manifold_lib.cc.

◆ ~TransfiniteInterpolationManifold()

template<int dim, int spacedim>
TransfiniteInterpolationManifold< dim, spacedim >::~TransfiniteInterpolationManifold ( )
overridevirtual

Destructor.

Definition at line 1615 of file manifold_lib.cc.

Member Function Documentation

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > TransfiniteInterpolationManifold< dim, spacedim >::clone ( ) const
overridevirtual

Make a clone of this Manifold object.

Implements Manifold< dim, spacedim >.

Definition at line 1625 of file manifold_lib.cc.

◆ initialize()

template<int dim, int spacedim>
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.

Note
The triangulation used to construct the manifold must not be destroyed during the usage of this object.

Definition at line 1637 of file manifold_lib.cc.

◆ get_new_point()

template<int dim, int spacedim>
Point< spacedim > TransfiniteInterpolationManifold< dim, spacedim >::get_new_point ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const ArrayView< const double > &  weights 
) const
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 2548 of file manifold_lib.cc.

◆ get_new_points()

template<int dim, int spacedim>
void TransfiniteInterpolationManifold< dim, spacedim >::get_new_points ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Table< 2, double > &  weights,
ArrayView< Point< spacedim > >  new_points 
) const
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 2568 of file manifold_lib.cc.

◆ get_possible_cells_around_points()

template<int dim, int spacedim>
std::array< unsigned int, 20 > TransfiniteInterpolationManifold< dim, spacedim >::get_possible_cells_around_points ( const ArrayView< const Point< spacedim > > &  surrounding_points) const
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 2201 of file manifold_lib.cc.

◆ compute_chart_points()

template<int dim, int spacedim>
Triangulation< dim, spacedim >::cell_iterator TransfiniteInterpolationManifold< dim, spacedim >::compute_chart_points ( const ArrayView< const Point< spacedim > > &  surrounding_points,
ArrayView< Point< dim > >  chart_points 
) const
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 2287 of file manifold_lib.cc.

◆ pull_back()

template<int dim, int spacedim>
Point< dim > TransfiniteInterpolationManifold< dim, spacedim >::pull_back ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< spacedim > &  p,
const Point< dim > &  initial_guess 
) const
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)

Note
This internal function is currently not compatible with the ChartManifold::pull_back() function because the given class represents an atlas of charts, not a single chart. Thus, the pull_back() operation is only valid with the additional information of the chart, given by a cell on the coarse grid. An alternative implementation could shift the index depending on the coarse cell for a 1-to-1 relation between the chart space and the image space.

Definition at line 2051 of file manifold_lib.cc.

◆ push_forward()

template<int dim, int spacedim>
Point< spacedim > TransfiniteInterpolationManifold< dim, spacedim >::push_forward ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< dim > &  chart_point 
) const
private

Push forward operation.

Note
This internal function is currently not compatible with the ChartManifold::push_forward() function because the given class represents an atlas of charts, not a single chart. Thus, the push_forward() operation is only valid with the additional information of the chart, given by a cell on the coarse grid. An alternative implementation could shift the index depending on the coarse cell for a 1-to-1 relation between the chart space and the image space.

Definition at line 2001 of file manifold_lib.cc.

◆ push_forward_gradient()

template<int dim, int spacedim>
DerivativeForm< 1, dim, spacedim > TransfiniteInterpolationManifold< 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

Gradient of the push_forward method.

Note
This internal function is not compatible with the ChartManifold::push_forward_gradient() function because the given class represents an atlas of charts, not a single chart. Furthermore, this private function also requires the user to provide the result of the push_forward() call on the chart point for the single use case of this function, namely inside a Newton iteration where the gradient is computed by finite differences.

Definition at line 2022 of file manifold_lib.cc.

◆ get_intermediate_point()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Manifold< dim, spacedim >::get_intermediate_point ( const Point< spacedim > &  p1,
const Point< spacedim > &  p2,
const double  w 
) const
virtualinherited

Return an intermediate point between two given points. Overloading this function allows the default pair-wise reduction implementation of the method get_new_point() that takes a Quadrature object as input to work properly.

An implementation of this function should returns a parametric curve on the manifold, joining the points p1 and p2, with parameter w in the interval [0,1]. In particular get_intermediate_point(p1, p2, 0.0) should return p1 and get_intermediate_point(p1, p2, 1.0) should return p2.

In its default implementation, this function calls the project_to_manifold() method with the convex combination of p1 and p2. User classes can get away by simply implementing the project_to_manifold() method.

Reimplemented in ChartManifold< dim, spacedim, chartdim >, ChartManifold< dim, 3, 3 >, ChartManifold< dim, dim, 3 >, ChartManifold< dim, dim, dim >, ChartManifold< dim, spacedim, 1 >, ChartManifold< dim, spacedim, 2 >, ChartManifold< dim, spacedim_A+spacedim_B, chartdim_A+chartdim_B >, ChartManifold< dim_A, spacedim_A, chartdim_A >, ChartManifold< dim_B, spacedim_B, chartdim_B >, and SphericalManifold< dim, spacedim >.

◆ project_to_manifold()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Manifold< dim, spacedim >::project_to_manifold ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Point< spacedim > &  candidate 
) const
virtualinherited

Given a point which lies close to the given manifold, it modifies it and projects it to manifold itself.

This class is used by the default implementation of the function get_new_point() and should be implemented by derived classes. The default implementation simply throws an exception if called.

If your manifold is simple, you could implement this function only, and the default behavior should work out of the box.

Reimplemented in FlatManifold< dim, spacedim >, FlatManifold< chartdim, chartdim >, FlatManifold< dim, dim >, FlatManifold< dim, spacedim >, OpenCASCADE::NormalProjectionManifold< dim, spacedim >, OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >, and OpenCASCADE::NormalToMeshProjectionManifold< dim, spacedim >.

◆ get_new_point_on_line()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Manifold< dim, spacedim >::get_new_point_on_line ( const typename Triangulation< dim, spacedim >::line_iterator &  line) const
virtualinherited

Backward compatibility interface. Return the point which shall become the new middle vertex of the two children of a regular line. In 2d, this line is a line at the boundary, while in 3d, it is bounding a face at the boundary (the lines therefore is also on the boundary).

The default implementation of this function passes its argument to the Manifolds::get_default_points_and_weights() function, and then calls the Manifold<dim,spacedim>::get_new_point() function. User derived classes can overload Manifold<dim,spacedim>::get_new_point() or Manifold<dim,spacedim>::project_to_manifold(), which is called by the default implementation of Manifold<dim,spacedim>::get_new_point().

◆ get_new_point_on_quad()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Manifold< dim, spacedim >::get_new_point_on_quad ( const typename Triangulation< dim, spacedim >::quad_iterator &  quad) const
virtualinherited

Backward compatibility interface. Return the point which shall become the common point of the four children of a quad at the boundary in three or more spatial dimensions. This function therefore is only useful in at least three dimensions and should not be called for lower dimensions.

This function is called after the four lines bounding the given quad are refined, so you may want to use the information provided by quad->line(i)->child(j), i=0...3, j=0,1.

The default implementation of this function passes its argument to the Manifolds::get_default_points_and_weights() function, and then calls the Manifold<dim,spacedim>::get_new_point() function. User derived classes can overload Manifold<dim,spacedim>::get_new_point() or Manifold<dim,spacedim>::project_to_manifold(), which is called by the default implementation of Manifold<dim,spacedim>::get_new_point().

◆ get_new_point_on_hex()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Manifold< dim, spacedim >::get_new_point_on_hex ( const typename Triangulation< dim, spacedim >::hex_iterator &  hex) const
virtualinherited

Backward compatibility interface. Return the point which shall become the common point of the eight children of a hex in three or spatial dimensions. This function therefore is only useful in at least three dimensions and should not be called for lower dimensions.

This function is called after the all the bounding objects of the given hex are refined, so you may want to use the information provided by hex->quad(i)->line(j)->child(k), i=0...5, j=0...3, k=0,1.

The default implementation of this function passes its argument to the Manifolds::get_default_points_and_weights() function, and then calls the Manifold<dim,spacedim>::get_new_point() function. User derived classes can overload Manifold<dim,spacedim>::get_new_point() or Manifold<dim,spacedim>::project_to_manifold(), which is called by the default implementation of Manifold<dim,spacedim>::get_new_point().

◆ get_new_point_on_face()

template<int dim, int spacedim = dim>
Point< spacedim > Manifold< dim, spacedim >::get_new_point_on_face ( const typename Triangulation< dim, spacedim >::face_iterator &  face) const
inherited

Backward compatibility interface. Depending on dim=2 or dim=3 this function calls the get_new_point_on_line or the get_new_point_on_quad function. It throws an exception for dim=1. This wrapper allows dimension independent programming.

◆ get_new_point_on_cell()

template<int dim, int spacedim = dim>
Point< spacedim > Manifold< dim, spacedim >::get_new_point_on_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
inherited

Backward compatibility interface. Depending on dim=1, dim=2 or dim=3 this function calls the get_new_point_on_line, get_new_point_on_quad or the get_new_point_on_hex function. This wrapper allows dimension independent programming.

◆ get_tangent_vector()

template<int dim, int spacedim = dim>
virtual Tensor< 1, spacedim > Manifold< dim, spacedim >::get_tangent_vector ( const Point< spacedim > &  x1,
const Point< spacedim > &  x2 
) const
virtualinherited

Return a vector that, at \(\mathbf x_1\), is tangential to the geodesic that connects two points \(\mathbf x_1,\mathbf x_2\). The geodesic is the shortest line between these two points, where "shortest" is defined via a metric specific to a particular implementation of this class in a derived class. For example, in the case of a FlatManifold, the shortest line between two points is just the straight line, and in this case the tangent vector is just the difference \(\mathbf d=\mathbf x_2-\mathbf x_1\). On the other hand, for a manifold that describes a surface embedded in a higher dimensional space (e.g., the surface of a sphere), then the tangent vector is tangential to the surface, and consequently may point in a different direction than the straight line that connects the two points.

While tangent vectors are often normalized to unit length, the vectors returned by this function are normalized as described in the introduction of this class. Specifically, if \(\mathbf s(t)\) traces out the geodesic between the two points where \(\mathbf x_1 = \mathbf s(0)\) and \(\mathbf x_2 = \mathbf s(1)\), then the returned vector must equal \(\mathbf s'(0)\). In other words, the norm of the returned vector also encodes, in some sense, the length of the geodesic because a curve \(\mathbf s(t)\) must move "faster" if the two points it connects between arguments \(t=0\) and \(t=1\) are farther apart.

The default implementation of this function approximates \(\mathbf s'(0) \approx \frac{\mathbf s(\epsilon)-\mathbf x_1}{\epsilon}\) for a small value of \(\epsilon\), and the evaluation of \(\mathbf s(\epsilon)\) is done by calling get_new_point(). If possible, derived classes should override this function by an implementation of the exact derivative.

Parameters
x1The first point that describes the geodesic, and the one at which the "direction" is to be evaluated.
x2The second point that describes the geodesic.
Returns
A "direction" vector tangential to the geodesic.

Reimplemented in FlatManifold< dim, spacedim >, FlatManifold< chartdim, chartdim >, FlatManifold< dim, dim >, FlatManifold< dim, spacedim >, ChartManifold< dim, spacedim, chartdim >, ChartManifold< dim, 3, 3 >, ChartManifold< dim, dim, 3 >, ChartManifold< dim, dim, dim >, ChartManifold< dim, spacedim, 1 >, ChartManifold< dim, spacedim, 2 >, ChartManifold< dim, spacedim_A+spacedim_B, chartdim_A+chartdim_B >, ChartManifold< dim_A, spacedim_A, chartdim_A >, ChartManifold< dim_B, spacedim_B, chartdim_B >, and SphericalManifold< dim, spacedim >.

◆ normal_vector()

template<int dim, int spacedim = dim>
virtual Tensor< 1, spacedim > Manifold< dim, spacedim >::normal_vector ( const typename Triangulation< dim, spacedim >::face_iterator &  face,
const Point< spacedim > &  p 
) const
virtualinherited

Return the normal vector to a face embedded in this manifold, at the point p. It is not required that the normals actually point outward from the domain even if the face iterator given points to a face on the boundary of the domain. If p is not in fact on the surface, but only close-by, try to return something reasonable, for example the normal vector at the surface point closest to p. (The point p will in fact not normally lie on the actual surface, but rather be a quadrature point mapped by some polynomial mapping; the mapped surface, however, will not usually coincide with the actual surface.)

This function only makes sense if dim==spacedim because otherwise there is no unique normal vector but in fact a (spacedim-dim+1)-dimensional tangent space of vectors that are all both normal to the face and normal to the dim-dimensional surface that lives in spacedim-dimensional space. For example, think of a two-dimensional mesh that covers a two-dimensional surface in three-dimensional space. In that case, each face (edge) is one-dimensional, and there are two linearly independent vectors that are both normal to the edge: one is normal to the edge and tangent to the surface (intuitively, that would be the one that points from the current cell to the neighboring one, if the surface was locally flat), and the other one is rooted in the edge but points perpendicular to the surface (which is also perpendicular to the edge that lives within the surface). Thus, because there are no obviously correct semantics for this function if spacedim is greater than dim, the function will simply throw an error in that situation.

The face iterator gives an indication which face this function is supposed to compute the normal vector for. This is useful if the boundary of the domain is composed of different nondifferential pieces (for example when using the FlatManifold class to approximate a geometry that is completely described by the coarse mesh, with piecewise (bi-)linear components between the vertices, but where the boundary may have a kink at the vertices itself).

Note
In 2d, the default implementation of this function computes the normal vector by taking the tangent direction from p to the further one of the two vertices that make up an edge, and then rotates it outward (with respect to the coordinate system of the edge) by 90 degrees. In 3d, the default implementation is more complicated, aiming at avoiding problems with numerical round-off for points close to one of the vertices, and avoiding tangent directions that are linearly dependent.

Reimplemented in FlatManifold< dim, spacedim >, FlatManifold< chartdim, chartdim >, FlatManifold< dim, dim >, FlatManifold< dim, spacedim >, PolarManifold< dim, spacedim >, and SphericalManifold< dim, spacedim >.

◆ get_normals_at_vertices()

template<int dim, int spacedim = dim>
virtual void Manifold< dim, spacedim >::get_normals_at_vertices ( const typename Triangulation< dim, spacedim >::face_iterator &  face,
FaceVertexNormals face_vertex_normals 
) const
virtualinherited

Compute the normal vectors to the boundary at each vertex of the given face embedded in the Manifold. It is not required that the normal vectors be normed somehow. Neither is it required that the normals actually point outward.

This function is needed to compute data for C1 mappings. The default implementation calls normal_vector() on each vertex.

Note that when computing normal vectors at a vertex where the boundary is not differentiable, you have to make sure that you compute the one-sided limits, i.e. limit with respect to points inside the given face.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ triangulation

template<int dim, int spacedim = dim>
const Triangulation<dim, spacedim>* TransfiniteInterpolationManifold< dim, spacedim >::triangulation
private

The underlying triangulation.

Definition at line 1124 of file manifold_lib.h.

◆ level_coarse

template<int dim, int spacedim = dim>
int TransfiniteInterpolationManifold< dim, spacedim >::level_coarse
private

The level of the mesh cells where the transfinite approximation is applied, usually level 0.

Definition at line 1130 of file manifold_lib.h.

◆ coarse_cell_is_flat

template<int dim, int spacedim = dim>
std::vector<bool> TransfiniteInterpolationManifold< dim, spacedim >::coarse_cell_is_flat
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 1137 of file manifold_lib.h.

◆ chart_manifold

template<int dim, int spacedim = dim>
FlatManifold<dim> TransfiniteInterpolationManifold< dim, spacedim >::chart_manifold
private

A flat manifold used to compute new points in the chart space where we use a FlatManifold description.

Definition at line 1143 of file manifold_lib.h.

◆ quadratic_approximation

template<int dim, int spacedim = dim>
std::vector<internal::MappingQImplementation:: InverseQuadraticApproximation<dim, spacedim> > TransfiniteInterpolationManifold< dim, spacedim >::quadratic_approximation
private

A vector of quadratic approximations to the inverse map from real points to chart points for each of the coarse mesh cells.

Definition at line 1151 of file manifold_lib.h.

◆ clear_signal

template<int dim, int spacedim = dim>
boost::signals2::connection TransfiniteInterpolationManifold< dim, spacedim >::clear_signal
private

The connection to Triangulation::signals::clear that must be reset once this class goes out of scope.

Definition at line 1157 of file manifold_lib.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


The documentation for this class was generated from the following files: