Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30: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
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members

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

Public Member Functions

constexpr ReferenceCell ()
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 0 > &node_indices, const std::array< unsigned, 0 > &nodes_per_direction, const bool legacy_format) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 1 > &node_indices, const std::array< unsigned, 1 > &nodes_per_direction, const bool legacy_format) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 2 > &node_indices, const std::array< unsigned, 2 > &nodes_per_direction, const bool legacy_format) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 3 > &node_indices, const std::array< unsigned, 3 > &nodes_per_direction, const bool legacy_format) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 0 > &, const std::array< unsigned, 0 > &, const bool) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 1 > &, const std::array< unsigned, 1 > &, const bool) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 2 > &node_indices, const std::array< unsigned, 2 > &nodes_per_direction, const bool) const
 
template<>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, 3 > &node_indices, const std::array< unsigned, 3 > &nodes_per_direction, const bool legacy_format) const
 
Querying information about the kind of reference cells
bool is_hyper_cube () const
 
bool is_simplex () const
 
unsigned int get_dimension () const
 
Shape functions, mappings, quadratures defined on a reference cell
template<int dim>
double d_linear_shape_function (const Point< dim > &xi, const unsigned int i) const
 
template<int dim>
Tensor< 1, dim > d_linear_shape_function_gradient (const Point< dim > &xi, const unsigned int i) const
 
template<int dim, int spacedim = dim>
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping (const unsigned int degree) const
 
template<int dim, int spacedim = dim>
const Mapping< dim, spacedim > & get_default_linear_mapping () const
 
template<int dim>
Quadrature< dim > get_gauss_type_quadrature (const unsigned n_points_1d) const
 
template<int dim>
Quadrature< dim > get_midpoint_quadrature () const
 
template<int dim>
const Quadrature< dim > & get_nodal_type_quadrature () const
 
Querying the number of building blocks of a reference cell
unsigned int n_vertices () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intvertex_indices () const
 
template<int dim>
Point< dim > vertex (const unsigned int v) const
 
unsigned int n_lines () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intline_indices () const
 
unsigned int n_faces () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intface_indices () const
 
unsigned int n_isotropic_children () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intisotropic_child_indices () const
 
ReferenceCell face_reference_cell (const unsigned int face_no) const
 
double volume () const
 
template<int dim>
Point< dim > barycenter () const
 
template<int dim>
bool contains_point (const Point< dim > &p, const double tolerance=0) const
 
template<int dim>
Point< dim > closest_point (const Point< dim > &p) const
 
template<int dim>
Tensor< 1, dim > unit_tangential_vectors (const unsigned int face_no, const unsigned int i) const
 
template<int dim>
Tensor< 1, dim > unit_normal_vectors (const unsigned int face_no) const
 
unsigned int n_face_orientations (const unsigned int face_no) const
 
template<typename T , std::size_t N>
unsigned char compute_orientation (const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1) const
 
template<typename T >
unsigned char get_combined_orientation (const ArrayView< const T > &vertices_0, const ArrayView< const T > &vertices_1) const
 
template<typename T , std::size_t N>
std::array< T, N > permute_according_orientation (const std::array< T, N > &vertices, const unsigned int orientation) const
 
template<typename T >
boost::container::small_vector< T, 8 > permute_by_combined_orientation (const ArrayView< const T > &vertices, const unsigned char orientation) const
 
unsigned char get_inverse_combined_orientation (const unsigned char orientation) const
 
ArrayView< const unsigned intfaces_for_given_vertex (const unsigned int vertex_index) const
 
Translating between deal.II indexing and formats used by other programs
unsigned int exodusii_vertex_to_deal_vertex (const unsigned int vertex_n) const
 
unsigned int exodusii_face_to_deal_face (const unsigned int face_n) const
 
unsigned int unv_vertex_to_deal_vertex (const unsigned int vertex_n) const
 
unsigned int vtk_linear_type () const
 
unsigned int vtk_quadratic_type () const
 
unsigned int vtk_lagrange_type () const
 
template<int dim>
unsigned int vtk_lexicographic_to_node_index (const std::array< unsigned, dim > &node_indices, const std::array< unsigned, dim > &nodes_per_direction, const bool legacy_format) const
 
unsigned int vtk_vertex_to_deal_vertex (const unsigned int vertex_index) const
 
unsigned int gmsh_element_type () const
 

Static Public Member Functions

static ReferenceCell n_vertices_to_type (const int dim, const unsigned int n_vertices)
 

Private Member Functions

constexpr ReferenceCell (const std::uint8_t kind)
 

Private Attributes

std::uint8_t kind
 

Static Private Attributes

static constexpr ndarray< unsigned int, 2, 2 > line_vertex_permutations
 
