Reference documentation for deal.II version 9.4.1
\(\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
Public Types | 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 Types

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

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
 
void get_normals_at_vertices (const Triangulation< 1 >::face_iterator &, Manifold< 1, 1 >::FaceVertexNormals &) const
 
void get_normals_at_vertices (const Triangulation< 1, 2 >::face_iterator &, Manifold< 1, 2 >::FaceVertexNormals &) const
 
Point< 3 > get_new_point (const ArrayView< const Tensor< 1, 3 > > &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights, const Point< 3 > &candidate_point) const
 
Point< 3 > get_new_point (const ArrayView< const Tensor< 1, 3 > > &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights, const Point< 3 > &candidate_point) const
 
Point< 3 > get_new_point (const ArrayView< const Tensor< 1, 3 > > &directions, const ArrayView< const double > &distances, const ArrayView< const double > &weights, const Point< 3 > &candidate_point) 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
 

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
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
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)
 
void check_no_subscribers () const noexcept
 
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)
 

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. (See also the extensive discussion in step-65.)

Definition at line 236 of file manifold_lib.h.

Member Typedef Documentation

◆ FaceVertexNormals

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

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

Member Function Documentation

◆ clone()

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

Make a clone of this Manifold object.

Implements Manifold< dim, spacedim >.

Definition at line 370 of file manifold_lib.cc.

◆ 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
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 379 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
overridevirtual

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

Reimplemented from Manifold< dim, spacedim >.

Definition at line 445 of file manifold_lib.cc.

◆ normal_vector() [1/3]

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

◆ get_normals_at_vertices() [1/6]

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

◆ get_new_points() [1/2]

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

◆ get_new_point() [1/5]

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

◆ guess_new_point()

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

◆ get_new_point() [2/5]

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

◆ get_new_points() [2/2]

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

◆ get_normals_at_vertices() [2/6]

void SphericalManifold< 1, 1 >::get_normals_at_vertices ( const Triangulation< 1 >::face_iterator &  ,
Manifold< 1, 1 >::FaceVertexNormals  
) const

Definition at line 526 of file manifold_lib.cc.

◆ get_normals_at_vertices() [3/6]

void SphericalManifold< 1, 2 >::get_normals_at_vertices ( const Triangulation< 1, 2 >::face_iterator &  ,
Manifold< 1, 2 >::FaceVertexNormals  
) const

Definition at line 537 of file manifold_lib.cc.

◆ get_new_point() [3/5]

Point< 3 > SphericalManifold< 1, 3 >::get_new_point ( const ArrayView< const Tensor< 1, 3 > > &  directions,
const ArrayView< const double > &  distances,
const ArrayView< const double > &  weights,
const Point< 3 > &  candidate_point 
) const

Definition at line 1004 of file manifold_lib.cc.

◆ get_new_point() [4/5]

Point< 3 > SphericalManifold< 2, 3 >::get_new_point ( const ArrayView< const Tensor< 1, 3 > > &  directions,
const ArrayView< const double > &  distances,
const ArrayView< const double > &  weights,
const Point< 3 > &  candidate_point 
) const

Definition at line 1017 of file manifold_lib.cc.

◆ get_new_point() [5/5]

Point< 3 > SphericalManifold< 3, 3 >::get_new_point ( const ArrayView< const Tensor< 1, 3 > > &  directions,
const ArrayView< const double > &  distances,
const ArrayView< const double > &  weights,
const Point< 3 > &  candidate_point 
) const

Definition at line 1030 of file manifold_lib.cc.

◆ project_to_manifold()

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

template<int dim, int spacedim = dim>
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]

template<int dim, int spacedim = dim>
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 418 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 429 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 440 of file manifold.cc.

◆ get_new_point_on_hex() [1/2]

template<int dim, int spacedim = dim>
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 462 of file manifold.cc.

◆ get_new_point_on_face() [1/4]

template<int dim, int spacedim = dim>
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 385 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 396 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 407 of file manifold.cc.

◆ get_new_point_on_cell()

template<int dim, int spacedim = dim>
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 144 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 165 of file manifold.cc.

◆ get_normals_at_vertices() [4/6]

template<int dim, int spacedim = dim>
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() [5/6]

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

Definition at line 250 of file manifold.cc.

◆ get_normals_at_vertices() [6/6]

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

Definition at line 272 of file manifold.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 317 of file manifold_lib.h.

◆ polar_manifold

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


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