Reference documentation for deal.II version 9.0.0
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
ConeBoundary< dim > Class Template Reference

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

Inheritance diagram for ConeBoundary< dim >:
[legend]

Public Member Functions

 ConeBoundary (const double radius_0, const double radius_1, const Point< dim > x_0, const Point< dim > x_1)
 
virtual std::unique_ptr< Manifold< dim, dim > > clone () const override
 
double get_radius (const Point< dim > x) const
 
virtual Point< dim > get_new_point_on_line (const typename Triangulation< dim >::line_iterator &line) const override
 
virtual Point< dim > get_new_point_on_quad (const typename Triangulation< dim >::quad_iterator &quad) const override
 
virtual void get_intermediate_points_on_line (const typename Triangulation< dim >::line_iterator &line, std::vector< Point< dim > > &points) const override
 
virtual void get_intermediate_points_on_quad (const typename Triangulation< dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const override
 
virtual void get_normals_at_vertices (const typename Triangulation< dim >::face_iterator &face, typename Boundary< dim >::FaceVertexNormals &face_vertex_normals) const override
 
virtual Tensor< 1, dim > normal_vector (const typename Triangulation< dim >::face_iterator &face, const Point< dim > &p) const override
 
- Public Member Functions inherited from StraightBoundary< dim >
 StraightBoundary ()
 
virtual Point< dim > get_new_point_on_line (const typename Triangulation< dim, dim >::line_iterator &line) const override
 
virtual Point< dim > get_new_point_on_quad (const typename Triangulation< dim, dim >::quad_iterator &quad) const override
 
virtual void get_intermediate_points_on_line (const typename Triangulation< dim, dim >::line_iterator &line, std::vector< Point< dim > > &points) const override
 
virtual void get_intermediate_points_on_quad (const typename Triangulation< dim, dim >::quad_iterator &quad, std::vector< Point< dim > > &points) const override
 
virtual Tensor< 1, dim > normal_vector (const typename Triangulation< dim, dim >::face_iterator &face, const Point< dim > &p) const override
 
virtual void get_normals_at_vertices (const typename Triangulation< dim, dim >::face_iterator &face, typename Boundary< dim, dim >::FaceVertexNormals &face_vertex_normals) const override
 
virtual Point< dim > project_to_surface (const typename Triangulation< dim, dim >::line_iterator &line, const Point< dim > &candidate) const override
 
virtual Point< dim > project_to_surface (const typename Triangulation< dim, dim >::quad_iterator &quad, const Point< dim > &candidate) const override
 
virtual Point< dim > project_to_surface (const typename Triangulation< dim, dim >::hex_iterator &hex, const Point< dim > &candidate) const override
 
- Public Member Functions inherited from Boundary< dim, dim >
virtual ~Boundary () override=default
 
virtual void get_intermediate_points_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const
 
virtual void get_intermediate_points_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
 
void get_intermediate_points_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face, std::vector< Point< spacedim > > &points) const
 
virtual Point< spacedim > project_to_surface (const typename Triangulation< dim, spacedim >::line_iterator &line, const Point< spacedim > &candidate) const
 
virtual Point< spacedim > project_to_surface (const typename Triangulation< dim, spacedim >::quad_iterator &quad, const Point< spacedim > &candidate) const
 
virtual Point< spacedim > project_to_surface (const typename Triangulation< dim, spacedim >::hex_iterator &hex, const Point< spacedim > &candidate) const
 
- Public Member Functions inherited from FlatManifold< dim, spacedim >
 FlatManifold (const Tensor< 1, spacedim > &periodicity=Tensor< 1, spacedim >(), const double tolerance=1e-10)
 
virtual Point< spacedim > get_new_point (const ArrayView< const Point< spacedim >> &surrounding_points, const ArrayView< const double > &weights) 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 > project_to_manifold (const ArrayView< const Point< spacedim >> &points, const Point< spacedim > &candidate) 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
 
const Tensor< 1, spacedim > & get_periodicity () const
 
- Public Member Functions inherited from Manifold< dim, spacedim >
virtual ~Manifold ()=default
 
virtual Point< spacedim > get_intermediate_point (const Point< spacedim > &p1, const Point< spacedim > &p2, const double w) 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
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Protected Attributes

const double radius_0
 
const double radius_1
 
const Point< dim > x_0
 
const Point< dim > x_1
 

Private Member Functions

void get_intermediate_points_between_points (const Point< dim > &p0, const Point< dim > &p1, std::vector< Point< dim > > &points) const
 

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, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Protected Member Functions inherited from Boundary< dim, dim >
const std::vector< Point< 1 > > & get_line_support_points (const unsigned int n_intermediate_points) const
 

Detailed Description

template<int dim>
class ConeBoundary< dim >

Boundary object for the hull of a (truncated) cone with two different radii at the two ends. If one radius is chosen to be 0 the object describes the boundary of a cone. In three dimensions, points are projected on an arbitrarily oriented (truncated) cone described by the two endpoints and the corresponding radii. Similar to HyperBallBoundary, new points are projected by dividing the straight line between the old two points and adjusting the radius from the axis.

