Reference documentation for deal.II version 9.3.3
\(\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\}}\)
Public Types | Public Member Functions | Public Attributes | Static Private Member Functions | Private Attributes | List of all members
PolarManifold< dim, spacedim > Class Template Referenceabstract

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

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

Public Types

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

Public Member Functions

 PolarManifold (const Point< spacedim > center=Point< spacedim >())
 
virtual std::unique_ptr< Manifold< dim, spacedim > > clone () const override
 
virtual Point< spacedim > pull_back (const Point< spacedim > &space_point) const override
 
virtual Point< spacedim > push_forward (const Point< spacedim > &chart_point) const override
 
virtual DerivativeForm< 1, spacedim, spacedim > push_forward_gradient (const Point< spacedim > &chart_point) const override
 
virtual Tensor< 1, spacedim > normal_vector (const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) 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 > push_forward (const Point< chartdim > &chart_point) const=0
 
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient (const Point< chartdim > &chart_point) const
 
virtual Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const override
 
const Tensor< 1, chartdim > & get_periodicity () const
 
Point< 1 > get_new_point_on_quad (const Triangulation< 1, 1 >::quad_iterator &) const
 
Point< 2 > get_new_point_on_quad (const Triangulation< 1, 2 >::quad_iterator &) const
 
Point< 3 > get_new_point_on_quad (const Triangulation< 1, 3 >::quad_iterator &) const
 
Point< 3 > get_new_point_on_hex (const Triangulation< 3, 3 >::hex_iterator &hex) const
 
Point< 1 > get_new_point_on_face (const Triangulation< 1, 1 >::face_iterator &) const
 
Point< 2 > get_new_point_on_face (const Triangulation< 1, 2 >::face_iterator &) const
 
Point< 3 > get_new_point_on_face (const Triangulation< 1, 3 >::face_iterator &) const
 
Tensor< 1, 2 > normal_vector (const Triangulation< 2, 2 >::face_iterator &face, const Point< 2 > &p) const
 
Tensor< 1, 3 > normal_vector (const Triangulation< 3, 3 >::face_iterator &face, const Point< 3 > &p) const
 
void get_normals_at_vertices (const Triangulation< 2, 2 >::face_iterator &face, FaceVertexNormals &n) const
 
void get_normals_at_vertices (const Triangulation< 3, 3 >::face_iterator &face, FaceVertexNormals &n) const
 
Computing the location of points.
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 normal vectors
virtual void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
 

Public Attributes

const Point< spacedim > center
 

Static Private Member Functions

static Tensor< 1, spacedim > get_periodicity ()
 

Private Attributes

const FlatManifold< chartdim, chartdim > sub_manifold
 

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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

Detailed Description

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

Manifold description for a polar coordinate system.

You can use this Manifold object to describe any sphere, circle, hypersphere or hyperdisc in two or three dimensions, both as a co-dimension one manifold descriptor or as co-dimension zero manifold descriptor, provided that the north and south poles (in three dimensions) and the center (in both two and three dimensions) are excluded from the Manifold (as they are singular points of the polar change of coordinates).

The two template arguments match the meaning of the two template arguments in Triangulation<dim, spacedim>, however this Manifold can be used to describe both thin and thick objects, and the behavior is identical when dim <= spacedim, i.e., the functionality of PolarManifold<2,3> is identical to PolarManifold<3,3>.

This class works by transforming points to polar coordinates (in both two and three dimensions), taking the average in that coordinate system, and then transforming the point back to Cartesian coordinates. In order for this manifold to work correctly, it cannot be attached to cells containing the center of the coordinate system or the north and south poles in three dimensions. These points are singular points of the coordinate transformation, and taking averages around these points does not make any sense.

Definition at line 70 of file manifold_lib.h.

Member Typedef Documentation

◆ FaceVertexNormals

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

Constructor & Destructor Documentation

◆ PolarManifold()

template<int dim, int spacedim>
PolarManifold< dim, spacedim >::PolarManifold ( const Point< spacedim >  center = Point<spacedim>())

The Constructor takes the center of the spherical coordinates system. This class uses the pull_back and push_forward mechanism to transform from Cartesian to spherical coordinate systems, taking into account the periodicity of base Manifold in two dimensions, while in three dimensions it takes the middle point, and project it along the radius using the average radius of the surrounding points.