static constexpr ndarray< unsigned int, 6, 3 > triangle_vertex_permutations
 
static constexpr ndarray< unsigned int, 8, 4 > quadrilateral_vertex_permutations
 

Friends

constexpr ReferenceCell internal::make_reference_cell_from_int (const std::uint8_t)
 
std::ostream & operator<< (std::ostream &out, const ReferenceCell &reference_cell)
 
std::istream & operator>> (std::istream &in, ReferenceCell &reference_cell)
 

Relationships between objects in the cell and on faces

static constexpr unsigned char default_combined_face_orientation ()
 
static constexpr unsigned char reversed_combined_line_orientation ()
 
unsigned int child_cell_on_face (const unsigned int face, const unsigned int subface, const unsigned char face_orientation=default_combined_face_orientation()) const
 
std::array< unsigned int, 2 > standard_vertex_to_face_and_vertex_index (const unsigned int vertex) const
 
std::array< unsigned int, 2 > standard_line_to_face_and_line_index (const unsigned int line) const
 
unsigned int line_to_cell_vertices (const unsigned int line, const unsigned int vertex) const
 
unsigned int face_to_cell_lines (const unsigned int face, const unsigned int line, const unsigned char face_orientation) const
 
unsigned int face_to_cell_vertices (const unsigned int face, const unsigned int vertex, const unsigned char face_orientation) const
 
template<int dim>
Point< dim > face_vertex_location (const unsigned int face, const unsigned int vertex) const
 
unsigned int standard_to_real_face_vertex (const unsigned int vertex, const unsigned int face, const unsigned char face_orientation) const
 
unsigned int standard_to_real_face_line (const unsigned int line, const unsigned int face, const unsigned char face_orientation) const
 
bool standard_vs_true_line_orientation (const unsigned int line, const unsigned int face, const unsigned char face_orientation, const bool line_orientation) const
 

Other functions

static constexpr std::size_t memory_consumption ()
 
std::string to_string () const
 
constexpr operator std::uint8_t () const
 
constexpr bool operator== (const ReferenceCell &type) const
 
constexpr bool operator!= (const ReferenceCell &type) const
 
template<class Archive >
void serialize (Archive &archive, const unsigned int)
 

Detailed Description

A type that describes the kinds of reference cells that can be used. This includes quadrilaterals and hexahedra (i.e., "hypercubes"), triangles and tetrahedra (simplices), and the pyramids and wedges necessary when using mixed 3d meshes. This class then describes geometric, topological, and other kinds of information about these kinds of reference cells. This includes how many vertices or faces a certain kind of reference cell has (topological information), where these vertices lie, what the cell's volume or center of mass is (geometric information), and how to output these cells in various output formats or what appropriate quadrature rules are. The documentation of this class is separated into a number of sections to group the many member functions into different categories such as those mentioned above.

Objects of this type should not be created in user code, and as a consequence the class does not have a user-accessible constructor other than the default constructor (which creates an invalid object). Rather, there is a finite number of specific reference cell objects defined in the ReferenceCells namespace that completely enumerate all of the possible values. User codes should therefore rely exclusively on assigning ReferenceCell objects from these special objects, and comparing against those special objects.

The purposes and intents of this class are described in the reference cell glossary entry.

Definition at line 98 of file reference_cell.h.

Constructor & Destructor Documentation

◆ ReferenceCell() [1/2]

constexpr ReferenceCell::ReferenceCell ( )
inlineconstexpr

Default constructor. Initialize this object as an invalid object. The end result is that the current object equals ReferenceCells::Invalid.

Generally, ReferenceCell objects are created by assignment from the special objects in namespace ReferenceCells, which is the only way to obtain a valid object.

Definition at line 1038 of file reference_cell.h.

◆ ReferenceCell() [2/2]

constexpr ReferenceCell::ReferenceCell ( const std::uint8_t  kind)
inlineconstexprprivate

Constructor. This is the constructor used to create the different static member variables of this class. It is private but can be called by a function in an internal namespace that is a friend of this class.

Definition at line 939 of file reference_cell.h.

Member Function Documentation

◆ n_vertices_to_type()

ReferenceCell ReferenceCell::n_vertices_to_type ( const int  dim,
const unsigned int  n_vertices 
)
inlinestatic

Return the correct ReferenceCell for a given structural dimension and number of vertices. For example, if dim==2 and n_vertices==4, this function will return ReferenceCells::Quadrilateral. But if dim==3 and n_vertices==4, it will return ReferenceCells::Tetrahedron.

Definition at line 2151 of file reference_cell.h.

◆ is_hyper_cube()

bool ReferenceCell::is_hyper_cube ( ) const
inline

