Reference documentation for deal.II version 9.2.0
|
#include <deal.II/matrix_free/shape_info.h>
Public Member Functions | |
ShapeInfo () | |
template<int dim> | |
ShapeInfo (const Quadrature< 1 > &quad, const FiniteElement< dim > &fe, const unsigned int base_element=0) | |
template<int dim> | |
void | reinit (const Quadrature< 1 > &quad, const FiniteElement< dim > &fe_dim, const unsigned int base_element=0) |
const UnivariateShapeData< Number > & | get_shape_data (const unsigned int dimension=0, const unsigned int component=0) const |
std::size_t | memory_consumption () const |
Public Attributes | |
ElementType | element_type |
std::vector< unsigned int > | lexicographic_numbering |
std::vector< UnivariateShapeData< Number > > | data |
::Table< 2, UnivariateShapeData< Number > * > | data_access |
unsigned int | n_dimensions |
unsigned int | n_components |
unsigned int | n_q_points |
unsigned int | dofs_per_component_on_cell |
unsigned int | n_q_points_face |
unsigned int | dofs_per_component_on_face |
::Table< 2, unsigned int > | face_to_cell_index_nodal |
::Table< 2, unsigned int > | face_to_cell_index_hermite |
Private Member Functions | |
bool | check_1d_shapes_symmetric (UnivariateShapeData< Number > &univariate_shape_data) |
bool | check_1d_shapes_collocation (const UnivariateShapeData< Number > &univariate_shape_data) const |
This struct stores a tensor (Kronecker) product view of the finite element and quadrature formula used for evaluation. It is based on a single or a collection of UnivariateShapeData object(s) that describe one-dimensional ingredients, plus some additional information about how these are combined and how indices are laid out in the multi-dimensional case such as the hierarchical -> lexicographic ordering of FE_Q.
Definition at line 290 of file shape_info.h.
internal::MatrixFreeFunctions::ShapeInfo< Number >::ShapeInfo | ( | ) |
Empty constructor. Does nothing.
|
inline |
Constructor that initializes the data fields using the reinit method.
Definition at line 500 of file shape_info.h.
void internal::MatrixFreeFunctions::ShapeInfo< Number >::reinit | ( | const Quadrature< 1 > & | quad, |
const FiniteElement< dim > & | fe_dim, | ||
const unsigned int | base_element = 0 |
||
) |
Initializes the data fields. Takes a one-dimensional quadrature formula and a finite element as arguments and evaluates the shape functions, gradients and Hessians on the one-dimensional unit cell. This function assumes that the finite element is derived from a one- dimensional element by a tensor product and that the zeroth shape function in zero evaluates to one.
|
inline |
Return data of univariate shape functions which defines the dimension dimension
of tensor product shape functions regarding vector component component
of the underlying finite element.
Definition at line 516 of file shape_info.h.
std::size_t internal::MatrixFreeFunctions::ShapeInfo< Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
|
private |
Check whether we have symmetries in the shape values. In that case, also fill the shape_???_eo fields.
|
private |
Check whether symmetric 1D basis functions are such that the shape values form a diagonal matrix, i.e., the nodal points are collocated with the quadrature points. This allows for specialized algorithms that save some operations in the evaluation.
ElementType internal::MatrixFreeFunctions::ShapeInfo< Number >::element_type |
Encodes the type of element detected at construction. FEEvaluation will select the most efficient algorithm based on the given element type.
Definition at line 297 of file shape_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::ShapeInfo< Number >::lexicographic_numbering |
Renumbering from deal.II's numbering of cell degrees of freedom to lexicographic numbering used inside the FEEvaluation schemes of the underlying element in the DoFHandler. For vector-valued elements, the renumbering starts with a lexicographic numbering of the first component, then everything of the second component, and so on.
Definition at line 349 of file shape_info.h.
std::vector<UnivariateShapeData<Number> > internal::MatrixFreeFunctions::ShapeInfo< Number >::data |
Stores data of univariate shape functions defining the underlying tensor product finite element.
Definition at line 355 of file shape_info.h.
::Table<2, UnivariateShapeData<Number> *> internal::MatrixFreeFunctions::ShapeInfo< Number >::data_access |
Grants access to univariate shape function data of given dimension and vector component. Rows identify dimensions and columns identify vector components.
Definition at line 362 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_dimensions |
Stores the number of space dimensions.
Definition at line 367 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_components |
Stores the number of vector components of the underlying vector-valued finite element.
Definition at line 373 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_q_points |
Stores the number of quadrature points in dim
dimensions for a cell.
Definition at line 379 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::dofs_per_component_on_cell |
Stores the number of DoFs per cell of the scalar element in dim
dimensions.
Definition at line 385 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_q_points_face |
Stores the number of quadrature points per face in dim
dimensions.
Definition at line 390 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::dofs_per_component_on_face |
Stores the number of DoFs per face in dim
dimensions.
Definition at line 395 of file shape_info.h.
::Table<2, unsigned int> internal::MatrixFreeFunctions::ShapeInfo< Number >::face_to_cell_index_nodal |
For nodal basis functions with nodes located at the boundary of the unit cell, face integrals that involve only the values of the shape functions (approximations of first derivatives in DG) do not need to load all degrees of freedom of the cell but rather only the degrees of freedom located on the face. While it would also be possible to compute these indices on the fly, we choose to simplify the code and store the indirect addressing in this class.
The first table index runs through the faces of a cell, and the second runs through the nodal degrees of freedom of the face, using dofs_per_face
entries.
The indices stored in this member variable are as follows. Consider for example a 2D element of degree 3 with the following degrees of freedom in lexicographic numbering:
The first row stores the indices on the face with index 0, i.e., the numbers 0, 4, 8, 12
, the second row holds the indices 3, 7, 11, 15
for face 1, the third row holds the indices 0, 1, 2, 3
for face 2, and the last (fourth) row holds the indices 12, 13, 14, 15
. Similarly, the indices are stored in 3D. (Note that the y faces in 3D use indices reversed in terms of the lexicographic numbers due to the orientation of the coordinate system.)
nodal_at_cell_boundaries
evaluates to true
. Definition at line 432 of file shape_info.h.
::Table<2, unsigned int> internal::MatrixFreeFunctions::ShapeInfo< Number >::face_to_cell_index_hermite |
The face_to_cell_index_nodal
provides a shortcut for the evaluation of values on the faces. For Hermite-type basis functions, one can go one step further and also use shortcuts to get derivatives more cheaply where only two layers of degrees of freedom contribute to the derivative on the face. In the lexicographic ordering, the respective indices is in the next "layer" of degrees of freedom as compared to the nodal values. This array stores the indirect addressing of both the values and the gradient.
The first table index runs through the faces of a cell, and the second runs through the pairs of the nodal degrees of freedom of the face and the derivatives, using 2*dofs_per_face
entries.
The indices stored in this member variable are as follows. Consider for example a 2D element of degree 3 with the following degrees of freedom in lexicographic numbering:
The first row stores the indices for values and gradients on the face with index 0, i.e., the numbers 0, 1, 5, 6, 10, 11, 15, 16, 20, 21
, the second row holds the indices 4, 3, 9, 8, 14, 13, 19, 18, 24, 23
for face 1, the third row holds the indices 0, 5, 1, 6, 2, 7, 3, 8, 4, 9
for face 2, and the last (fourth) row holds the indices 20, 15, 21, 16, 22, 17, 23, 18, 24, 19
. Similarly, the indices are stored in 3D. (Note that the y faces in 3D use indices reversed in terms of the lexicographic numbers due to the orientation of the coordinate system.)
element_type
evaluates to tensor_symmetric_hermite
. Definition at line 472 of file shape_info.h.