Definition at line 127 of file manifold_lib.cc.

Member Function Documentation

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< Manifold< dim, spacedim > > PolarManifold< dim, spacedim >::clone
overridevirtual

Make a clone of this Manifold object.

Implements Manifold< dim, spacedim >.

Definition at line 137 of file manifold_lib.cc.

◆ pull_back()

template<int dim, int spacedim>
Point< spacedim > PolarManifold< dim, spacedim >::pull_back ( const Point< spacedim > &  space_point) const
overridevirtual

Pull back the given point from the Euclidean space. Will return the polar coordinates associated with the point space_point. Only used when spacedim = 2.

Implements ChartManifold< dim, dim, dim >.

Definition at line 197 of file manifold_lib.cc.

◆ push_forward() [1/2]

template<int dim, int spacedim>
Point< spacedim > PolarManifold< dim, spacedim >::push_forward ( const Point< spacedim > &  chart_point) const
overridevirtual

Given a point in the spherical coordinate system, this method returns the Euclidean coordinates associated to the polar coordinates chart_point. Only used when spacedim = 3.

Definition at line 163 of file manifold_lib.cc.

◆ push_forward_gradient() [1/2]

template<int dim, int spacedim>
DerivativeForm< 1, spacedim, spacedim > PolarManifold< dim, spacedim >::push_forward_gradient ( const Point< spacedim > &  chart_point) const
overridevirtual

Given a point in the spacedim dimensional Euclidean space, this method returns the derivatives of the function \(F\) that maps from the polar coordinate system to the Euclidean coordinate system. In other words, it is a matrix of size \(\text{spacedim}\times\text{spacedim}\).

This function is used in the computations required by the get_tangent_vector() function.

Refer to the general documentation of this class for more information.

Definition at line 236 of file manifold_lib.cc.

◆ normal_vector() [1/3]

template<int dim, int spacedim>
Tensor< 1, spacedim > PolarManifold< dim, spacedim >::normal_vector ( const typename Triangulation< dim, spacedim >::face_iterator &  face,
const Point< spacedim > &  p 
) const
overridevirtual

Return the normal vector to a face embedded in this manifold, at the point p. 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 from Manifold< dim, spacedim >.

Definition at line 325 of file manifold_lib.cc.

◆ get_periodicity() [1/2]

template<int dim, int spacedim>
Tensor< 1, spacedim > PolarManifold< dim, spacedim >::get_periodicity
staticprivate

Helper function which returns the periodicity associated with this coordinate system, according to dim, chartdim, and spacedim.

Definition at line 146 of file manifold_lib.cc.

◆ get_intermediate_point()

Point< spacedim > ChartManifold< dim, spacedim, chartdim >::get_intermediate_point ( const Point< spacedim > &  p1,
const Point< spacedim > &  p2,
const double  w 
) const
overridevirtualinherited

Refer to the general documentation of this class and the documentation of the base class for more information.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 934 of file manifold.cc.

◆ get_new_point()

Point< spacedim > ChartManifold< dim, spacedim, chartdim >::get_new_point ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const ArrayView< const double > &  weights 
) const
overridevirtualinherited

Refer to the general documentation of this class and the documentation of the base class for more information.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 943 of file manifold.cc.

◆ get_new_points()

void ChartManifold< dim, spacedim, chartdim >::get_new_points ( const ArrayView< const Point< spacedim > > &  surrounding_points,
const Table< 2, double > &  weights,
ArrayView< Point< spacedim > >  new_points 
) const
overridevirtualinherited

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.

The implementation of this function first transforms the surrounding_points to the chart space by calling pull_back(). Then, new points are computed on the chart by usual interpolation according to the given weights, which are finally transformed to the image space by push_forward().

This implementation can be much more efficient for computing multiple new points from the same surrounding points than separate calls to get_new_point() in case the pull_back() operation is expensive. This is because pull_back() is only called once for the surrounding points and the interpolation is done for all given weights using this set of points. Often, pull_back() is also more expensive than push_forward() because the former might involve some kind of Newton iteration in non-trivial manifolds.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 968 of file manifold.cc.

◆ push_forward() [2/2]

virtual Point< spacedim > ChartManifold< dim, spacedim, chartdim >::push_forward ( const Point< chartdim > &  chart_point) const
pure virtualinherited