◆ is_simplex()

bool ReferenceCell::is_simplex ( ) const
inline

Return true if the object is a Vertex, Line, Triangle, or Tetrahedron.

Definition at line 1130 of file reference_cell.h.

◆ get_dimension()

unsigned int ReferenceCell::get_dimension ( ) const
inline

Return the dimension of the reference cell represented by the current object.

Definition at line 1140 of file reference_cell.h.

◆ d_linear_shape_function()

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

Compute the value of the \(i\)-th linear shape function at location \(\xi\) for the current reference-cell type.

Definition at line 2193 of file reference_cell.h.

◆ d_linear_shape_function_gradient()

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

Compute the gradient of the \(i\)-th linear shape function at location \(\xi\) for the current reference-cell type.

Definition at line 2282 of file reference_cell.h.

◆ get_default_mapping()

template<int dim, int spacedim>
std::unique_ptr< Mapping< dim, spacedim > > ReferenceCell::get_default_mapping ( const unsigned int  degree) const

Return a default mapping of degree degree matching the current reference cell. If this reference cell is a hypercube, then the returned mapping is a MappingQ; otherwise, it is an object of type MappingFE initialized with FE_SimplexP (if the reference cell is a triangle or tetrahedron), with FE_PyramidP (if the reference cell is a pyramid), or with FE_WedgeP (if the reference cell is a wedge).

Definition at line 121 of file reference_cell.cc.

◆ get_default_linear_mapping()

template<int dim, int spacedim>
const Mapping< dim, spacedim > & ReferenceCell::get_default_linear_mapping ( ) const

Return a default linear mapping matching the current reference cell. If this reference cell is a hypercube, then the returned mapping is a MappingQ1; otherwise, it is an object of type MappingFE initialized with FE_SimplexP (if the reference cell is a triangle or tetrahedron), with FE_PyramidP (if the reference cell is a pyramid), or with FE_WedgeP (if the reference cell is a wedge). In other words, the term "linear" in the name of the function has to be understood as \(d\)-linear (i.e., bilinear or trilinear) for some of the coordinate directions.

Definition at line 148 of file reference_cell.cc.

◆ get_gauss_type_quadrature()

template<int dim>
template Quadrature< 0 > ReferenceCell::get_gauss_type_quadrature ( const unsigned  n_points_1d) const

Return a Gauss-type quadrature matching the given reference cell (QGauss, QGaussSimplex, QGaussPyramid, QGaussWedge).

Parameters
[in]n_points_1dThe number of quadrature points in each direction (QGauss) or an indication of what polynomial degree needs to be integrated exactly for the other types.

Definition at line 186 of file reference_cell.cc.

◆ get_midpoint_quadrature()

template<int dim>
Quadrature< dim > ReferenceCell::get_midpoint_quadrature ( ) const

Return a quadrature object that has a single quadrature point at the barycenter of the cell with quadrature weight equal to the volume of the reference cell. This quadrature formula is exact for integrals of constant and linear integrands.

The object returned by this function generalizes what the QMidpoint class represents to other reference cells.

Definition at line 1167 of file reference_cell.h.

◆ get_nodal_type_quadrature()

template<int dim>
template const Quadrature< 0 > & ReferenceCell::get_nodal_type_quadrature ( ) const

Return a quadrature rule whose quadrature points are the vertices of the given reference cell. For 1d line segments, this corresponds to the quadrature points of the trapezoidal rule, which by taking tensor products easily generalizes also to other hypercube elements (see also QTrapezoid). For all reference cell shapes, the quadrature points are ordered in the same order as the vertices of the reference cell.

Note
The weights of the quadrature object are left unfilled and consequently the object cannot usefully be used for actually computing integrals. This is in contrast to, for example, the QTrapezoid class that correctly sets quadrature weights.

Definition at line 208 of file reference_cell.cc.

◆ n_vertices()

unsigned int ReferenceCell::n_vertices ( ) const
inline

Return the number of vertices that make up the reference cell in question. A vertex is a "corner" (a zero-dimensional object) of the reference cell.

Definition at line 1176 of file reference_cell.h.

◆ vertex_indices()

std_cxx20::ranges::iota_view< unsigned int, unsigned int > ReferenceCell::vertex_indices ( ) const
inline

Return an object that can be thought of as an array containing all indices from zero to n_vertices().

Definition at line 1435 of file reference_cell.h.

◆ vertex()

template<int dim>
Point< dim > ReferenceCell::vertex ( const unsigned int  v) const

Return the location of the vth vertex of the reference cell that corresponds to the current object.

Because the ReferenceCell class does not have a dim argument, it has to be explicitly specified in the call to this function.

Definition at line 1237 of file reference_cell.h.

