Reference documentation for deal.II version GIT relicensing-439-g5fda5c893d 2024-04-20 06:50:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Static Public Member Functions | Static Public Attributes | List of all members
GeometryInfo< dim > Struct Template Reference

#include <deal.II/base/geometry_info.h>

Static Public Member Functions

static std_cxx20::ranges::iota_view< unsigned int, unsigned intface_indices ()
 
static std_cxx20::ranges::iota_view< unsigned int, unsigned intvertex_indices ()
 
static unsigned int n_children (const RefinementCase< dim > &refinement_case)
 
static unsigned int n_subfaces (const internal::SubfaceCase< dim > &subface_case)
 
static double subface_ratio (const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
 
static RefinementCase< dim - 1 > face_refinement_case (const RefinementCase< dim > &cell_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static RefinementCase< dim > min_cell_refinement_case_for_face_refinement (const RefinementCase< dim - 1 > &face_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static RefinementCase< 1 > line_refinement_case (const RefinementCase< dim > &cell_refinement_case, const unsigned int line_no)
 
static RefinementCase< dim > min_cell_refinement_case_for_line_refinement (const unsigned int line_no)
 
static unsigned int child_cell_on_face (const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
 
static unsigned int line_to_cell_vertices (const unsigned int line, const unsigned int vertex)
 
static unsigned int face_to_cell_vertices (const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static unsigned int face_to_cell_lines (const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static unsigned int standard_to_real_face_vertex (const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static unsigned int real_to_standard_face_vertex (const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static unsigned int standard_to_real_face_line (const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static unsigned int standard_to_real_line_vertex (const unsigned int vertex, const bool line_orientation=true)
 
static std::array< unsigned int, 2 > standard_quad_vertex_to_line_vertex_index (const unsigned int vertex)
 
static std::array< unsigned int, 2 > standard_hex_vertex_to_quad_vertex_index (const unsigned int vertex)
 
static std::array< unsigned int, 2 > standard_hex_line_to_quad_line_index (const unsigned int line)
 
static unsigned int real_to_standard_face_line (const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
 
static Point< dim > unit_cell_vertex (const unsigned int vertex)
 
static unsigned int child_cell_from_point (const Point< dim > &p)
 
static Point< dim > cell_to_child_coordinates (const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
 
static Point< dim > child_to_cell_coordinates (const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
 
static bool is_inside_unit_cell (const Point< dim > &p)
 
static bool is_inside_unit_cell (const Point< dim > &p, const double eps)
 
template<typename Number = double>
static Point< dim, Number > project_to_unit_cell (const Point< dim, Number > &p)
 
static double distance_to_unit_cell (const Point< dim > &p)
 
static double d_linear_shape_function (const Point< dim > &xi, const unsigned int i)
 
static Tensor< 1, dim > d_linear_shape_function_gradient (const Point< dim > &xi, const unsigned int i)
 
template<int spacedim>
static void alternating_form_at_vertices (const Point< spacedim >(&vertices)[vertices_per_cell], Tensor< spacedim - dim, spacedim >(&forms)[vertices_per_cell])
 
static ::ExceptionBaseExcInvalidCoordinate (double arg1)
 
static ::ExceptionBaseExcInvalidSubface (int arg1, int arg2, int arg3)
 

Static Public Attributes

static constexpr unsigned int max_children_per_cell = 1 << dim
 
static constexpr unsigned int faces_per_cell = 2 * dim
 
static constexpr unsigned int max_children_per_face
 
static constexpr unsigned int vertices_per_cell = 1 << dim
 
static constexpr unsigned int vertices_per_face
 
static constexpr unsigned int lines_per_face
 
static constexpr unsigned int quads_per_face
 
static constexpr unsigned int lines_per_cell
 
static constexpr unsigned int quads_per_cell
 
static constexpr unsigned int hexes_per_cell
 
static constexpr std::array< unsigned int, vertices_per_cellucd_to_deal
 
static constexpr std::array< unsigned int, vertices_per_celldx_to_deal
 
static constexpr ndarray< unsigned int, vertices_per_cell, dim > vertex_to_face
 
static constexpr std::array< unsigned int, faces_per_cellunit_normal_direction
 
static constexpr std::array< int, faces_per_cellunit_normal_orientation
 
static constexpr std::array< Tensor< 1, dim >, faces_per_cellunit_normal_vector
 
static constexpr ndarray< Tensor< 1, dim >, faces_per_cell, dim - 1 > unit_tangential_vectors
 
static constexpr std::array< unsigned int, faces_per_cellopposite_face
 

Detailed Description

template<int dim>
struct GeometryInfo< dim >

This class provides dimension independent information to all topological structures that make up the unit, or reference cell. That said, this class only describes information about hypercube reference cells (i.e., lines, quadrilaterals, or hexahedra), which historically were the only kinds of cells supported by deal.II. This is no longer the case today, and consequently this class has been superseded by the ReferenceCell class – see there for more information. The rest of this class's documentation is therefore partly historical.

It is the one central point in the library where information about the numbering of vertices, lines, or faces of the reference cell is collected. Consequently, the information of this class is used extensively in the geometric description of Triangulation objects, as well as in various other parts of the code. In particular, it also serves as the focus of writing code in a dimension independent way; for example, instead of writing a loop over vertices 0<=v<4 in 2d, one would write it as 0<=v<GeometryInfo<dim>::vertices_per_cell, thus allowing the code to work in 3d as well without changes.

The most frequently used parts of the class are its static members like vertices_per_cell, faces_per_cell, etc. However, the class also offers information about more abstract questions like the orientation of faces, etc. The following documentation gives a textual description of many of these concepts.

Implementation conventions for two spatial dimensions

From version 5.2 onwards deal.II is based on a numbering scheme that uses a lexicographic ordering (with x running fastest) wherever possible, hence trying to adopt a kind of 'canonical' ordering.

The ordering of vertices and faces (lines) in 2d is defined by

The resulting numbering of vertices and faces (lines) in 2d as well as the directions of lines is shown in the following.

*       3
*    2-->--3
*    |     |
*   0^     ^1
*    |     |
*    0-->--1
*        2
* 

Note that the orientation of lines has to be correct upon construction of a grid; however, it is automatically preserved upon refinement.

Further we define that child lines have the same direction as their parent, i.e. that line->child(0)->vertex(0)==line->vertex(0) and line->child(1)->vertex(1)==line->vertex(1). This also implies, that the first sub-line (line->child(0)) is the one at vertex(0) of the old line.

Similarly we define, that the four children of a quad are adjacent to the vertex with the same number of the old quad.

Note that information about several of these conventions can be extracted at run- or compile-time from the member functions and variables of the present class.

Coordinate systems

When explicit coordinates are required for points in a cell (e.g for quadrature formulae or the point of definition of trial functions), we define the following coordinate system for the unit cell:

*  y^   2-----3
*   |   |     |
*   |   |     |
*   |   |     |
*   |   0-----1
*   *------------>x
* 

Here, vertex 0 is the origin of the coordinate system, vertex 1 has coordinates (1,0), vertex 2 at (0,1) and vertex 3 at (1,1). The GeometryInfo<dim>::unit_cell_vertex() function can be used to query this information at run-time.

Implementation conventions for three spatial dimensions

By convention, we will use the following numbering conventions for vertices, lines and faces of hexahedra in three space dimensions. Before giving these conventions we declare the following sketch to be the standard way of drawing 3d pictures of hexahedra:

*                       *-------*        *-------*
*                      /|       |       /       /|
*                     / |       |      /       / |
*  z                 /  |       |     /       /  |
*  ^                *   |       |    *-------*   |
*  |   ^y           |   *-------*    |       |   *
*  |  /             |  /       /     |       |  /
*  | /              | /       /      |       | /
*  |/               |/       /       |       |/
*  *------>x        *-------*        *-------*
* 

The left part of the picture shows the left, bottom and back face of the cube, while the right one shall be the top, right and front face. You may recover the whole cube by moving the two parts together into one.

Note again that information about several of the following conventions can be extracted at run- or compile-time from the member functions and variables of the present class.

Vertices

The ordering of vertices in 3d is defined by the same rules as in the 2d case. In particular, the following is still true:

Hence, the vertices are numbered as follows

*       6-------7        6-------7
*      /|       |       /       /|
*     / |       |      /       / |
*    /  |       |     /       /  |
*   4   |       |    4-------5   |
*   |   2-------3    |       |   3
*   |  /       /     |       |  /
*   | /       /      |       | /
*   |/       /       |       |/
*   0-------1        0-------1
* 

We note, that first the vertices on the bottom face (z=0) are numbered exactly the same way as are the vertices on a quadrilateral. Then the vertices on the top face (z=1) are numbered similarly by moving the bottom face to the top. Again, the GeometryInfo<dim>::unit_cell_vertex() function can be used to query this information at run-time.

Lines

Here, the same holds as for the vertices:

The fact that edges (just as vertices and faces) are entities that are stored in their own right rather than constructed from cells each time they are needed, means that adjacent cells actually have pointers to edges that are thus shared between them. This implies that the convention that sets of parallel edges have parallel directions is not only a local condition. Before a list of cells is passed to an object of the Triangulation class for creation of a triangulation, you therefore have to make sure that cells are oriented in a compatible fashion, so that edge directions are globally according to above convention. However, the GridTools::consistently_order_cells() function can do this for you, by reorienting cells and edges of an arbitrary list of input cells that need not be already sorted.

Faces

The numbering of faces in 3d is defined by a rule analogous to 2d:

Therefore, the faces are numbered in the ordering: left, right, front, back, bottom and top face:

*       *-------*        *-------*
*      /|       |       /       /|
*     / |   3   |      /   5   / |
*    /  |       |     /       /  |
*   *   |       |    *-------*   |
*   | 0 *-------*    |       | 1 *
*   |  /       /     |       |  /
*   | /   4   /      |   2   | /
*   |/       /       |       |/
*   *-------*        *-------*
* 

The standard direction of the faces is such, that the induced 2d local coordinate system (x,y) implies (right hand rule) a normal in face normal direction, see N2a). In the following we show the local coordinate system and the numbering of face lines:

The face line numbers (0,1,2,3) correspond to following cell line numbers.

You can get these numbers using the GeometryInfo<3>::face_to_cell_lines() function.

The face normals can be deduced from the face orientation by applying the right hand side rule (x,y -> normal). We note, that in the standard orientation of faces in 2d, faces 0 and 2 have normals that point into the cell, and faces 1 and 3 have normals pointing outward. In 3d, faces 0, 2, and 4 have normals that point into the cell, while the normals of faces 1, 3, and 5 point outward. This information, again, can be queried from GeometryInfo<dim>::unit_normal_orientation.

However, it turns out that a significant number of 3d meshes cannot satisfy this convention. This is due to the fact that the face convention for one cell already implies something for the neighbor, since they share a common face and fixing it for the first cell also fixes the normal vectors of the opposite faces of both cells. It is easy to construct cases of loops of cells for which this leads to cases where we cannot find orientations for all faces that are consistent with this convention.

For this reason, above convention is only what we call the standard orientation. deal.II actually allows faces in 3d to have either the standard direction, or its opposite, in which case the lines that make up a cell would have reverted orders, and the above line equivalences would not hold any more. You can ask a cell whether a given face has standard orientation by calling cell->face_orientation(face_no): if the result is true, then the face has standard orientation, otherwise its normal vector is pointing the other direction. There are not very many places in application programs where you need this information actually, but a few places in the library make use of this. Note that in 2d, the result is always true. More information on the topic can be found in this glossary article.

In order to allow all kinds of meshes in 3d, including Moebius-loops, a face might even be rotated looking from one cell, whereas it is according to the standard looking at it from the neighboring cell sharing that particular face. In order to cope with this, two flags face_flip and face_rotation are available, to represent rotations by 180 and 90 degree, respectively. Setting both flags amounts to a rotation of 270 degrees (all counterclockwise). You can ask the cell for these flags like for the face_orientation. In order to enable rotated faces, even lines can deviate from their standard direction in 3d. This information is available as the line_orientation flag for cells and faces in 3d. Again, this is something that should be internal to the library and application program will probably never have to bother about it. For more information on this see also this glossary entry.

Children

The eight children of an isotropically refined cell are numbered according to the vertices they are adjacent to:

*       *----*----*        *----*----*
*      /| 6  |  7 |       / 6  /  7 /|
*     *6|    |    |      *----*----*7|
*    /| *----*----*     / 4  /  5 /| *
*   * |/|    |    |    *----*----* |/|
*   |4* | 2  |  3 |    | 4  |  5 |5*3|
*   |/|2*----*----*    |    |    |/| *
*   * |/ 2  /  3 /     *----*----* |/
*   |0*----*----*      |    |    |1*
*   |/0   /  1 /       | 0  |  1 |/
*   *----*----*        *----*----*
* 

Taking into account the orientation of the faces, the following children are adjacent to the respective faces:

You can get these numbers using the GeometryInfo<3>::child_cell_on_face() function. As each child is adjacent to the vertex with the same number these numbers are also given by the GeometryInfo<3>::face_to_cell_vertices() function.

Note that, again, the above list only holds for faces in their standard orientation. If a face is not in standard orientation, then the children at positions 1 and 2 (counting from 0 to 3) would be swapped. In fact, this is what the child_cell_on_face and the face_to_cell_vertices functions of GeometryInfo<3> do, when invoked with a face_orientation=false argument.

The information which child cell is at which position of which face is most often used when computing jump terms across faces with hanging nodes, using objects of type FESubfaceValues. Sitting on one cell, you would look at a face and figure out which child of the neighbor is sitting on a given subface between the present and the neighboring cell. To avoid having to query the standard orientation of the faces of the two cells every time in such cases, you should use a function call like cell->neighbor_child_on_subface(face_no,subface_no), which returns the correct result both in 2d (where face orientations are immaterial) and 3d (where it is necessary to use the face orientation as additional argument to GeometryInfo<3>::child_cell_on_face).

For anisotropic refinement, the child cells can not be numbered according to adjacent vertices, thus the following conventions are used:

*            RefinementCase<3>::cut_x
*
*       *----*----*        *----*----*
*      /|    |    |       /    /    /|
*     / |    |    |      / 0  /  1 / |
*    /  | 0  |  1 |     /    /    /  |
*   *   |    |    |    *----*----*   |
*   | 0 |    |    |    |    |    | 1 |
*   |   *----*----*    |    |    |   *
*   |  /    /    /     | 0  | 1  |  /
*   | / 0  /  1 /      |    |    | /
*   |/    /    /       |    |    |/
*   *----*----*        *----*----*
* 
*            RefinementCase<3>::cut_y
*
*       *---------*        *---------*
*      /|         |       /    1    /|
*     * |         |      *---------* |
*    /| |    1    |     /    0    /| |
*   * |1|         |    *---------* |1|
*   | | |         |    |         | | |
*   |0| *---------*    |         |0| *
*   | |/    1    /     |    0    | |/
*   | *---------*      |         | *
*   |/    0    /       |         |/
*   *---------*        *---------*
* 
*            RefinementCase<3>::cut_z
*
*       *---------*        *---------*
*      /|    1    |       /         /|
*     / |         |      /    1    / |
*    /  *---------*     /         /  *
*   * 1/|         |    *---------* 1/|
*   | / |    0    |    |    1    | / |
*   |/  *---------*    |         |/  *
*   * 0/         /     *---------* 0/
*   | /    0    /      |         | /
*   |/         /       |    0    |/
*   *---------*        *---------*
* 
*            RefinementCase<3>::cut_xy
*
*       *----*----*        *----*----*
*      /|    |    |       / 2  /  3 /|
*     * |    |    |      *----*----* |
*    /| | 2  |  3 |     / 0  /  1 /| |
*   * |2|    |    |    *----*----* |3|
*   | | |    |    |    |    |    | | |
*   |0| *----*----*    |    |    |1| *
*   | |/ 2  /  3 /     | 0  |  1 | |/
*   | *----*----*      |    |    | *
*   |/ 0  /  1 /       |    |    |/
*   *----*----*        *----*----*
* 
*            RefinementCase<3>::cut_xz
*
*       *----*----*        *----*----*
*      /| 1  |  3 |       /    /    /|
*     / |    |    |      / 1  /  3 / |
*    /  *----*----*     /    /    /  *
*   * 1/|    |    |    *----*----* 3/|
*   | / | 0  |  2 |    | 1  |  3 | / |
*   |/  *----*----*    |    |    |/  *
*   * 0/    /    /     *----*----* 2/
*   | / 0  /  2 /      |    |    | /
*   |/    /    /       | 0  |  2 |/
*   *----*----*        *----*----*
* 
*            RefinementCase<3>::cut_yz
*
*       *---------*        *---------*
*      /|    3    |       /    3    /|
*     * |         |      *---------* |
*    /|3*---------*     /    2    /|3*
*   * |/|         |    *---------* |/|
*   |2* |    1    |    |    2    |2* |
*   |/|1*---------*    |         |/|1*
*   * |/    1    /     *---------* |/
*   |0*---------*      |         |0*
*   |/    0    /       |    0    |/
*   *---------*        *---------*
* 

This information can also be obtained by the GeometryInfo<3>::child_cell_on_face function.

Coordinate systems

We define the following coordinate system for the explicit coordinates of the vertices of the unit cell:

*                       6-------7        6-------7
*                      /|       |       /       /|
*                     / |       |      /       / |
*  z                 /  |       |     /       /  |
*  ^                4   |       |    4-------5   |
*  |   ^y           |   2-------3    |       |   3
*  |  /             |  /       /     |       |  /
*  | /              | /       /      |       | /
*  |/               |/       /       |       |/
*  *------>x        0-------1        0-------1
* 

By the convention laid down as above, the vertices have the following coordinates (lexicographic, with x running fastest):

Note
Instantiations for this template are provided for dimensions 1,2,3,4, and there is a specialization for dim=0 (see the section on Template instantiations in the manual).

Definition at line 1985 of file geometry_info.h.

Member Function Documentation

◆ face_indices()

template<int dim>
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > GeometryInfo< dim >::face_indices ( )
static

Return an object that can be thought of as an array containing all indices from zero to faces_per_cell. This allows to write code using range-based for loops of the following kind:

for (auto &cell : triangulation.active_cell_iterators())
for (auto face_index : GeometryInfo<dim>::face_indices())
if (cell->face(face_index)->at_boundary())
... do something ...
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()

Here, we are looping over all faces of all cells, with face_index taking on all valid indices for faces (zero and one in 1d, zero through three in 2d, and zero through 5 in 3d).

See also
deal.II and Modern C++ standards

◆ vertex_indices()

template<int dim>
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > GeometryInfo< dim >::vertex_indices ( )
static

Return an object that can be thought of as an array containing all indices from zero to vertices_per_cell. This allows to write code using range-based for loops of the following kind:

for (auto &cell : triangulation.active_cell_iterators())
for (auto vertex_index : GeometryInfo<dim>::vertex_indices())
if (cell->vertex(vertex_index) satisfies some condition)
... do something ...
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()

Here, we are looping over all vertices of all cells, with vertex_index taking on all valid indices.

See also
deal.II and Modern C++ standards

◆ n_children()

template<int dim>
static unsigned int GeometryInfo< dim >::n_children ( const RefinementCase< dim > &  refinement_case)
static

Return the number of children of a cell (or face) refined with ref_case.

◆ n_subfaces()

template<int dim>
static unsigned int GeometryInfo< dim >::n_subfaces ( const internal::SubfaceCase< dim > &  subface_case)
static

Return the number of subfaces of a face refined according to internal::SubfaceCase face_ref_case.

◆ subface_ratio()

template<int dim>
static double GeometryInfo< dim >::subface_ratio ( const internal::SubfaceCase< dim > &  subface_case,
const unsigned int  subface_no 
)
static

Given a face on the reference element with a internal::SubfaceCase<dim> face_refinement_case this function returns the ratio between the area of the subface_no th subface and the area(=1) of the face.

E.g. for internal::SubfaceCase::cut_xy the ratio is 1/4 for each of the subfaces.

◆ face_refinement_case()

template<int dim>
static RefinementCase< dim - 1 > GeometryInfo< dim >::face_refinement_case ( const RefinementCase< dim > &  cell_refinement_case,
const unsigned int  face_no,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Given a cell refined with the RefinementCase cell_refinement_case return the SubfaceCase of the face_no th face.

◆ min_cell_refinement_case_for_face_refinement()

template<int dim>
static RefinementCase< dim > GeometryInfo< dim >::min_cell_refinement_case_for_face_refinement ( const RefinementCase< dim - 1 > &  face_refinement_case,
const unsigned int  face_no,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Given the SubfaceCase face_refinement_case of the face_no th face, return the smallest RefinementCase of the cell, which corresponds to that refinement of the face.

◆ line_refinement_case()

template<int dim>
static RefinementCase< 1 > GeometryInfo< dim >::line_refinement_case ( const RefinementCase< dim > &  cell_refinement_case,
const unsigned int  line_no 
)
static

Given a cell refined with the RefinementCase cell_refinement_case return the RefinementCase of the line_no th face.

◆ min_cell_refinement_case_for_line_refinement()

template<int dim>
static RefinementCase< dim > GeometryInfo< dim >::min_cell_refinement_case_for_line_refinement ( const unsigned int  line_no)
static

Return the minimal / smallest RefinementCase of the cell, which ensures refinement of line line_no.

◆ child_cell_on_face()

template<int dim>
static unsigned int GeometryInfo< dim >::child_cell_on_face ( const RefinementCase< dim > &  ref_case,
const unsigned int  face,
const unsigned int  subface,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false,
const RefinementCase< dim - 1 > &  face_refinement_case = RefinementCase< dim - 1 >::isotropic_refinement 
)
static

This field stores which child cells are adjacent to a certain face of the mother cell.

For example, in 2d the layout of a cell is as follows:

* .      3
* .   2-->--3
* .   |     |
* . 0 ^     ^ 1
* .   |     |
* .   0-->--1
* .      2
* 

Vertices and faces are indicated with their numbers, faces also with their directions.

Now, when refined, the layout is like this:

* *--*--*
* | 2|3 |
* *--*--*
* | 0|1 |
* *--*--*
* 

Thus, the child cells on face 0 are (ordered in the direction of the face) 0 and 2, on face 3 they are 2 and 3, etc.

For three spatial dimensions, the exact order of the children is laid down in the general documentation of this class.

Through the face_orientation, face_flip and face_rotation arguments this function handles faces oriented in the standard and non-standard orientation. face_orientation defaults to true, face_flip and face_rotation default to false (standard orientation) and has no effect in 2d. The concept of face orientations is explained in this glossary entry.

In the case of anisotropically refined cells and faces, the RefinementCase of the face, face_ref_case, might have an influence on which child is behind which given subface, thus this is an additional argument, defaulting to isotropic refinement of the face.

◆ line_to_cell_vertices()

template<int dim>
static unsigned int GeometryInfo< dim >::line_to_cell_vertices ( const unsigned int  line,
const unsigned int  vertex 
)
static

Map line vertex number to cell vertex number, i.e. give the cell vertex number of the vertexth vertex of line line, e.g. GeometryInfo<2>::line_to_cell_vertices(3,0)=2.

The order of the lines, as well as their direction (which in turn determines which is the first and which the second vertex on a line) is the canonical one in deal.II, as described in the general documentation of this class.

For dim=2 this call is simply passed down to the face_to_cell_vertices() function.

◆ face_to_cell_vertices()

template<int dim>
static unsigned int GeometryInfo< dim >::face_to_cell_vertices ( const unsigned int  face,
const unsigned int  vertex,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map face vertex number to cell vertex number, i.e. give the cell vertex number of the vertexth vertex of face face, e.g. GeometryInfo<2>::face_to_cell_vertices(3,0)=2, see the image under point N4 in the 2d section of this class's documentation.

Through the face_orientation, face_flip and face_rotation arguments this function handles faces oriented in the standard and non-standard orientation. face_orientation defaults to true, face_flip and face_rotation default to false (standard orientation). In 2d only face_flip is considered. See this glossary article for more information.

As the children of a cell are ordered according to the vertices of the cell, this call is passed down to the child_cell_on_face() function. Hence this function is simply a wrapper of child_cell_on_face() giving it a suggestive name.

◆ face_to_cell_lines()

template<int dim>
static unsigned int GeometryInfo< dim >::face_to_cell_lines ( const unsigned int  face,
const unsigned int  line,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map face line number to cell line number, i.e. give the cell line number of the lineth line of face face, e.g. GeometryInfo<3>::face_to_cell_lines(5,0)=4.

Through the face_orientation, face_flip and face_rotation arguments this function handles faces oriented in the standard and non-standard orientation. face_orientation defaults to true, face_flip and face_rotation default to false (standard orientation) and has no effect in 2d.

◆ standard_to_real_face_vertex()

template<int dim>
static unsigned int GeometryInfo< dim >::standard_to_real_face_vertex ( const unsigned int  vertex,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map the vertex index vertex of a face in standard orientation to one of a face with arbitrary face_orientation, face_flip and face_rotation. The values of these three flags default to true, false and false, respectively. this combination describes a face in standard orientation.

This function is only implemented in 3d.

◆ real_to_standard_face_vertex()

template<int dim>
static unsigned int GeometryInfo< dim >::real_to_standard_face_vertex ( const unsigned int  vertex,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map the vertex index vertex of a face with arbitrary face_orientation, face_flip and face_rotation to a face in standard orientation. The values of these three flags default to true, false and false, respectively. this combination describes a face in standard orientation.

This function is only implemented in 3d.

◆ standard_to_real_face_line()

template<int dim>
static unsigned int GeometryInfo< dim >::standard_to_real_face_line ( const unsigned int  line,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map the line index line of a face in standard orientation to one of a face with arbitrary face_orientation, face_flip and face_rotation. The values of these three flags default to true, false and false, respectively. this combination describes a face in standard orientation.

This function is only implemented in 3d.

◆ standard_to_real_line_vertex()

template<int dim>
static unsigned int GeometryInfo< dim >::standard_to_real_line_vertex ( const unsigned int  vertex,
const bool  line_orientation = true 
)
static

Map the vertex index vertex of a line in standard orientation to one of a face with arbitrary line_orientation. The value of this flag default to true.

◆ standard_quad_vertex_to_line_vertex_index()

template<int dim>
static std::array< unsigned int, 2 > GeometryInfo< dim >::standard_quad_vertex_to_line_vertex_index ( const unsigned int  vertex)
static

Decompose the vertex index in a quad into a pair of a line index and a vertex index within this line.

Note
Which line is selected is not of importance (and not exposed on purpose).

◆ standard_hex_vertex_to_quad_vertex_index()

template<int dim>
static std::array< unsigned int, 2 > GeometryInfo< dim >::standard_hex_vertex_to_quad_vertex_index ( const unsigned int  vertex)
static

Decompose the vertex index in a hex into a pair of a quad index and a vertex index within this quad.

Note
Which quad is selected is not of importance (and not exposed on purpose).

◆ standard_hex_line_to_quad_line_index()

template<int dim>
static std::array< unsigned int, 2 > GeometryInfo< dim >::standard_hex_line_to_quad_line_index ( const unsigned int  line)
static

Decompose the line index in a hex into a pair of a quad index and a line index within this quad.

Note
Which quad is selected is not of importance (and not exposed on purpose).

◆ real_to_standard_face_line()

template<int dim>
static unsigned int GeometryInfo< dim >::real_to_standard_face_line ( const unsigned int  line,
const bool  face_orientation = true,
const bool  face_flip = false,
const bool  face_rotation = false 
)
static

Map the line index line of a face with arbitrary face_orientation, face_flip and face_rotation to a face in standard orientation. The values of these three flags default to true, false and false, respectively. this combination describes a face in standard orientation.

This function is only implemented in 3d.

◆ unit_cell_vertex()

template<int dim>
static Point< dim > GeometryInfo< dim >::unit_cell_vertex ( const unsigned int  vertex)
static

Return the position of the ith vertex on the unit cell. The order of vertices is the canonical one in deal.II, as described in the general documentation of this class.

◆ child_cell_from_point()

template<int dim>
static unsigned int GeometryInfo< dim >::child_cell_from_point ( const Point< dim > &  p)
static

Given a point p in unit coordinates, return the number of the child cell in which it would lie in. If the point lies on the interface of two children, return any one of their indices. The result is always less than GeometryInfo<dimension>::max_children_per_cell.

The order of child cells is described the general documentation of this class.

◆ cell_to_child_coordinates()

template<int dim>
static Point< dim > GeometryInfo< dim >::cell_to_child_coordinates ( const Point< dim > &  p,
const unsigned int  child_index,
const RefinementCase< dim >  refine_case = RefinementCase< dim >::isotropic_refinement 
)
static

Given coordinates p on the unit cell, return the values of the coordinates of this point in the coordinate system of the given child. Neither original nor returned coordinates need actually be inside the cell, we simply perform a scale-and-shift operation with a shift that depends on the number of the child.

◆ child_to_cell_coordinates()

template<int dim>
static Point< dim > GeometryInfo< dim >::child_to_cell_coordinates ( const Point< dim > &  p,
const unsigned int  child_index,
const RefinementCase< dim >  refine_case = RefinementCase< dim >::isotropic_refinement 
)
static

The reverse function to the one above: take a point in the coordinate system of the child, and transform it to the coordinate system of the mother cell.

◆ is_inside_unit_cell() [1/2]

template<int dim>
static bool GeometryInfo< dim >::is_inside_unit_cell ( const Point< dim > &  p)
static

Return true if the given point is inside the unit cell of the present space dimension.

◆ is_inside_unit_cell() [2/2]

template<int dim>
static bool GeometryInfo< dim >::is_inside_unit_cell ( const Point< dim > &  p,
const double  eps 
)
static

Return true if the given point is inside the unit cell of the present space dimension. This function accepts an additional parameter which specifies how much the point position may actually be outside the true unit cell. This is useful because in practice we may often not be able to compute the coordinates of a point in reference coordinates exactly, but only up to numerical roundoff.

The tolerance parameter may be less than zero, indicating that the point should be safely inside the cell.

◆ project_to_unit_cell()

template<int dim>
template<typename Number = double>
static Point< dim, Number > GeometryInfo< dim >::project_to_unit_cell ( const Point< dim, Number > &  p)
static

Projects a given point onto the unit cell, i.e. each coordinate outside [0..1] is modified to lie within that interval.

◆ distance_to_unit_cell()

template<int dim>
static double GeometryInfo< dim >::distance_to_unit_cell ( const Point< dim > &  p)
static

Return the infinity norm of the vector between a given point p outside the unit cell to the closest unit cell boundary. For points inside the cell, this is defined as zero.

◆ d_linear_shape_function()

template<int dim>
static double GeometryInfo< dim >::d_linear_shape_function ( const Point< dim > &  xi,
const unsigned int  i 
)
static

Compute the value of the \(i\)-th \(d\)-linear (i.e. (bi-,tri-)linear) shape function at location \(\xi\).

◆ d_linear_shape_function_gradient()

template<int dim>
static Tensor< 1, dim > GeometryInfo< dim >::d_linear_shape_function_gradient ( const Point< dim > &  xi,
const unsigned int  i 
)
static

Compute the gradient of the \(i\)-th \(d\)-linear (i.e. (bi-,tri-)linear) shape function at location \(\xi\).

◆ alternating_form_at_vertices()

template<int dim>
template<int spacedim>
template void GeometryInfo< dim >::alternating_form_at_vertices ( const Point< spacedim >(&)  vertices[vertices_per_cell],
Tensor< spacedim - dim, spacedim >(&)  forms[vertices_per_cell] 
)
static

For a (bi-, tri-)linear mapping from the reference cell, face, or edge to the object specified by the given vertices, compute the alternating form of the transformed unit vectors vertices. For an object of dimensionality dim, there are dim vectors with spacedim components each, and the alternating form is a tensor of rank spacedim-dim that corresponds to the wedge product of the dim unit vectors, and it corresponds to the volume and normal vectors of the mapping from reference element to the element described by the vertices.

For example, if dim==spacedim==2, then the alternating form is a scalar (because spacedim-dim=0) and its value equals \(\mathbf v_1\wedge \mathbf v_2=\mathbf v_1^\perp \cdot\mathbf v_2\), where \(\mathbf v_1^\perp\) is a vector that is rotated to the right by 90 degrees from \(\mathbf v_1\). If dim==spacedim==3, then the result is again a scalar with value \(\mathbf v_1\wedge \mathbf v_2 \wedge \mathbf v_3 = (\mathbf v_1\times \mathbf v_2)\cdot \mathbf v_3\), where \(\mathbf v_1, \mathbf v_2, \mathbf v_3\) are the images of the unit vectors at a vertex of the unit dim-dimensional cell under transformation to the dim-dimensional cell in spacedim-dimensional space. In both cases, i.e. for dim==2 or 3, the result happens to equal the determinant of the Jacobian of the mapping from reference cell to cell in real space. Note that it is the actual determinant, not its absolute value as often used in transforming integrals from one coordinate system to another. In particular, if the object specified by the vertices is a parallelogram (i.e. a linear transformation of the reference cell) then the computed values are the same at all vertices and equal the (signed) area of the cell; similarly, for parallel-epipeds, it is the volume of the cell.

Likewise, if we have dim==spacedim-1 (e.g. we have a quad in 3d space, or a line in 2d), then the alternating product denotes the normal vector (i.e. a rank-1 tensor, since spacedim-dim=1) to the object at each vertex, where the normal vector's magnitude denotes the area element of the transformation from the reference object to the object given by the vertices. In particular, if again the mapping from reference object to the object under consideration here is linear (not bi- or trilinear), then the returned vectors are all parallel, perpendicular to the mapped object described by the vertices, and have a magnitude equal to the area/volume of the mapped object. If dim=1, spacedim=2, then the returned value is \(\mathbf v_1^\perp\), where \(\mathbf v_1\) is the image of the sole unit vector of a line mapped to the line in 2d given by the vertices; if dim=2, spacedim=3, then the returned values are \(\mathbf v_1 \wedge \mathbf v_2=\mathbf v_1 \times \mathbf v_2\) where \(\mathbf v_1,\mathbf v_2\) are the two three-dimensional vectors that are tangential to the quad mapped into three-dimensional space.

This function is used in order to determine how distorted a cell is (see the entry on distorted cells in the glossary).

Member Data Documentation

◆ max_children_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::max_children_per_cell = 1 << dim
staticconstexpr

Maximum number of children of a refined cell, i.e. the number of children of an isotropically refined cell.

If a cell is refined anisotropically, the actual number of children may be less than the value given here.

Definition at line 1994 of file geometry_info.h.

◆ faces_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::faces_per_cell = 2 * dim
staticconstexpr

Number of faces of a cell.

Definition at line 1999 of file geometry_info.h.

◆ max_children_per_face

template<int dim>
constexpr unsigned int GeometryInfo< dim >::max_children_per_face
staticconstexpr
Initial value:
=
static constexpr unsigned int max_children_per_cell

Maximum number of children of a refined face, i.e. the number of children of an isotropically refined face.

If a cell is refined anisotropically, the actual number of children may be less than the value given here.

Definition at line 2027 of file geometry_info.h.

◆ vertices_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::vertices_per_cell = 1 << dim
staticconstexpr

Number of vertices of a cell.

Definition at line 2033 of file geometry_info.h.

◆ vertices_per_face

template<int dim>
constexpr unsigned int GeometryInfo< dim >::vertices_per_face
staticconstexpr
Initial value:
=
static constexpr unsigned int vertices_per_cell

Number of vertices on each face.

Definition at line 2056 of file geometry_info.h.

◆ lines_per_face

template<int dim>
constexpr unsigned int GeometryInfo< dim >::lines_per_face
staticconstexpr
Initial value:
=
static constexpr unsigned int lines_per_cell

Number of lines on each face.

Definition at line 2062 of file geometry_info.h.

◆ quads_per_face

template<int dim>
constexpr unsigned int GeometryInfo< dim >::quads_per_face
staticconstexpr
Initial value:
=
static constexpr unsigned int quads_per_cell

Number of quads on each face.

Definition at line 2068 of file geometry_info.h.

◆ lines_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::lines_per_cell
staticconstexpr
Initial value:

Number of lines of a cell.

The formula to compute this makes use of the fact that when going from one dimension to the next, the object of the lower dimension is copied once (thus twice the old number of lines) and then a new line is inserted between each vertex of the old object and the corresponding one in the copy.

Definition at line 2080 of file geometry_info.h.

◆ quads_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::quads_per_cell
staticconstexpr
Initial value:

Number of quadrilaterals of a cell.

This number is computed recursively just as the previous one, with the exception that new quads result from connecting an original line and its copy.

Definition at line 2091 of file geometry_info.h.

◆ hexes_per_cell

template<int dim>
constexpr unsigned int GeometryInfo< dim >::hexes_per_cell
staticconstexpr
Initial value:
=
static constexpr unsigned int hexes_per_cell

Number of hexahedra of a cell.

Definition at line 2098 of file geometry_info.h.

◆ ucd_to_deal

template<int dim>
constexpr std::array< unsigned int, GeometryInfo< dim >::vertices_per_cell > GeometryInfo< dim >::ucd_to_deal
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::ucd_to_deal()

Rearrange vertices for UCD output. For a cell being written in UCD format, each entry in this field contains the number of a vertex in deal.II that corresponds to the UCD numbering at this location.

Typical example: write a cell and arrange the vertices, such that UCD understands them.

for (i=0; i< n_vertices; ++i)
out << cell->vertex(ucd_to_deal[i]);
static constexpr std::array< unsigned int, vertices_per_cell > ucd_to_deal

As the vertex numbering in deal.II versions <= 5.1 happened to coincide with the UCD numbering, this field can also be used like a old_to_lexicographic mapping.

Definition at line 2119 of file geometry_info.h.

◆ dx_to_deal

template<int dim>
constexpr std::array< unsigned int, GeometryInfo< dim >::vertices_per_cell > GeometryInfo< dim >::dx_to_deal
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::dx_to_deal()

Rearrange vertices for OpenDX output. For a cell being written in OpenDX format, each entry in this field contains the number of a vertex in deal.II that corresponds to the DX numbering at this location.

Typical example: write a cell and arrange the vertices, such that OpenDX understands them.

for (i=0; i< n_vertices; ++i)
out << cell->vertex(dx_to_deal[i]);
static constexpr std::array< unsigned int, vertices_per_cell > dx_to_deal

Definition at line 2135 of file geometry_info.h.

◆ vertex_to_face

template<int dim>
constexpr ndarray< unsigned int, GeometryInfo< dim >::vertices_per_cell, dim > GeometryInfo< dim >::vertex_to_face
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::vertex_to_face()

This field stores for each vertex to which faces it belongs. In any given dimension, the number of faces is equal to the dimension. The first index in this 2d-array runs over all vertices, the second index over dim faces to which the vertex belongs.

The order of the faces for each vertex is such that the first listed face bounds the reference cell in x direction, the second in y direction, and so on.

Definition at line 2149 of file geometry_info.h.

◆ unit_normal_direction

template<int dim>
constexpr std::array< unsigned int, GeometryInfo< dim >::faces_per_cell > GeometryInfo< dim >::unit_normal_direction
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::unit_normal_direction()

For each face of the reference cell, this field stores the coordinate direction in which its normal vector points. In dim dimension these are the 2*dim first entries of {0,0,1,1,2,2,3,3}.

Note that this is only the coordinate number. The actual direction of the normal vector is obtained by multiplying the unit vector in this direction with unit_normal_orientation.

Definition at line 2601 of file geometry_info.h.

◆ unit_normal_orientation

template<int dim>
constexpr std::array< int, GeometryInfo< dim >::faces_per_cell > GeometryInfo< dim >::unit_normal_orientation
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::unit_normal_orientation()

Orientation of the unit normal vector of a face of the reference cell. In dim dimension these are the 2*dim first entries of {-1,1,-1,1,-1,1,-1,1}.

Each value is either 1 or -1, corresponding to a normal vector pointing in the positive or negative coordinate direction, respectively.

Note that this is only the standard orientation of faces. At least in 3d, actual faces of cells in a triangulation can also have the opposite orientation, depending on a flag that one can query from the cell it belongs to. For more information, see the glossary entry on face orientation.

Definition at line 2620 of file geometry_info.h.

◆ unit_normal_vector

template<int dim>
constexpr std::array< Tensor< 1, dim >, GeometryInfo< dim >::faces_per_cell > GeometryInfo< dim >::unit_normal_vector
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::unit_normal_vector()

Unit normal vector (Point<dim>) of a face of the reference cell.

Note that this is only the standard orientation of faces. At least in 3d, actual faces of cells in a triangulation can also have the opposite orientation, depending on a flag that one can query from the cell it belongs to. For more information, see the glossary entry on face orientation.

Definition at line 2634 of file geometry_info.h.

◆ unit_tangential_vectors

template<int dim>
constexpr ndarray< Tensor< 1, dim >, GeometryInfo< dim >::faces_per_cell, dim - 1 > GeometryInfo< dim >::unit_tangential_vectors
staticconstexpr
Initial value:
= internal::GeometryInfoHelper::Initializers<
static constexpr ndarray< Tensor< 1, dim >, faces_per_cell, dim - 1 > unit_tangential_vectors

Unit tangential vectors (array of dim-1 elements of Point<dim>) of a face of the reference cell, arranged in a right-hand coordinate system such that the cross product between the two vectors return the unit normal vector.

Note that this is only the standard orientation of faces. At least in 3d, actual faces of cells in a triangulation can also have the opposite orientation, depending on a flag that one can query from the cell it belongs to. For more information, see the glossary entry on face orientation.

Definition at line 2651 of file geometry_info.h.

◆ opposite_face

template<int dim>
constexpr std::array< unsigned int, GeometryInfo< dim >::faces_per_cell > GeometryInfo< dim >::opposite_face
staticconstexpr
Initial value:
=
internal::GeometryInfoHelper::Initializers<dim>::opposite_face()

List of numbers which denotes which face is opposite to a given face. Its entries are the first 2*dim entries of { 1, 0, 3, 2, 5, 4, 7, 6}.

Definition at line 2659 of file geometry_info.h.


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