Reference documentation for deal.II version 8.5.1
Public Member Functions | Public Attributes | List of all members
SphericalManifold< dim, spacedim > Class Template Reference

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

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

Public Member Functions

 SphericalManifold (const Point< spacedim > center=Point< spacedim >())
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) const
 
virtual Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const
 
virtual Point< spacedim > get_new_point (const ::Quadrature< spacedim > &quadrature) const 1
 
virtual Point< spacedim > get_new_point (const std::vector< Point< spacedim > > &vertices, const std::vector< double > &weights) const
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold ()
 
virtual Point< spacedim > get_new_point (const Quadrature< spacedim > &quad) const 1
 
virtual void add_new_points (const std::vector< Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, std::vector< Point< spacedim > > &new_points) const
 
virtual Point< spacedim > project_to_manifold (const std::vector< Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const
 
virtual Point< spacedim > get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const
 
virtual Point< spacedim > get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
 
virtual Point< spacedim > get_new_point_on_hex (const typename Triangulation< dim, spacedim >::hex_iterator &hex) const
 
Point< spacedim > get_new_point_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face) const
 
Point< spacedim > get_new_point_on_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
virtual Tensor< 1, spacedim > normal_vector (const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const
 
virtual void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Public Attributes

const Point< spacedim > center
 

Additional Inherited Members

- Public Types inherited from Manifold< dim, spacedim >
typedef Tensor< 1, spacedim > FaceVertexNormals[GeometryInfo< dim >::vertices_per_face]
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Detailed Description

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

Manifold description for a spherical space coordinate system.

You can use this Manifold object to describe any sphere, circle, hypersphere or hyperdisc in two or three dimensions. This manifold can be used as a co-dimension one manifold descriptor of a spherical surface embedded in a higher dimensional space, or as a co-dimension zero manifold descriptor for a body with positive volume, provided that the center of the spherical space is excluded from the domain.

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 SphericalManifold<2,3> is identical to SphericalManifold<3,3>.

While PolarManifold reflects the usual notion of polar coordinates, it may not be suitable for domains that contain either the north or south poles. Consider for instance the pair of points \(x_1=(1,\pi/3,0)\) and \(x_2=(1,\pi/3,\pi)\) in polar coordinates (lying on the surface of a sphere with radius one, on a parallel at at height \(\pi/3\)). In this case connecting the points with a straight line in polar coordinates would take the long road around the globe, without passing through the north pole.

These two points would be connented (using a PolarManifold) by the curve

\begin{align*} s: [0,1] & \rightarrow & \mathbb S^3 \\ t & \mapsto & (1,\pi/3,0) + (0,0,t\pi) \end{align*}

This curve is not a geodesic on the sphere, and it is not how we would connect those two points. A better curve, would be the one passing through the North pole:

\[ s(t) = x_1 \cos(\alpha(t)) + \kappa \times x_1 \sin(\alpha(t)) + \kappa ( \kappa \cdot x_1) (1-\cos(\alpha(t))). \]

where \(\kappa = \frac{x_1 \times x_2}{\Vert x_1 \times x_2 \Vert}\) and \(\alpha(t) = t \cdot \arccos(x_1 \cdot x_2)\) for \(t\in[0,1]\). Indeed, this is a geodesic, and it is the natural choice when connecting points on the surface of the sphere. In the examples above, the PolarManifold class implements the first way of connecting two points on the surface of a sphere, while SphericalManifold implements the second way, i.e., if the codimension of the Manifold is one, than this Manifold connects points using geodesics. In all other cases it is a continuus extension of the codimension one case.

In particular, this class implements a Manifold that joins any two points in space by first projecting them onto the surface of a sphere with unit radius, then connecting them with a geodesic, and finally rescaling the final radius so that the resulting one is the weighted average of the starting radii. This Manifold is identical to PolarManifold in dimension two, while for dimension three it returns points that are more uniformly distributed on the sphere, and it is invariant with respect to rotations of the coordinate system, therefore avoiding the problems that PolarManifold has at the poles. Notice, in particular, that computing tangent vectors at the poles with a PolarManifold is not well defined, while it is perfectly fine with this class.

For mathematical reasons, it is impossible to construct a unique map of a sphere using only geodesic curves, and therefore, using this class with MappingManifold is discouraged. If you use this Manifold to describe the geometry of a sphere, you should use MappingQ as the underlying mapping, and not MappingManifold.

This Manifold can be used only on geometries where a ball with finite radius is removed from the center. Indeed, the center is a singular point for this manifold, and if you try to connect two points across the center, they would travel on spherical coordinates, avoiding the center.

The ideal geometry for this Manifold is an HyperShell. If you plan to use this Manifold on a HyperBall, you have to make sure you do not attach this Manifold to the cell containing the center.

Author
Mauro Bardelloni, Luca Heltai, 2016

Definition at line 200 of file manifold_lib.h.

Constructor & Destructor Documentation

◆ SphericalManifold()

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

The Constructor takes the center of the spherical coordinates.

Definition at line 167 of file manifold_lib.cc.

Member Function Documentation

◆ get_intermediate_point()

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

Given any two points in space, first project them on the surface of a sphere with unit radius, then connect them with a geodesic and find the intermediate point, and finally rescale the final radius so that the resulting one is the convex combination of the starting radii.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 174 of file manifold_lib.cc.

◆ get_tangent_vector()

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

Compute the derivative of the get_intermediate_point() function with parameter w equal to zero.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 233 of file manifold_lib.cc.

◆ get_new_point() [1/2]

template<int dim, int spacedim = dim>
virtual Point<spacedim> SphericalManifold< dim, spacedim >::get_new_point ( const ::Quadrature< spacedim > &  quadrature) const
virtual

Return a point on the spherical manifold which is intermediate with respect to the surrounding points.

Deprecated:
Use the other function that takes points and weights separately instead.

◆ get_new_point() [2/2]

template<int dim, int spacedim>
Point< spacedim > SphericalManifold< dim, spacedim >::get_new_point ( const std::vector< Point< spacedim > > &  vertices,
const std::vector< double > &  weights 
) const
virtual

Return a point on the spherical manifold which is intermediate with respect to the surrounding points.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 275 of file manifold_lib.cc.

Member Data Documentation

◆ center

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

The center of the spherical coordinate system.

Definition at line 252 of file manifold_lib.h.


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