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

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

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

Public Member Functions

 FunctionManifold (const Function< chartdim > &push_forward_function, const Function< spacedim > &pull_back_function, const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >(), const double tolerance=1e-10)
 
 FunctionManifold (const std::string push_forward_expression, const std::string pull_back_expression, const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >(), const typename FunctionParser< spacedim >::ConstMap=typename FunctionParser< spacedim >::ConstMap(), const std::string chart_vars=FunctionParser< chartdim >::default_variable_names(), const std::string space_vars=FunctionParser< spacedim >::default_variable_names(), const double tolerance=1e-10, const double h=1e-8)
 
 ~FunctionManifold ()
 
virtual Point< spacedim > push_forward (const Point< chartdim > &chart_point) const
 
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient (const Point< chartdim > &chart_point) const
 
virtual Point< chartdim > pull_back (const Point< spacedim > &space_point) const
 
- Public Member Functions inherited from ChartManifold< dim, spacedim, chartdim >
 ChartManifold (const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
 
virtual ~ChartManifold ()
 
virtual Point< spacedim > get_new_point (const Quadrature< spacedim > &quad) const 1
 
virtual Point< spacedim > get_new_point (const std::vector< Point< spacedim > > &surrounding_points, const std::vector< double > &weights) const
 
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 Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const
 
const Tensor< 1, chartdim > & get_periodicity () const
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold ()
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) 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)
 

Private Attributes

const FunctionParser< spacedim >::ConstMap const_map
 
SmartPointer< const Function< chartdim >, FunctionManifold< dim, spacedim, chartdim > > push_forward_function
 
SmartPointer< const Function< spacedim >, FunctionManifold< dim, spacedim, chartdim > > pull_back_function
 
const double tolerance
 
const bool owns_pointers
 

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, int chartdim = dim>
class FunctionManifold< dim, spacedim, chartdim >

Manifold description derived from ChartManifold, based on explicit Function<spacedim> and Function<chartdim> objects describing the push_forward() and pull_back() functions.

You can use this Manifold object to describe any arbitrary shape domain, as long as you can express it in terms of an invertible map, for which you provide both the forward expression, and the inverse expression.

In debug mode, a check is performed to verify that the transformations are actually one the inverse of the other.

Author
Luca Heltai, 2014

Definition at line 351 of file manifold_lib.h.

Constructor & Destructor Documentation

◆ FunctionManifold() [1/2]

template<int dim, int spacedim, int chartdim>
FunctionManifold< dim, spacedim, chartdim >::FunctionManifold ( const Function< chartdim > &  push_forward_function,
const Function< spacedim > &  pull_back_function,
const Tensor< 1, chartdim > &  periodicity = Tensor<1,chartdim>(),
const double  tolerance = 1e-10 
)

Explicit functions constructor. Takes a push_forward function of spacedim components, and a pull_back function of chartdim components. See the documentation of the base class ChartManifold for the meaning of the optional periodicity argument.

The tolerance argument is used in debug mode to actually check that the two functions are one the inverse of the other.

Definition at line 377 of file manifold_lib.cc.

◆ FunctionManifold() [2/2]

template<int dim, int spacedim, int chartdim>
FunctionManifold< dim, spacedim, chartdim >::FunctionManifold ( const std::string  push_forward_expression,
const std::string  pull_back_expression,
const Tensor< 1, chartdim > &  periodicity = Tensor<1,chartdim>(),
const typename FunctionParser< spacedim >::ConstMap  const_map = typename FunctionParser<spacedim>::ConstMap(),
const std::string  chart_vars = FunctionParser<chartdim>::default_variable_names(),
const std::string  space_vars = FunctionParser<spacedim>::default_variable_names(),
const double  tolerance = 1e-10,
const double  h = 1e-8 
)

Expressions constructor. Takes the expressions of the push_forward function of spacedim components, and of the pull_back function of chartdim components. See the documentation of the base class ChartManifold for the meaning of the optional periodicity argument.

The strings should be the readable by the default constructor of the FunctionParser classes. You can specify custom variable expressions with the last two optional arguments. If you don't, the default names are used, i.e., "x,y,z".

The tolerance argument is used in debug mode to actually check that the two functions are one the inverse of the other.

