ASPECT
Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
aspect::SphericalManifold< dim, spacedim > Class Template Reference

Inherits Manifold< dim, spacedim >.

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 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
 
void do_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
 

Detailed Description

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

The deal.II class SphericalManifold has a design flaw that made it impossible to derive from the class. This is fixed post-9.5, see https://github.com/dealii/dealii/pull/16242 and https://github.com/dealii/dealii/pull/16248, but we can't use deal.II 9.5 and earlier for this class. The current class here is therefore a copy of the fixed class.

Definition at line 1107 of file compat.h.

Constructor & Destructor Documentation

§ SphericalManifold()

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

Constructor.

Parameters
[in]centerThe center of the coordinate system. Defaults to the origin.

Member Function Documentation

§ clone()

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

Make a clone of this Manifold object.

Reimplemented in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ get_intermediate_point()

template<int dim, int spacedim = dim>
virtual Point<spacedim> aspect::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 in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ get_tangent_vector()

template<int dim, int spacedim = dim>
virtual Tensor<1, spacedim> aspect::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 in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ normal_vector()

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

§ get_normals_at_vertices()

template<int dim, int spacedim = dim>
virtual void aspect::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.

Reimplemented in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ get_new_points()

template<int dim, int spacedim = dim>
virtual void aspect::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 in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ get_new_point()

template<int dim, int spacedim = dim>
virtual Point<spacedim> aspect::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 in aspect::GeometryModel::internal::SphericalManifoldWithTopography< dim >.

§ guess_new_point()

template<int dim, int spacedim = dim>
std::pair<double, Tensor<1, spacedim> > aspect::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.

§ do_get_new_points()

template<int dim, int spacedim = dim>
void aspect::SphericalManifold< dim, spacedim >::do_get_new_points ( const ArrayView< const Point< spacedim >> &  surrounding_points,
const ArrayView< const double > &  weights,
ArrayView< Point< spacedim >>  new_points 
) const
private

This function provides an internal implementation of the get_new_points() interface.

It computes 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.

Member Data Documentation

§ center

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

The center of the spherical coordinate system.

Definition at line 1191 of file compat.h.

§ polar_manifold

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

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

Definition at line 1230 of file compat.h.


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