◆ n_lines()

unsigned int ReferenceCell::n_lines ( ) const
inline

Return the number of lines that make up the reference cell in question. A line is an "edge" (a one-dimensional object) of the reference cell.

Definition at line 1206 of file reference_cell.h.

◆ line_indices()

std_cxx20::ranges::iota_view< unsigned int, unsigned int > ReferenceCell::line_indices ( ) const
inline

Return an object that can be thought of as an array containing all indices from zero to n_lines().

Definition at line 1443 of file reference_cell.h.

◆ n_faces()

unsigned int ReferenceCell::n_faces ( ) const
inline

Return the number of faces that make up the reference cell in question. A face is a (dim-1)-dimensional object bounding the reference cell.

Definition at line 1356 of file reference_cell.h.

◆ face_indices()

std_cxx20::ranges::iota_view< unsigned int, unsigned int > ReferenceCell::face_indices ( ) const
inline

Return an object that can be thought of as an array containing all indices from zero to n_faces().

Definition at line 1386 of file reference_cell.h.

◆ n_isotropic_children()

unsigned int ReferenceCell::n_isotropic_children ( ) const
inline

Return the number of cells one would get by isotropically refining the current cell. Here, "isotropic refinement" means that we subdivide in each "direction" of a cell. For example, a square would be refined into four children by introducing new vertices along each edge and a new vertex in the cell center. For triangles, one would introduce new vertices at the center of each edge, and connect them to obtain four children. Similar constructions can be done for the other reference cell types.

Definition at line 1394 of file reference_cell.h.

◆ isotropic_child_indices()

std_cxx20::ranges::iota_view< unsigned int, unsigned int > ReferenceCell::isotropic_child_indices ( ) const
inline

Return an object that can be thought of as an array containing all indices from zero to n_isotropic_children().

Definition at line 1427 of file reference_cell.h.

◆ face_reference_cell()

ReferenceCell ReferenceCell::face_reference_cell ( const unsigned int  face_no) const
inline

Return the reference-cell type of face face_no of the current object. For example, if the current object is ReferenceCells::Tetrahedron, then face_no must be between in the interval \([0,4)\) and the function will always return ReferenceCells::Triangle. If the current object is ReferenceCells::Hexahedron, then face_no must be between in the interval \([0,6)\) and the function will always return ReferenceCells::Quadrilateral. For wedges and pyramids, the returned object may be either ReferenceCells::Triangle or ReferenceCells::Quadrilateral, depending on the given index.

Definition at line 1451 of file reference_cell.h.

◆ default_combined_face_orientation()

constexpr unsigned char ReferenceCell::default_combined_face_orientation ( )
inlinestaticconstexpr

Return the default combined face orientation flag (i.e., the default set of orientations, defined by orientation, rotate, and flip for a face in 3d).

Definition at line 1488 of file reference_cell.h.

◆ reversed_combined_line_orientation()

constexpr unsigned char ReferenceCell::reversed_combined_line_orientation ( )
inlinestaticconstexpr

Return the reversed (non-default orientation) line orientation flag. As lines only have two possible orientations, this function and ReferenceCell::default_combined_face_orientation() encode all of its possible orientation states.

Note
Line orientations are typically stored as booleans, but to better enable dimension-independent programming relevant functions typically present these values as unsigned chars.

Definition at line 1499 of file reference_cell.h.

◆ child_cell_on_face()

unsigned int ReferenceCell::child_cell_on_face ( const unsigned int  face,
const unsigned int  subface,
const unsigned char  face_orientation = default_combined_face_orientation() 
) const
inline

Return which child cells are adjacent to a certain face of the mother cell.

For example, in 2d the layout of a quadrilateral 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.

The face_orientation argument is meant exclusively for quadrilaterals and hexahedra at the moment. It determines how this function handles faces oriented in the standard and non-standard orientation. It represents a bit-code for the overall face_orientation, face_flip and face_rotation and defaults to the standard orientation. The concept of face orientations is explained in this glossary entry.

Definition at line 1509 of file reference_cell.h.

◆ standard_vertex_to_face_and_vertex_index()

std::array< unsigned int, 2 > ReferenceCell::standard_vertex_to_face_and_vertex_index ( const unsigned int  vertex) const
inline

For a given vertex in a cell, return a pair of a face index and a vertex index within this face.

Note
In practice, a vertex is of course generally part of more than one face, and one could return different faces and the corresponding index within. Which face this function chooses is often not of importance (and not exposed by this function on purpose).

Definition at line 1575 of file reference_cell.h.

◆ standard_line_to_face_and_line_index()

std::array< unsigned int, 2 > ReferenceCell::standard_line_to_face_and_line_index ( const unsigned int  line) const
inline

