Boundary< dim, spacedim > Class Template Reference
[Boundary description for triangulations]

Inheritance diagram for Boundary< dim, spacedim >:

Inheritance graph
[legend]

List of all members.

Public Types

typedef Tensor< 1, spacedim > FaceVertexNormals [GeometryInfo< dim >::vertices_per_face]

Public Member Functions

virtual ~Boundary ()
virtual Point< spacedim > get_new_point_on_line (const typename Triangulation< dim, spacedim >::line_iterator &line) const =0
virtual Point< spacedim > get_new_point_on_quad (const typename Triangulation< dim, spacedim >::quad_iterator &quad) const
Point< spacedim > get_new_point_on_face (const typename Triangulation< dim, spacedim >::face_iterator &face) const
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 void get_normals_at_vertices (const typename Triangulation< dim, spacedim >::face_iterator &face, FaceVertexNormals &face_vertex_normals) const


Detailed Description

template<int dim, int spacedim = dim>
class Boundary< dim, spacedim >

This class is used to represent a boundary to a triangulation. When a triangulation creates a new vertex on the boundary of the domain, it determines the new vertex' coordinates through the following code (here in two dimensions):
 *     ...
 *     Point<2> new_vertex = boundary.get_new_point_on_line (line);
 *     ...
 *   
line denotes the line at the boundary that shall be refined and for which we seek the common point of the two child lines.

In 3D, a new vertex may be placed on the middle of a line or on the middle of a side. Respectively, the library calls

 *     ...
 *     Point<3> new_line_vertices[4]
 *       = { boundary.get_new_point_on_line (face->line(0)),
 *           boundary.get_new_point_on_line (face->line(1)),
 *           boundary.get_new_point_on_line (face->line(2)),
 *           boundary.get_new_point_on_line (face->line(3))  };
 *     ...
 *   
to get the four midpoints of the lines bounding the quad at the boundary, and after that
 *     ...
 *     Point<3> new_quad_vertex = boundary.get_new_point_on_quad (face);
 *     ...
 *   
to get the midpoint of the face. It is guaranteed that this order (first lines, then faces) holds, so you can use information from the children of the four lines of a face, since these already exist at the time the midpoint of the face is to be computed.

Since iterators are passed to the functions, you may use information about boundary indicators and the like, as well as all other information provided by these objects.

There are specialisations, StraightBoundary, which places the new point right into the middle of the given points, and HyperBallBoundary creating a hyperball with given radius around a given center point.

Author:
Wolfgang Bangerth, 1999, 2001, Ralf Hartmann, 2001, 2008

Member Typedef Documentation

template<int dim, int spacedim = dim>
typedef Tensor<1,spacedim> Boundary< dim, spacedim >::FaceVertexNormals[GeometryInfo< dim >::vertices_per_face]

Type keeping information about the normals at the vertices of a face of a cell. Thus, there are GeometryInfo<dim>::vertices_per_face normal vectors, that define the tangent spaces of the boundary at the vertices. Note that the vectors stored in this object are not required to be normalized, nor to actually point outward, as one often will only want to check for orthogonality to define the tangent plane; if a function requires the normals to be normalized, then it must do so itself.

For obvious reasons, this type is not useful in 1d.


Constructor & Destructor Documentation

template<int dim, int spacedim = dim>
virtual Boundary< dim, spacedim >::~Boundary (  )  [virtual]

Destructor. Does nothing here, but needs to be declared to make it virtual.


Member Function Documentation

template<int dim, int spacedim = dim>
virtual Point<spacedim> Boundary< dim, spacedim >::get_new_point_on_line ( const typename Triangulation< dim, spacedim >::line_iterator &  line  )  const [pure virtual]

Return the point which shall become the new middle vertex of the two children of a regular line. In 2D, this line is a line at the boundary, while in 3d, it is bounding a face at the boundary (the lines therefore is also on the boundary).

Implemented in StraightBoundary< dim, spacedim >, HyperBallBoundary< dim, spacedim >, StraightBoundary< dim, spacedim >, StraightBoundary< spacedim, spacedim >, and StraightBoundary< dim, dim >.

template<int dim, int spacedim = dim>
virtual Point<spacedim> Boundary< dim, spacedim >::get_new_point_on_quad ( const typename Triangulation< dim, spacedim >::quad_iterator &  quad  )  const [virtual]