Definition at line 393 of file manifold_lib.cc.

◆ ~FunctionManifold()

template<int dim, int spacedim, int chartdim>
FunctionManifold< dim, spacedim, chartdim >::~FunctionManifold ( )

If needed, we delete the pointers we own.

Definition at line 415 of file manifold_lib.cc.

Member Function Documentation

◆ push_forward()

template<int dim, int spacedim, int chartdim>
Point< spacedim > FunctionManifold< dim, spacedim, chartdim >::push_forward ( const Point< chartdim > &  chart_point) const
virtual

Given a point in the chartdim coordinate system, uses the push_forward_function to compute the push_forward of points in chartdim space dimensions to spacedim space dimensions.

Implements ChartManifold< dim, spacedim, chartdim >.

Definition at line 431 of file manifold_lib.cc.

◆ push_forward_gradient()

template<int dim, int spacedim, int chartdim>
DerivativeForm< 1, chartdim, spacedim > FunctionManifold< dim, spacedim, chartdim >::push_forward_gradient ( const Point< chartdim > &  chart_point) const
virtual

Given a point in the chartdim dimensional Euclidean space, this method returns the derivatives of the function \(F\) that maps from the sub_manifold coordinate system to the Euclidean coordinate system. In other words, it is a matrix of size \(\text{spacedim}\times\text{chartdim}\).

This function is used in the computations required by the get_tangent_vector() function. The default implementation calls the get_gradient() method of the FunctionManifold::push_forward_function() member class. If you construct this object using the constructor that takes two string expression, then the default implementation of this method uses a finite difference scheme to compute the gradients(see the AutoDerivativeFunction() class for details), and you can specify the size of the spatial step size at construction time with the h parameter.

Refer to the general documentation of this class for more information.

Reimplemented from ChartManifold< dim, spacedim, chartdim >.

Definition at line 455 of file manifold_lib.cc.

◆ pull_back()

template<int dim, int spacedim, int chartdim>
Point< chartdim > FunctionManifold< dim, spacedim, chartdim >::pull_back ( const Point< spacedim > &  space_point) const
virtual

Given a point in the spacedim coordinate system, uses the pull_back_function to compute the pull_back of points in spacedim space dimensions to chartdim space dimensions.

Implements ChartManifold< dim, spacedim, chartdim >.

Definition at line 469 of file manifold_lib.cc.

Member Data Documentation

◆ const_map

template<int dim, int spacedim = dim, int chartdim = dim>
const FunctionParser<spacedim>::ConstMap FunctionManifold< dim, spacedim, chartdim >::const_map
private

Constants for the FunctionParser classes.

Definition at line 440 of file manifold_lib.h.

◆ push_forward_function

template<int dim, int spacedim = dim, int chartdim = dim>
SmartPointer<const Function<chartdim>, FunctionManifold<dim,spacedim,chartdim> > FunctionManifold< dim, spacedim, chartdim >::push_forward_function
private

Pointer to the push_forward function.

Definition at line 446 of file manifold_lib.h.

◆ pull_back_function

template<int dim, int spacedim = dim, int chartdim = dim>
SmartPointer<const Function<spacedim>, FunctionManifold<dim,spacedim,chartdim> > FunctionManifold< dim, spacedim, chartdim >::pull_back_function
private

Pointer to the pull_back function.

Definition at line 452 of file manifold_lib.h.

◆ tolerance

template<int dim, int spacedim = dim, int chartdim = dim>
const double FunctionManifold< dim, spacedim, chartdim >::tolerance
private

Relative tolerance. In debug mode, we check that the two functions provided at construction time are actually one the inverse of the other. This value is used as relative tolerance in this check.

Definition at line 459 of file manifold_lib.h.

◆ owns_pointers

template<int dim, int spacedim = dim, int chartdim = dim>
const bool FunctionManifold< dim, spacedim, chartdim >::owns_pointers
private

Check ownership of the smart pointers. Indicates whether this class is the owner of the objects pointed to by the previous two member variables. This value is set in the constructor of the class. If true, then the destructor will delete the function objects pointed to be the two pointers.

Definition at line 468 of file manifold_lib.h.


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