For a given line in a cell, return a pair of a face index and a line index within this face.

Note
In practice, a line is of course generally part of more than one face, and one could return different faces and the corresponding index within. Which face this function chooses is often not of importance (and not exposed by this function on purpose).

Definition at line 1637 of file reference_cell.h.

◆ line_to_cell_vertices()

unsigned int ReferenceCell::line_to_cell_vertices ( const unsigned int  line,
const unsigned int  vertex 
) const
inline

Map line vertex number to cell vertex number, i.e., return the cell vertex number of the vertexth vertex of line line.

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

Definition at line 1699 of file reference_cell.h.

◆ face_to_cell_lines()

unsigned int ReferenceCell::face_to_cell_lines ( const unsigned int  face,
const unsigned int  line,
const unsigned char  face_orientation 
) const
inline

Map face line number to cell line number.

Definition at line 1781 of file reference_cell.h.

◆ face_to_cell_vertices()

unsigned int ReferenceCell::face_to_cell_vertices ( const unsigned int  face,
const unsigned int  vertex,
const unsigned char  face_orientation 
) const
inline

Map face vertex number to cell vertex number.

Definition at line 1868 of file reference_cell.h.

◆ face_vertex_location()

template<int dim>
Point< dim > ReferenceCell::face_vertex_location ( const unsigned int  face,
const unsigned int  vertex 
) const

For a given face, in standard orientation, return the location of one of its vertices in the ambient space of the cell. For example, for a square or triangular 2d cell, the zeroth vertex of its zeroth face is located at \((0,0)\) – a location in 2d space.

Parameters
[in]faceThe number of face. This number must be between zero and n_faces().
[in]vertexThe number of the vertex within the face. This number must be between zero and face_reference_cell(face).n_vertices().
Returns
The location of the vertex so identified in dim-dimensional space.
Postcondition
The output of calling reference_cell.face_vertex_location<dim>(f,v) is identical to calling reference_cell.vertex<dim>( reference_cell.face_to_cell_vertices( f, v, ReferenceCell::default_combined_face_orientation())).

Definition at line 1970 of file reference_cell.h.

◆ standard_to_real_face_vertex()

unsigned int ReferenceCell::standard_to_real_face_vertex ( const unsigned int  vertex,
const unsigned int  face,
const unsigned char  face_orientation 
) const
inline

Correct vertex index depending on face orientation.

Definition at line 1980 of file reference_cell.h.

◆ standard_to_real_face_line()

unsigned int ReferenceCell::standard_to_real_face_line ( const unsigned int  line,
const unsigned int  face,
const unsigned char  face_orientation 
) const
inline

Correct line index depending on face orientation.

Definition at line 2028 of file reference_cell.h.

◆ standard_vs_true_line_orientation()

bool ReferenceCell::standard_vs_true_line_orientation ( const unsigned int  line,
const unsigned int  face,
const unsigned char  face_orientation,
const bool  line_orientation 
) const
inline

Return whether the line with index line is oriented in standard direction within a cell, given the face_orientation of the face within the current cell, and line_orientation flag for the line within that face. true indicates that the line is oriented from vertex 0 to vertex 1, whereas it is the other way around otherwise. In 1d and 2d, this is always true, but in 3d it may be different, see the respective discussion in the documentation of the GeometryInfo class.

Definition at line 2619 of file reference_cell.h.

◆ volume()

double ReferenceCell::volume ( ) const
inline

Return the \(d\)-dimensional volume of the reference cell that corresponds to the current object, where \(d\) is the dimension of the space it lives in. For example, since the quadrilateral reference cell is \([0,1]^2\), its volume is one, whereas the volume of the reference triangle is 0.5 because it occupies the area \(\{0 \le x,y \le 1, x+y\le 1\}\).

For ReferenceCells::Vertex, the reference cell is a zero-dimensional point in a zero-dimensional space. As a consequence, one cannot meaningfully define a volume for it. The function returns one for this case, because this makes it possible to define useful quadrature rules based on the center of a reference cell and its volume.

Definition at line 2316 of file reference_cell.h.

◆ barycenter()

template<int dim>
Point< dim > ReferenceCell::barycenter ( ) const
inline

Return the barycenter (i.e., the center of mass) of the reference cell that corresponds to the current object. The function is not called center() because one can define the center of an object in a number of different ways whereas the barycenter of a reference cell \(K\) is unambiguously defined as

\[ \mathbf x_K = \frac{1}{V} \int_K \mathbf x \; dx \]

where \(V\) is the volume of the reference cell (see also the volume() function).

Definition at line 2347 of file reference_cell.h.

◆ contains_point()

