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

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

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

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
 
- Public Member Functions inherited from ChartManifold< dim, dim, dim >
 ChartManifold (const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
 
virtual ~ChartManifold () override=default
 
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
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold () override=default
 
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
 
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
 
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
 
virtual Point< spacedim > project_to_manifold (const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const
 
virtual Point< spacedim > get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const
 
virtual Point< spacedim > get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
 
virtual Point< spacedim > get_new_point_on_hex (const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
 
Point< spacedim > get_new_point_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face) const
 
Point< spacedim > get_new_point_on_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
virtual void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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)
 

Public Attributes

const Point< spacedim > center
 

Static Private Member Functions

static Tensor< 1, spacedim > get_periodicity ()
 

Additional Inherited Members

- Public Types inherited from Manifold< dim, spacedim >
using FaceVertexNormals = std::array< Tensor< 1, spacedim >, GeometryInfo< dim >::vertices_per_face >
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

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.

Author
Luca Heltai, Mauro Bardelloni, 2014-2016

Definition at line 63 of file manifold_lib.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 122 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 132 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 192 of file manifold_lib.cc.

◆ push_forward()

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 158 of file manifold_lib.cc.

◆ push_forward_gradient()

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 231 of file manifold_lib.cc.

◆ normal_vector()

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 320 of file manifold_lib.cc.

◆ get_periodicity()

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 141 of file manifold_lib.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 124 of file manifold_lib.h.


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