Reference documentation for deal.II version 9.0.0
|
#include <deal.II/grid/tria_boundary_lib.h>
Public Member Functions | |
CylinderBoundary (const double radius=1.0, const unsigned int axis=0) | |
CylinderBoundary (const double radius, const Point< spacedim > &direction, const Point< spacedim > &point_on_axis) | |
virtual std::unique_ptr< Manifold< dim, spacedim > > | clone () const override |
virtual Point< spacedim > | get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const override |
virtual Point< spacedim > | get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const override |
virtual void | get_intermediate_points_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line, std::vector< Point< spacedim > > &points) const override |
virtual void | get_intermediate_points_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const override |
virtual void | get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, typename Boundary< dim, spacedim >::FaceVertexNormals &face_vertex_normals) const override |
double | get_radius () const |
Public Member Functions inherited from StraightBoundary< dim, spacedim > | |
StraightBoundary () | |
virtual Tensor< 1, spacedim > | normal_vector (const typename Triangulation< dim, spacedim >::face_iterator &face, const Point< spacedim > &p) const override |
virtual Point< spacedim > | project_to_surface (const typename Triangulation< dim, spacedim >::line_iterator &line, const Point< spacedim > &candidate) const override |
virtual Point< spacedim > | project_to_surface (const typename Triangulation< dim, spacedim >::quad_iterator &quad, const Point< spacedim > &candidate) const override |
virtual Point< spacedim > | project_to_surface (const typename Triangulation< dim, spacedim >::hex_iterator &hex, const Point< spacedim > &candidate) const override |
Public Member Functions inherited from Boundary< dim, spacedim > | |
virtual | ~Boundary () override=default |
void | get_intermediate_points_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face, std::vector< Point< spacedim > > &points) 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 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_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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcRadiusNotSet () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Attributes | |
const double | radius |
const Point< spacedim > | direction |
const Point< spacedim > | point_on_axis |
Private Member Functions | |
void | get_intermediate_points_between_points (const Point< spacedim > &p0, const Point< spacedim > &p1, std::vector< Point< spacedim > > &points) const |
Static Private Member Functions | |
static Point< spacedim > | get_axis_vector (const unsigned int axis) |
Additional Inherited Members | |
Public Types inherited from Manifold< dim, spacedim > | |
typedef Tensor< 1, spacedim > | FaceVertexNormals[GeometryInfo< dim >::vertices_per_face] |
Protected Member Functions inherited from Boundary< dim, spacedim > | |
const std::vector< Point< 1 > > & | get_line_support_points (const unsigned int n_intermediate_points) const |
Boundary object for the hull of a cylinder. In three dimensions, points are projected on a circular tube 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. The radius of the tube can be given independently. 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 was developed to be used in conjunction with the cylinder
function of GridGenerator. It should be used for the hull of the cylinder only (boundary indicator 0). Its use is discussed in detail in the results section of step-49.
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.
Definition at line 54 of file tria_boundary_lib.h.
CylinderBoundary< dim, spacedim >::CylinderBoundary | ( | const double | radius = 1.0 , |
const unsigned int | axis = 0 |
||
) |
Constructor. Using default values for the constructor arguments yields a circular tube along the x-axis (axis=0
). Choose axis=1
or axis=2
for a tube along the y- or z-axis, respectively.
Definition at line 31 of file tria_boundary_lib.cc.
CylinderBoundary< dim, spacedim >::CylinderBoundary | ( | const double | radius, |
const Point< spacedim > & | direction, | ||
const Point< spacedim > & | point_on_axis | ||
) |
Constructor. If constructed with this constructor, the boundary 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.
Definition at line 41 of file tria_boundary_lib.cc.
|
overridevirtual |
Clone this Boundary object.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 54 of file tria_boundary_lib.cc.
|
overridevirtual |
Refer to the general documentation of this class and the documentation of the base class.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 77 of file tria_boundary_lib.cc.
|
overridevirtual |
Refer to the general documentation of this class and the documentation of the base class.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 147 of file tria_boundary_lib.cc.
|
overridevirtual |
Refer to the general documentation of this class and the documentation of the base class.
Calls get_intermediate_points_between_points
.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 157 of file tria_boundary_lib.cc.
|
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
.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 234 of file tria_boundary_lib.cc.
|
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.
Reimplemented from StraightBoundary< dim, spacedim >.
Definition at line 259 of file tria_boundary_lib.cc.
double CylinderBoundary< dim, spacedim >::get_radius | ( | ) | const |
Return the radius of the cylinder.
Definition at line 277 of file tria_boundary_lib.cc.
|
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 170 of file tria_boundary_lib.cc.
|
staticprivate |
Given a number for the axis, return a vector that denotes this direction.
Definition at line 63 of file tria_boundary_lib.cc.
|
protected |
Radius of the cylinder.
Definition at line 140 of file tria_boundary_lib.h.
|
protected |
The direction vector of the axis.
Definition at line 145 of file tria_boundary_lib.h.
|
protected |
An arbitrary point on the axis.
Definition at line 150 of file tria_boundary_lib.h.