template<int dim>
bool ReferenceCell::contains_point ( const Point< dim > &  p,
const double  tolerance = 0 
) const
inline

Return true if the given point is inside the reference cell of the present space dimension up to some tolerance. This function accepts an additional parameter (which defaults to zero) which specifies by how much the point position may actually be outside the true reference 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. For example, strictly speaking one would expect that for points on the boundary of the reference cell, the function would return true if the tolerance was zero. But in practice, this may or may not actually be true; for example, the point \((1/3, 2/3)\) is on the boundary of the reference triangle because \(1/3+2/3 \le 1\), but since neither of its coordinates are exactly representable in floating point arithmetic, the floating point representations of \(1/3\) and \(2/3\) may or may not add up to anything that is less than or equal to one.

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

Definition at line 2380 of file reference_cell.h.

◆ closest_point()

template<int dim>
Point< dim > ReferenceCell::closest_point ( const Point< dim > &  p) const

Return the point on the surface of the reference cell closest (in the Euclidean norm) to p.

Definition at line 950 of file reference_cell.cc.

◆ unit_tangential_vectors()

template<int dim>
Tensor< 1, dim > ReferenceCell::unit_tangential_vectors ( const unsigned int  face_no,
const unsigned int  i 
) const
inline

Return \(i\)-th unit tangential vector of a face of the reference cell. The vectors are arranged such that the cross product between the two vectors returns the unit normal vector.

Precondition
\(i\) must be between zero and dim-1.

Definition at line 2489 of file reference_cell.h.

◆ unit_normal_vectors()

template<int dim>
Tensor< 1, dim > ReferenceCell::unit_normal_vectors ( const unsigned int  face_no) const
inline

Return the unit normal vector of a face of the reference cell.

Definition at line 2570 of file reference_cell.h.

◆ n_face_orientations()

unsigned int ReferenceCell::n_face_orientations ( const unsigned int  face_no) const
inline

Return the number of orientations for a face in the ReferenceCell. For example, for hexahedra this is 8 for every face since quadrilaterals have 8 possible orientations.

Definition at line 2600 of file reference_cell.h.

◆ compute_orientation()

template<typename T , std::size_t N>
unsigned char ReferenceCell::compute_orientation ( const std::array< T, N > &  vertices_0,
const std::array< T, N > &  vertices_1 
) const
inline

Determine the orientation of the current entity described by its vertices vertices_1 relative to an entity described by vertices_0. The two arrays given as arguments can be arrays of global vertex indices or local vertex indices, arrays of vertex locations, or arrays of any other objects identifying the vertices and the order in which they are encountered in a cell.

The size of the arrays, i.e., the template argument N, must be equal to or larger than the number of vertices of the current entity. If it is larger, only those elements of the input and output arrays are read from or written to that correspond to valid vertex indices.

Deprecated:
Use get_combined_orientation() instead.

Definition at line 2756 of file reference_cell.h.

◆ get_combined_orientation()

template<typename T >
unsigned char ReferenceCell::get_combined_orientation ( const ArrayView< const T > &  vertices_0,
const ArrayView< const T > &  vertices_1 
) const

Determine the relative orientation of the current entity described by its vertices vertices_1 relative to an entity described by vertices_0. Relative orientations are special cases of permutations since every vertex has to appear in the list of vertices of a reoriented cell as well; however, not every permutation can denote the same cell: For example, a square's vertices can be rotated by 90, 180, or 270 degrees, and the cell can be inverted (in essence looking at it from the other side), but one can't just exchange the order of two adjacent vertices because then the resulting cell is no longer a square but an object with two edges that cross each other.

The two arrays given as arguments can be arrays of global vertex indices or local vertex indices, arrays of vertex locations, or arrays of any other objects identifying the vertices and the order in which they are encountered in a cell.

The size of the input arrays must be equal to the number of vertices of the current entity.

Returns
A number that describes a relative orientation. How exactly this index is defined is not important, but it is consistent with the understanding the permute_by_combined_orientation() has of these orientation indices.

Definition at line 2777 of file reference_cell.h.

◆ permute_according_orientation()

template<typename T , std::size_t N>
std::array< T, N > ReferenceCell::permute_according_orientation ( const std::array< T, N > &  vertices,
const unsigned int  orientation 
) const
inline

Inverse function of compute_orientation(): Given a set of vertex-associated objects (such as vertex indices, locations, etc.) and a desired orientation permutation, return the permuted vertex information.

The size of the input and output arrays, i.e., the template argument N, must be equal to or larger than the number of vertices of the current entity. If it is larger, only those elements of the input and output arrays are read from or written to that correspond to valid vertex indices.

Deprecated:
Use permute_by_combined_orientation() instead.

Definition at line 2890 of file reference_cell.h.