Return the point which shall become the common point of the four children of a quad at the boundary in three or more spatial dimensions. This function therefore is only useful in at least three dimensions and should not be called for lower dimensions.

This function is called after the four lines bounding the given quad are refined, so you may want to use the information provided by quad->line(i)->child(j), i=0...3, j=0,1.

Because in 2D, this function is not needed, it is not made pure virtual, to avoid the need to overload it. The default implementation throws an error in any case, however.

Reimplemented in StraightBoundary< dim, spacedim >, HyperBallBoundary< dim, spacedim >, StraightBoundary< dim, spacedim >, StraightBoundary< spacedim, spacedim >, and StraightBoundary< dim, dim >.

template<int dim, int spacedim = dim>
Point<spacedim> Boundary< dim, spacedim >::get_new_point_on_face ( const typename Triangulation< dim, spacedim >::face_iterator &  face  )  const

Depending on dim=2 or dim=3 this function calls the get_new_point_on_line or the get_new_point_on_quad function. It throws an exception for dim=1. This wrapper allows dimension independent programming.

template<int dim, int spacedim = dim>
virtual void Boundary< dim, spacedim >::get_intermediate_points_on_line ( const typename Triangulation< dim, spacedim >::line_iterator &  line,
std::vector< Point< spacedim > > &  points 
) const [virtual]

Return equally spaced intermediate points on a line.

The number of points requested is given by the size of the vector points. It is the task of the derived classes to arrange the points in approximately equal distances.

This function is called by the MappingQ class. This happens on each face line of a cells that has got at least one boundary line.

As this function is not needed for MappingQ1, it is not made pure virtual, to avoid the need to overload it. The default implementation throws an error in any case, however.

Reimplemented in StraightBoundary< dim, spacedim >, HyperBallBoundary< dim, spacedim >, StraightBoundary< dim, spacedim >, StraightBoundary< spacedim, spacedim >, and StraightBoundary< dim, dim >.

template<int dim, int spacedim = dim>
virtual void Boundary< dim, spacedim >::get_intermediate_points_on_quad ( const typename Triangulation< dim, spacedim >::quad_iterator &  quad,
std::vector< Point< spacedim > > &  points 
) const [virtual]

Return equally spaced intermediate points on a boundary quad.

The number of points requested is given by the size of the vector points. It is required that this number is a square of another integer, i.e. n=points.size()=m*m. It is the task of the derived classes to arrange the points such they split the quad into (m+1)(m+1) approximately equal-sized subquads.

This function is called by the MappingQ<3> class. This happens each face quad of cells in 3d that has got at least one boundary face quad.

As this function is not needed for MappingQ1, it is not made pure virtual, to avoid the need to overload it. The default implementation throws an error in any case, however.

Reimplemented in StraightBoundary< dim, spacedim >, HyperBallBoundary< dim, spacedim >, StraightBoundary< dim, spacedim >, StraightBoundary< spacedim, spacedim >, and StraightBoundary< dim, dim >.

template<int dim, int spacedim = dim>
void Boundary< dim, spacedim >::get_intermediate_points_on_face ( const typename Triangulation< dim, spacedim >::face_iterator &  face,
std::vector< Point< spacedim > > &  points 
) const

Depending on dim=2 or dim=3 this function calls the get_intermediate_points_on_line or the get_intermediate_points_on_quad function. It throws an exception for dim=1. This wrapper allows dimension independent programming.

template<int dim, int spacedim = dim>
virtual void Boundary< dim, spacedim >::get_normals_at_vertices ( const typename Triangulation< dim, spacedim >::face_iterator &  face,
FaceVertexNormals face_vertex_normals 
) const [virtual]

Compute the normal vectors to the boundary at each vertex of the given face. It is not required that the normal vectors be normed somehow. Neither is it required that the normals actually point outward.

This function is needed to compute data for C1 mappings. The default implementation is to throw an error, so you need not overload this function in case you do not intend to use C1 mappings.

Note that when computing normal vectors at a vertex where the boundary is not differentiable, you have to make sure that you compute the one-sided limits, i.e. limit with respect to points inside the given face.


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

deal.II documentation generated on Sat Aug 15 16:51:42 2009 by doxygen 1.5.9