Given a point in the chartdim dimensional Euclidean space, this method returns a point on the manifold embedded in the spacedim Euclidean space.

Refer to the general documentation of this class for more information.

Implemented in CompositionManifold< dim, spacedim, chartdim, intermediate_dim, dim1, dim2 >, and FunctionManifold< dim, spacedim, chartdim >.

◆ push_forward_gradient() [2/2]

DerivativeForm< 1, chartdim, spacedim > ChartManifold< dim, spacedim, chartdim >::push_forward_gradient ( const Point< chartdim > &  chart_point) const
virtualinherited

Given a point in the chartdim dimensional Euclidean space, this method returns the derivatives of the function \(F\) that maps from the chartdim-dimensional to the spacedim-dimensional space. In other words, it is a matrix of size \(\text{spacedim}\times\text{chartdim}\).

This function is used in the computations required by the get_tangent_vector() function. Since not all users of the Manifold class interface will require calling that function, the current function is implemented but will trigger an exception whenever called. This allows derived classes to avoid implementing the push_forward_gradient function if this functionality is not needed in the user program.

Refer to the general documentation of this class for more information.

Reimplemented in CompositionManifold< dim, spacedim, chartdim, intermediate_dim, dim1, dim2 >, and FunctionManifold< dim, spacedim, chartdim >.

Definition at line 1006 of file manifold.cc.

◆ get_tangent_vector()

Tensor< 1, spacedim > ChartManifold< dim, spacedim, chartdim >::get_tangent_vector ( const Point< spacedim > &  x1,
const Point< spacedim > &  x2 
) const
overridevirtualinherited

Return a vector that, at \(\mathbf x_1\), is tangential to the geodesic that connects two points \(\mathbf x_1,\mathbf x_2\). See the documentation of the Manifold class and of Manifold::get_tangent_vector() for a more detailed description.

For the current class, we assume that this geodesic is the image under the push_forward() operation of a straight line of the pre-images of x1 and x2 (where pre-images are computed by pulling back the locations x1 and x2). In other words, if these preimages are \(\xi_1=F^{-1}(\mathbf x_1), \xi_2=F^{-1}(\mathbf x_2)\), then the geodesic in preimage (the chartdim-dimensional Euclidean) space is

\begin{align*} \zeta(t) &= \xi_1 + t (\xi_2-\xi_1) \\ &= F^{-1}(\mathbf x_1) + t\left[F^{-1}(\mathbf x_2) -F^{-1}(\mathbf x_1)\right] \end{align*}

In image space, i.e., in the space in which we operate, this leads to the curve

\begin{align*} \mathbf s(t) &= F(\zeta(t)) \\ &= F(\xi_1 + t (\xi_2-\xi_1)) \\ &= F\left(F^{-1}(\mathbf x_1) + t\left[F^{-1}(\mathbf x_2) -F^{-1}(\mathbf x_1)\right]\right). \end{align*}