◆ permute_by_combined_orientation()

template<typename T >
boost::container::small_vector< T, 8 > ReferenceCell::permute_by_combined_orientation ( const ArrayView< const T > &  vertices,
const unsigned char  orientation 
) const

This is the inverse function to get_combined_orientation(): Given a set of vertex-associated objects (such as vertex indices, locations, etc.) and a desired orientation permutation, return the permuted vertex information.

The size of the input array must be equal to the number of vertices of the current entity. The output is an array or permuted quantities of the same size. It is a vector that can store up to and including as many elements as cells can have vertices (namely eight, as in the case of hexahedra in 3d).

Definition at line 2917 of file reference_cell.h.

◆ get_inverse_combined_orientation()

unsigned char ReferenceCell::get_inverse_combined_orientation ( const unsigned char  orientation) const
inline

Return the inverse orientation. This is the value such that calling permute_by_combined_orientation() with o and then calling it again with get_inverse_combined_orientation(o) is the identity operation.

Definition at line 2991 of file reference_cell.h.

◆ faces_for_given_vertex()

ArrayView< const unsigned int > ReferenceCell::faces_for_given_vertex ( const unsigned int  vertex_index) const
inline

Return a vector of faces a given vertex_index belongs to.

Definition at line 1062 of file reference_cell.h.

◆ exodusii_vertex_to_deal_vertex()

unsigned int ReferenceCell::exodusii_vertex_to_deal_vertex ( const unsigned int  vertex_n) const

Map an ExodusII vertex number to a deal.II vertex number.

Definition at line 253 of file reference_cell.cc.

◆ exodusii_face_to_deal_face()

unsigned int ReferenceCell::exodusii_face_to_deal_face ( const unsigned int  face_n) const

Map an ExodusII face number to a deal.II face number.

Definition at line 296 of file reference_cell.cc.

◆ unv_vertex_to_deal_vertex()

unsigned int ReferenceCell::unv_vertex_to_deal_vertex ( const unsigned int  vertex_n) const

Map a UNV vertex number to a deal.II vertex number.

Definition at line 343 of file reference_cell.cc.

◆ vtk_linear_type()

unsigned int ReferenceCell::vtk_linear_type ( ) const

Return a VTK linear shape constant that corresponds to the reference cell.

Definition at line 379 of file reference_cell.cc.

◆ vtk_quadratic_type()

unsigned int ReferenceCell::vtk_quadratic_type ( ) const

Return a VTK quadratic shape constant that corresponds to the reference cell.

Definition at line 411 of file reference_cell.cc.

◆ vtk_lagrange_type()

unsigned int ReferenceCell::vtk_lagrange_type ( ) const

Return a VTK Lagrange shape constant that corresponds to the reference cell.

Definition at line 443 of file reference_cell.cc.

◆ vtk_lexicographic_to_node_index() [1/9]

template<int dim>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, dim > &  node_indices,
const std::array< unsigned, dim > &  nodes_per_direction,
const bool  legacy_format 
) const

Given a set of node indices of the form \((i)\) or \((i,j)\) or \((i,j,k)\) (depending on whether the reference cell is in 1d, 2d, or 3d), return the index the VTK format uses for this node for cells that are subdivided as many times in each of the coordinate directions as described by the second argument. For a uniformly subdivided cell, the second argument is an array whose elements will all be equal.

The last argument, legacy_format, indicates whether to use the old, VTK legacy format (when true) or the new, VTU format (when false).

◆ vtk_vertex_to_deal_vertex()

unsigned int ReferenceCell::vtk_vertex_to_deal_vertex ( const unsigned int  vertex_index) const

Map a VTK vertex number to a deal.II vertex number.

Definition at line 672 of file reference_cell.cc.

◆ gmsh_element_type()

unsigned int ReferenceCell::gmsh_element_type ( ) const

Return the GMSH element type code that corresponds to the reference cell.

Definition at line 727 of file reference_cell.cc.

◆ to_string()

std::string ReferenceCell::to_string ( ) const

Return a text representation of the reference cell represented by the current object.

Definition at line 88 of file reference_cell.cc.

◆ operator std::uint8_t()

constexpr ReferenceCell::operator std::uint8_t ( ) const
inlineconstexpr

Conversion operator to an integer.

Definition at line 945 of file reference_cell.h.

◆ operator==()

constexpr bool ReferenceCell::operator== ( const ReferenceCell type) const
inlineconstexpr

Operator for equality comparison.

Definition at line 953 of file reference_cell.h.

◆ operator!=()

constexpr bool ReferenceCell::operator!= ( const ReferenceCell type) const
inlineconstexpr

Operator for inequality comparison.

Definition at line 961 of file reference_cell.h.

