Loading [MathJax]/extensions/TeX/newcommand.js
 Reference documentation for deal.II version 9.1.1
\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Public Member Functions | Protected Attributes | Private Attributes | List of all members
CylindricalManifold< dim, spacedim > Class Template Reference

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

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

Public Member Functions

 CylindricalManifold (const unsigned int axis=0, const double tolerance=1e-10)
 
 CylindricalManifold (const Tensor< 1, spacedim > &direction, const Point< spacedim > &point_on_axis, const double tolerance=1e-10)
 
virtual std::unique_ptr< Manifold< dim, spacedim > > clone () const override
 
virtual Point< 3 > pull_back (const Point< spacedim > &space_point) const override
 
virtual Point< spacedim > push_forward (const Point< 3 > &chart_point) const override
 
virtual DerivativeForm< 1, 3, spacedim > push_forward_gradient (const Point< 3 > &chart_point) const override
 
virtual Point< spacedim > get_new_point (const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) const override
 
- Public Member Functions inherited from ChartManifold< dim, spacedim, 3 >
 ChartManifold (const Tensor< 1, chartdim > &periodicity=Tensor< 1, chartdim >())
 
virtual ~ChartManifold () override=default
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) 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 > push_forward (const Point< chartdim > &chart_point) const=0
 
virtual DerivativeForm< 1, chartdim, spacedim > push_forward_gradient (const Point< chartdim > &chart_point) const
 
virtual Tensor< 1, spacedim > get_tangent_vector (const Point< spacedim > &x1, const Point< spacedim > &x2) const override
 
const Tensor< 1, chartdim > & get_periodicity () const
 
- 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 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 &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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)
 

Protected Attributes

const Tensor< 1, spacedim > normal_direction
 
const Tensor< 1, spacedim > direction
 
const Point< spacedim > point_on_axis
 

Private Attributes

double tolerance
 

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 ::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 CylindricalManifold< dim, spacedim >

Cylindrical Manifold description. In three dimensions, points are transformed using a cylindrical coordinate system along the x-, y- or z-axis (when using the first constructor of this class), or an arbitrarily oriented cylinder described by the direction of its axis and a point located on the axis.

This class was developed to be used in conjunction with the cylinder or cylinder_shell functions of GridGenerator. This function will throw a run time exception whenever spacedim is not equal to three.

Author
Luca Heltai, Daniel Arndt, 2014, 2017

Definition at line 379 of file manifold_lib.h.

Constructor & Destructor Documentation

◆ CylindricalManifold() [1/2]

template<int dim, int spacedim>
CylindricalManifold< dim, spacedim >::CylindricalManifold ( const unsigned int  axis = 0,
const double  tolerance = 1e-10 
)

Constructor. Using default values for the constructor arguments yields a cylinder along the x-axis (axis=0). Choose axis=1 or axis=2 for a tube along the y- or z-axis, respectively. The tolerance value is used to determine if a point is on the axis.

Definition at line 988 of file manifold_lib.cc.

◆ CylindricalManifold() [2/2]

template<int dim, int spacedim>
CylindricalManifold< dim, spacedim >::CylindricalManifold ( const Tensor< 1, spacedim > &  direction,
const Point< spacedim > &  point_on_axis,
const double  tolerance = 1e-10 
)

Constructor. If constructed with this constructor, the manifold described is a cylinder with an axis that points in direction direction and goes through the given point_on_axis. The direction may be arbitrarily scaled, and the given point may be any point on the axis. The tolerance value is used to determine if a point is on the axis.

Definition at line 1003 of file manifold_lib.cc.

Member Function Documentation

◆ clone()

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

Make a clone of this Manifold object.

Implements Manifold< dim, spacedim >.

Definition at line 1023 of file manifold_lib.cc.

◆ pull_back()

template<int dim, int spacedim>
Point< 3 > CylindricalManifold< dim, spacedim >::pull_back ( const Point< spacedim > &  space_point) const
overridevirtual

Compute the cylindrical coordinates (r, \phi, \lambda) for the given space point where r denotes the distance from the axis, \phi the angle between the given point and the computed normal direction, and \lambda the axial position.

Implements ChartManifold< dim, spacedim, 3 >.

Definition at line 1062 of file manifold_lib.cc.

◆ push_forward()

template<int dim, int spacedim>
Point< spacedim > CylindricalManifold< dim, spacedim >::push_forward ( const Point< 3 > &  chart_point) const
overridevirtual

Compute the Cartesian coordinates for a chart point given in cylindrical coordinates (r, \phi, \lambda), where r denotes the distance from the axis, \phi the angle between the given point and the computed normal direction, and \lambda the axial position.

Definition at line 1088 of file manifold_lib.cc.

◆ push_forward_gradient()

template<int dim, int spacedim>
DerivativeForm< 1, 3, spacedim > CylindricalManifold< dim, spacedim >::push_forward_gradient ( const Point< 3 > &  chart_point) const
overridevirtual

Compute the derivatives of the mapping from cylindrical coordinates (r, \phi, \lambda) to cartesian coordinates where r denotes the distance from the axis, \phi the angle between the given point and the computed normal direction, and \lambda the axial position.

Definition at line 1109 of file manifold_lib.cc.

◆ get_new_point()

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

Compute new points on the CylindricalManifold. See the documentation of the base class for a detailed description of what this function does.

Reimplemented from ChartManifold< dim, spacedim, 3 >.

Definition at line 1033 of file manifold_lib.cc.

Member Data Documentation

◆ normal_direction

template<int dim, int spacedim = dim>
const Tensor<1, spacedim> CylindricalManifold< dim, spacedim >::normal_direction
protected

A vector orthogonal to the normal direction.

Definition at line 447 of file manifold_lib.h.

◆ direction

template<int dim, int spacedim = dim>
const Tensor<1, spacedim> CylindricalManifold< dim, spacedim >::direction
protected

The direction vector of the axis.

Definition at line 452 of file manifold_lib.h.

◆ point_on_axis

template<int dim, int spacedim = dim>
const Point<spacedim> CylindricalManifold< dim, spacedim >::point_on_axis
protected

An arbitrary point on the axis.

Definition at line 457 of file manifold_lib.h.

◆ tolerance

template<int dim, int spacedim = dim>
double CylindricalManifold< dim, spacedim >::tolerance
private

Relative tolerance to measure zero distances.

Definition at line 463 of file manifold_lib.h.


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