What the current function is supposed to return is \(\mathbf s'(0)\). By the chain rule, this is equal to

\begin{align*} \mathbf s'(0) &= \frac{d}{dt}\left. F\left(F^{-1}(\mathbf x_1) + t\left[F^{-1}(\mathbf x_2) -F^{-1}(\mathbf x_1)\right]\right) \right|_{t=0} \\ &= \nabla_\xi F\left(F^{-1}(\mathbf x_1)\right) \left[F^{-1}(\mathbf x_2) -F^{-1}(\mathbf x_1)\right]. \end{align*}

This formula may then have to be slightly modified by considering any periodicity that was assumed in the call to the constructor.

Thus, the computation of tangent vectors also requires the implementation of derivatives \(\nabla_\xi F(\xi)\) of the push-forward mapping. Here, \(F^{-1}(\mathbf x_2)-F^{-1}(\mathbf x_1)\) is a chartdim-dimensional vector, and \(\nabla_\xi F\left(F^{-1}(\mathbf x_1)\right) = \nabla_\xi F\left(\xi_1\right)\) is a spacedim-times-chartdim-dimensional matrix. Consequently, and as desired, the operation results in a spacedim-dimensional vector.

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 from Manifold< dim, spacedim >.

Definition at line 1064 of file manifold.cc.

◆ get_periodicity() [2/2]

const Tensor< 1, chartdim > & ChartManifold< dim, spacedim, chartdim >::get_periodicity
inherited

Return the periodicity associated with the submanifold.

Definition at line 1071 of file manifold.cc.

◆ project_to_manifold()

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< chartdim, chartdim >, FlatManifold< dim, spacedim >, FlatManifold< dim, dim >, OpenCASCADE::NormalProjectionManifold< dim, spacedim >, OpenCASCADE::DirectionalProjectionManifold< dim, spacedim >, and OpenCASCADE::NormalToMeshProjectionManifold< 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() [1/4]

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

Point< 1 > Manifold< 1, 1 >::get_new_point_on_quad ( const Triangulation< 1, 1 >::quad_iterator &  ) const
inherited

Definition at line 419 of file manifold.cc.

◆ get_new_point_on_quad() [3/4]

Point< 2 > Manifold< 1, 2 >::get_new_point_on_quad ( const Triangulation< 1, 2 >::quad_iterator &  ) const
inherited

Definition at line 430 of file manifold.cc.

◆ get_new_point_on_quad() [4/4]

Point< 3 > Manifold< 1, 3 >::get_new_point_on_quad ( const Triangulation< 1, 3 >::quad_iterator &  ) const
inherited

Definition at line 441 of file manifold.cc.

◆ get_new_point_on_hex() [1/2]

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

Point< 3 > Manifold< 3, 3 >::get_new_point_on_hex ( const Triangulation< 3, 3 >::hex_iterator &  hex) const
inherited

Definition at line 463 of file manifold.cc.

◆ get_new_point_on_face() [1/4]

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

Point< 1 > Manifold< 1, 1 >::get_new_point_on_face ( const Triangulation< 1, 1 >::face_iterator &  ) const
inherited

Definition at line 386 of file manifold.cc.

◆ get_new_point_on_face() [3/4]

Point< 2 > Manifold< 1, 2 >::get_new_point_on_face ( const Triangulation< 1, 2 >::face_iterator &  ) const
inherited

Definition at line 397 of file manifold.cc.

◆ get_new_point_on_face() [4/4]

Point< 3 > Manifold< 1, 3 >::get_new_point_on_face ( const Triangulation< 1, 3 >::face_iterator &  ) const
inherited

Definition at line 408 of file manifold.cc.

◆ 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.

◆ normal_vector() [2/3]

Tensor< 1, 2 > Manifold< 2, 2 >::normal_vector ( const Triangulation< 2, 2 >::face_iterator &  face,
const Point< 2 > &  p 
) const
inherited

Definition at line 145 of file manifold.cc.

◆ normal_vector() [3/3]

Tensor< 1, 3 > Manifold< 3, 3 >::normal_vector ( const Triangulation< 3, 3 >::face_iterator &  face,
const Point< 3 > &  p 
) const
inherited

Definition at line 166 of file manifold.cc.

◆ get_normals_at_vertices() [1/3]

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.

◆ get_normals_at_vertices() [2/3]

void Manifold< 2, 2 >::get_normals_at_vertices ( const Triangulation< 2, 2 >::face_iterator &  face,
FaceVertexNormals n 
) const
inherited

Definition at line 251 of file manifold.cc.

◆ get_normals_at_vertices() [3/3]

void Manifold< 3, 3 >::get_normals_at_vertices ( const Triangulation< 3, 3 >::face_iterator &  face,
FaceVertexNormals n 
) const
inherited

Definition at line 273 of file manifold.cc.

Member Data Documentation

◆ center

template<int dim, int spacedim = dim>
const Point<spacedim> PolarManifold< dim, spacedim >::center

The center of the spherical coordinate system.

Definition at line 131 of file manifold_lib.h.

◆ sub_manifold

const FlatManifold<chartdim, chartdim> ChartManifold< dim, spacedim, chartdim >::sub_manifold
privateinherited

The sub_manifold object is used to compute the average of the points in the chart coordinates system.

In an ideal world, it would have type FlatManifold<dim,chartdim>. However, this would instantiate cases where dim>spacedim, which leads to invalid situations. We instead use <chartdim,chartdim>, which is (i) always valid, and (ii) does not matter at all since the first (dim) argument of manifolds is, in fact, ignored as far as manifold functionality is concerned.

Definition at line 1085 of file manifold.h.


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