◆ serialize()

template<class Archive >
void ReferenceCell::serialize ( Archive &  archive,
const unsigned int   
)
inline

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 1046 of file reference_cell.h.

◆ memory_consumption()

constexpr std::size_t ReferenceCell::memory_consumption ( )
inlinestaticconstexpr

Return the number of bytes used by an instance of this class.

Definition at line 1054 of file reference_cell.h.

◆ vtk_lexicographic_to_node_index() [2/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 0 > &  node_indices,
const std::array< unsigned, 0 > &  nodes_per_direction,
const bool  legacy_format 
) const

◆ vtk_lexicographic_to_node_index() [3/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 1 > &  node_indices,
const std::array< unsigned, 1 > &  nodes_per_direction,
const bool  legacy_format 
) const

◆ vtk_lexicographic_to_node_index() [4/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 2 > &  node_indices,
const std::array< unsigned, 2 > &  nodes_per_direction,
const bool  legacy_format 
) const

◆ vtk_lexicographic_to_node_index() [5/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 3 > &  node_indices,
const std::array< unsigned, 3 > &  nodes_per_direction,
const bool  legacy_format 
) const

◆ vtk_lexicographic_to_node_index() [6/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 0 > &  ,
const std::array< unsigned, 0 > &  ,
const bool   
) const

Definition at line 475 of file reference_cell.cc.

◆ vtk_lexicographic_to_node_index() [7/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 1 > &  ,
const std::array< unsigned, 1 > &  ,
const bool   
) const

Definition at line 488 of file reference_cell.cc.

◆ vtk_lexicographic_to_node_index() [8/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 2 > &  node_indices,
const std::array< unsigned, 2 > &  nodes_per_direction,
const bool   
) const

◆ vtk_lexicographic_to_node_index() [9/9]

template<>
unsigned int ReferenceCell::vtk_lexicographic_to_node_index ( const std::array< unsigned, 3 > &  node_indices,
const std::array< unsigned, 3 > &  nodes_per_direction,
const bool  legacy_format 
) const

Friends And Related Symbol Documentation

◆ internal::make_reference_cell_from_int

constexpr ReferenceCell internal::make_reference_cell_from_int ( const std::uint8_t  )
friend

A kind of constructor – not quite private because it can be called by anyone, but at least hidden in an internal namespace.

◆ operator<<

std::ostream & operator<< ( std::ostream &  out,
const ReferenceCell reference_cell 
)
friend

Output operator that writes the reference_cell object to the stream in a text format in which the object is represented by an integer. The details of which integer value represents each kind of reference cell is unimportant and consequently not specified. If you want a string representation of what a ReferenceCell is, use ReferenceCell::to_string().

Definition at line 1131 of file reference_cell.cc.

◆ operator>>

std::istream & operator>> ( std::istream &  in,
ReferenceCell reference_cell 
)
friend

Input operator that reads the reference_cell object from the stream in a text format in which the object is represented by an integer. Which specific integer value represents which reference cell is unspecified, but the function uses the same translation as the corresponding output operator<<.

Definition at line 1145 of file reference_cell.cc.

Member Data Documentation

◆ kind

std::uint8_t ReferenceCell::kind
private

The variable that stores what this object actually corresponds to.

Definition at line 862 of file reference_cell.h.

◆ line_vertex_permutations

constexpr ndarray< unsigned int, 2, 2 > ReferenceCell::line_vertex_permutations
staticconstexprprivate
Initial value:
= {
{{{1, 0}}, {{0, 1}}}}

Table containing all vertex permutations for a line.

Definition at line 875 of file reference_cell.h.

◆ triangle_vertex_permutations

constexpr ndarray< unsigned int, 6, 3 > ReferenceCell::triangle_vertex_permutations
staticconstexprprivate
Initial value:
= {
{{{0, 2, 1}},
{{0, 1, 2}},
{{2, 1, 0}},
{{2, 0, 1}},
{{1, 0, 2}},
{{1, 2, 0}}}}

Table containing all vertex permutations for a triangle.

Definition at line 882 of file reference_cell.h.

◆ quadrilateral_vertex_permutations

constexpr ndarray< unsigned int, 8, 4 > ReferenceCell::quadrilateral_vertex_permutations
staticconstexprprivate
Initial value:
= {{{{0, 2, 1, 3}},
{{0, 1, 2, 3}},
{{2, 3, 0, 1}},
{{2, 0, 3, 1}},
{{3, 1, 2, 0}},
{{3, 2, 1, 0}},
{{1, 0, 3, 2}},
{{1, 3, 0, 2}}}}

Table containing all vertex permutations for a quadrilateral.

Definition at line 894 of file reference_cell.h.


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