This class is derived from StraightBoundary rather than from Boundary, which would seem natural, since this way we can use the StraightBoundary::in_between() function.

As an example of use, consider the following code snippet:

Triangulation<dim> triangulation;
Point<dim> p1, p2;
p1[0] = -1;
p2[0] = 1;
const ConeBoundary<dim> boundary (1, 0.5, p1, p2);
triangulation.set_boundary (0, boundary);
triangulation.refine_global (2);

This will produce the following meshes after the two refinements we perform, in 2d and 3d, respectively:

cone_2d.png
cone_3d.png
Deprecated:
The boundary classes have been deprecated in favor of the more general Manifold classes. For this class CylindricalManifold is an appropriate replacement for spacedim==3, else StraightManifold can be used.
Author
Markus Bürg, 2009

Definition at line 209 of file tria_boundary_lib.h.

Constructor & Destructor Documentation

◆ ConeBoundary()

template<int dim>
ConeBoundary< dim >::ConeBoundary ( const double  radius_0,
const double  radius_1,
const Point< dim >  x_0,
const Point< dim >  x_1 
)

Constructor. Here the boundary object is constructed. The points x_0 and x_1 describe the starting and ending points of the axis of the (truncated) cone. radius_0 denotes the radius corresponding to x_0 and radius_1 the one corresponding to x_1.

Definition at line 286 of file tria_boundary_lib.cc.

Member Function Documentation

◆ clone()

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

Clone this Boundary object.

Reimplemented from StraightBoundary< dim >.

Definition at line 300 of file tria_boundary_lib.cc.

◆ get_radius()

template<int dim>
double ConeBoundary< dim >::get_radius ( const Point< dim >  x) const

Return the radius of the (truncated) cone at given point x on the axis.

Definition at line 308 of file tria_boundary_lib.cc.

◆ get_new_point_on_line()

template<int dim>
Point< dim > ConeBoundary< dim >::get_new_point_on_line ( const typename Triangulation< dim >::line_iterator &  line) const
overridevirtual

Refer to the general documentation of this class and the documentation of the base class.

Definition at line 352 of file tria_boundary_lib.cc.

◆ get_new_point_on_quad()

template<int dim>
Point< dim > ConeBoundary< dim >::get_new_point_on_quad ( const typename Triangulation< dim >::quad_iterator &  quad) const
overridevirtual

Refer to the general documentation of this class and the documentation of the base class.

Definition at line 391 of file tria_boundary_lib.cc.

◆ get_intermediate_points_on_line()

template<int dim>
void ConeBoundary< dim >::get_intermediate_points_on_line ( const typename Triangulation< dim >::line_iterator &  line,
std::vector< Point< dim > > &  points 
) const
overridevirtual

Refer to the general documentation of this class and the documentation of the base class.

Calls get_intermediate_points_between_points.

Definition at line 403 of file tria_boundary_lib.cc.

◆ get_intermediate_points_on_quad()

template<int dim>
void ConeBoundary< dim >::get_intermediate_points_on_quad ( const typename Triangulation< dim >::quad_iterator &  quad,
std::vector< Point< dim > > &  points 
) const
overridevirtual

Refer to the general documentation of this class and the documentation of the base class.

Only implemented for dim=3 and for points.size()==1.

Definition at line 454 of file tria_boundary_lib.cc.

◆ get_normals_at_vertices()

template<int dim>
void ConeBoundary< dim >::get_normals_at_vertices ( const typename Triangulation< dim >::face_iterator &  face,
typename Boundary< dim >::FaceVertexNormals face_vertex_normals 
) const
overridevirtual

Compute the normals to the boundary at the vertices of the given face.

Refer to the general documentation of this class and the documentation of the base class.

Definition at line 477 of file tria_boundary_lib.cc.

◆ normal_vector()

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

Implementation of the function declared in the base class.

Refer to the general documentation of this class and the documentation of the base class.

Definition at line 501 of file tria_boundary_lib.cc.

◆ get_intermediate_points_between_points()

template<int dim>
void ConeBoundary< dim >::get_intermediate_points_between_points ( const Point< dim > &  p0,
const Point< dim > &  p1,
std::vector< Point< dim > > &  points 
) const
private

Called by get_intermediate_points_on_line and by get_intermediate_points_on_quad.

Refer to the general documentation of get_intermediate_points_on_line in the documentation of the base class.

Definition at line 322 of file tria_boundary_lib.cc.

Member Data Documentation

◆ radius_0

template<int dim>
const double ConeBoundary< dim >::radius_0
protected

First radius of the (truncated) cone.

Definition at line 300 of file tria_boundary_lib.h.

◆ radius_1

template<int dim>
const double ConeBoundary< dim >::radius_1
protected

Second radius of the (truncated) cone.

Definition at line 305 of file tria_boundary_lib.h.

◆ x_0

template<int dim>
const Point<dim> ConeBoundary< dim >::x_0
protected

Starting point of the axis.

Definition at line 310 of file tria_boundary_lib.h.

◆ x_1

template<int dim>
const Point<dim> ConeBoundary< dim >::x_1
protected

Ending point of the axis.

Definition at line 315 of file tria_boundary_lib.h.


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