Reference documentation for deal.II version 9.6.0
\(\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
FlatManifold< dim, spacedim > Class Template Reference

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

Inheritance diagram for FlatManifold< dim, spacedim >:

Public Types

using FaceVertexNormals
 

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
 
double get_tolerance () const
 
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 > 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 normal vectors
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

void check_no_subscribers () const noexcept
 

Static Private Member Functions

static ::ExceptionBaseExcPeriodicBox (int arg1, Point< spacedim > arg2, double arg3)
 

Private Attributes

const Tensor< 1, spacedim > periodicity
 
const double tolerance
 
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 FlatManifold< dim, spacedim >

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 634 of file manifold.h.

Member Typedef Documentation

◆ FaceVertexNormals

using Manifold< dim, spacedim >::FaceVertexNormals
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

◆ FlatManifold()

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

Member Function Documentation

◆ clone()

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

◆ get_new_point()

template<int dim, int spacedim = dim>
virtual Point< spacedim > FlatManifold< dim, spacedim >::get_new_point ( const ArrayView< const Point< spacedim > > & surrounding_points,
const ArrayView< const double > & weights ) const
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, dim >.

◆ get_new_points()

template<int dim, int spacedim = dim>
virtual void FlatManifold< 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 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, dim >.

◆ project_to_manifold()

template<int dim, int spacedim = dim>
virtual Point< spacedim > FlatManifold< dim, spacedim >::project_to_manifold ( const ArrayView< const Point< spacedim > > & points,
const Point< spacedim > & candidate ) const
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, dim >.

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

◆ get_tangent_vector()

template<int dim, int spacedim = dim>
virtual Tensor< 1, spacedim > FlatManifold< dim, spacedim >::get_tangent_vector ( const Point< spacedim > & x1,
const Point< spacedim > & x2 ) const
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().

Note
If you use this class as a stepping stone to build a manifold that only "slightly" deviates from a flat manifold, by overloading the project_to_manifold() function.
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. Here, this is \(\mathbf x_2-\mathbf x_1\), possibly modified by the periodicity of the domain as set in the constructor, to use the "shortest" connection between the points through the periodic boundary as necessary.

Reimplemented from Manifold< dim, dim >.

◆ normal_vector()

template<int dim, int spacedim = dim>
virtual Tensor< 1, spacedim > FlatManifold< dim, spacedim >::normal_vector ( const typename Triangulation< dim, spacedim >::face_iterator & face,
const Point< spacedim > & p ) const
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, dim >.

◆ get_normals_at_vertices() [1/2]

template<int dim, int spacedim = dim>
virtual void FlatManifold< dim, spacedim >::get_normals_at_vertices ( const typename Triangulation< dim, spacedim >::face_iterator & face,
typename Manifold< dim, spacedim >::FaceVertexNormals & face_vertex_normals ) const
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.

◆ get_periodicity()

template<int dim, int spacedim = dim>
const Tensor< 1, spacedim > & FlatManifold< dim, spacedim >::get_periodicity ( ) const

Return the periodicity of this Manifold.

◆ get_tolerance()

template<int dim, int spacedim = dim>
double FlatManifold< dim, spacedim >::get_tolerance ( ) const

Return the relative tolerance set in the constructor.

◆ get_intermediate_point()

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 >, and SphericalManifold< dim, spacedim >.

◆ get_new_point_on_line()

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()

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()

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()

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()

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_normals_at_vertices() [2/2]

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

◆ periodicity

template<int dim, int spacedim = dim>
const Tensor<1, spacedim> FlatManifold< dim, spacedim >::periodicity
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 802 of file manifold.h.

◆ tolerance

template<int dim, int spacedim = dim>
const double FlatManifold< dim, spacedim >::tolerance
private

Relative tolerance. This tolerance is used to compute distances in double precision.

Definition at line 816 of file manifold.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 file: