Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Public Attributes | Private Member Functions | Private 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 std::unique_ptr< Manifold< dim, spacedim > > clone () const override
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) 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
 
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 > get_new_point (const ArrayView< const Point< spacedim >> &vertices, const ArrayView< const double > &weights) const override
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold () override=default
 
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 (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) 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
 

Private Member Functions

std::pair< double, Tensor< 1, spacedim > > guess_new_point (const ArrayView< const Tensor< 1, spacedim >> &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights) const
 
Point< spacedim > get_new_point (const ArrayView< const Tensor< 1, spacedim >> &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights, const Point< spacedim > &candidate_point) const
 
virtual void get_new_points (const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights, ArrayView< Point< spacedim >> new_points) const
 

Private Attributes

const PolarManifold< spacedim > polar_manifold
 

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::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string 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. An example for the use of this function would be in the description of a hyper-shell or hyper-ball geometry, for example after creating a coarse mesh using GridGenerator::hyper_ball(). (However, it is worth mentioning that generating a good mesh for a disk or ball is complicated and requires addition steps. See the "Possibilities for extensions" section of step-6 for an extensive discussion of how one would construct such meshes and what one needs to do for it.)

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 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 connected (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., this Manifold connects points using geodesics. If more than two points are involved through a SphericalManifold::get_new_points() call, a so-called spherical average is used where the final point minimizes the weighted distance to all other points via geodesics.

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. It is advisable to combine this class with TransfiniteInterpolationManifold to ensure a smooth transition from a curved shape to the straight coordinate system in the center of the ball.

Author
Mauro Bardelloni, Luca Heltai, Daniel Arndt, 2016, 2017

Definition at line 223 of file manifold_lib.h.

Constructor & Destructor Documentation

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

Member Function Documentation

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

Make a clone of this Manifold object.

Implements Manifold< dim, spacedim >.

Definition at line 291 of file manifold_lib.cc.

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
overridevirtual

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

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

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

Reimplemented from Manifold< dim, spacedim >.

Definition at line 366 of file manifold_lib.cc.

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

Return the (normalized) normal vector at the point p.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 417 of file manifold_lib.cc.

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

Definition at line 478 of file manifold_lib.cc.

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

This function is optimized to perform on a collection of new points, by collecting operations that are not dependent on the weights outside of the loop over all new points.

The implementation does not allow for surrounding_points and new_points to point to the same array, so make sure to pass different objects into the function.

Reimplemented from Manifold< dim, spacedim >.

Definition at line 517 of file manifold_lib.cc.

template<int dim, int spacedim>
Point< spacedim > SphericalManifold< dim, spacedim >::get_new_point ( const ArrayView< const Point< spacedim >> &  vertices,
const ArrayView< const double > &  weights 
) const
overridevirtual

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

Reimplemented from Manifold< dim, spacedim >.

Definition at line 534 of file manifold_lib.cc.

template<int dim, int spacedim>
std::pair< double, Tensor< 1, spacedim > > SphericalManifold< dim, spacedim >::guess_new_point ( const ArrayView< const Tensor< 1, spacedim >> &  directions,
const ArrayView< const double > &  distances,
const ArrayView< const double > &  weights 
) const
private

Return a point on the spherical manifold which is intermediate with respect to the surrounding points. This function uses a linear average of the directions to find an estimated point. It returns a pair of radius and direction from the center point to the candidate point.

Definition at line 761 of file manifold_lib.cc.

template<int dim, int spacedim>
Point< spacedim > SphericalManifold< dim, spacedim >::get_new_point ( const ArrayView< const Tensor< 1, spacedim >> &  directions,
const ArrayView< const double > &  distances,
const ArrayView< const double > &  weights,
const Point< spacedim > &  candidate_point 
) const
private

Return a point on the spherical manifold which is intermediate with respect to the surrounding points. This function uses a candidate point as guess, and performs a Newton-style iteration to compute the correct point.

The main part of the implementation uses the ideas in the publication

Buss, Samuel R., and Jay P. Fillmore. "Spherical averages and applications to spherical splines and interpolation." ACM Transactions on Graphics (TOG) 20.2 (2001): 95-126.

and in particular the implementation provided at http://math.ucsd.edu/~sbuss/ResearchWeb/spheremean/

Definition at line 932 of file manifold_lib.cc.

template<int dim, int spacedim>
void SphericalManifold< dim, spacedim >::get_new_points ( const ArrayView< const Point< spacedim >> &  surrounding_points,
const ArrayView< const double > &  weights,
ArrayView< Point< spacedim >>  new_points 
) const
privatevirtual

Compute a new set of points that interpolate between the given points surrounding_points. weights is an array view with as many entries as surrounding_points.size() times new_points.size().

This function is optimized to perform on a collection of new points, by collecting operations that are not dependent on the weights outside of the loop over all new points.

The implementation does not allow for surrounding_points and new_points to point to the same array, so make sure to pass different objects into the function.

Definition at line 552 of file manifold_lib.cc.

Member Data Documentation

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

The center of the spherical coordinate system.

Definition at line 304 of file manifold_lib.h.

template<int dim, int spacedim = dim>
const PolarManifold<spacedim> SphericalManifold< dim, spacedim >::polar_manifold
private

A manifold description to be used for get_new_point in 2D.

Definition at line 360 of file manifold_lib.h.


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