Reference documentation for deal.II version 9.1.1
|
#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) |
std::size_t | memory_consumption () const |
bool | check_1d_shapes_symmetric (const unsigned int n_q_points_1d) |
bool | check_1d_shapes_collocation () |
Public Attributes | |
ElementType | element_type |
AlignedVector< Number > | shape_values |
AlignedVector< Number > | shape_gradients |
AlignedVector< Number > | shape_hessians |
AlignedVector< Number > | shape_values_eo |
AlignedVector< Number > | shape_gradients_eo |
AlignedVector< Number > | shape_hessians_eo |
AlignedVector< Number > | shape_gradients_collocation_eo |
AlignedVector< Number > | shape_hessians_collocation_eo |
AlignedVector< Number > | shape_data_on_face [2] |
AlignedVector< Number > | values_within_subface [2] |
AlignedVector< Number > | gradients_within_subface [2] |
AlignedVector< Number > | hessians_within_subface [2] |
std::vector< unsigned int > | lexicographic_numbering |
unsigned int | fe_degree |
unsigned int | n_q_points_1d |
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 |
bool | nodal_at_cell_boundaries |
::Table< 2, unsigned int > | face_to_cell_index_nodal |
::Table< 2, unsigned int > | face_to_cell_index_hermite |
The class that stores the shape functions, gradients and Hessians evaluated for a tensor product finite element and tensor product quadrature formula on the unit cell. Because of this structure, only one-dimensional data is stored.
Definition at line 106 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 382 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.
std::size_t internal::MatrixFreeFunctions::ShapeInfo< Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
bool internal::MatrixFreeFunctions::ShapeInfo< Number >::check_1d_shapes_symmetric | ( | const unsigned int | n_q_points_1d | ) |
Check whether we have symmetries in the shape values. In that case, also fill the shape_???_eo fields.
bool internal::MatrixFreeFunctions::ShapeInfo< Number >::check_1d_shapes_collocation | ( | ) |
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 146 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_values |
Stores the shape values of the 1D finite element evaluated on all 1D quadrature points in vectorized format, i.e., as an array of VectorizedArray<dim>::n_array_elements equal elements. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 155 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_gradients |
Stores the shape gradients of the 1D finite element evaluated on all 1D quadrature points in vectorized format, i.e., as an array of VectorizedArray<dim>::n_array_elements equal elements. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 164 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_hessians |
Stores the shape Hessians of the 1D finite element evaluated on all 1D quadrature points in vectorized format, i.e., as an array of VectorizedArray<dim>::n_array_elements equal elements. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 173 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_values_eo |
Stores the shape values in a different format, namely the so-called even-odd scheme where the symmetries in shape_values are used for faster evaluation.
Definition at line 180 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_gradients_eo |
Stores the shape gradients in a different format, namely the so- called even-odd scheme where the symmetries in shape_gradients are used for faster evaluation.
Definition at line 187 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_hessians_eo |
Stores the shape second derivatives in a different format, namely the so-called even-odd scheme where the symmetries in shape_hessians are used for faster evaluation.
Definition at line 194 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_gradients_collocation_eo |
Stores the shape gradients of the shape function space associated to the quadrature (collocation), given by FE_DGQ<1>(Quadrature<1>). For faster evaluation only the even-odd format is necessary.
Definition at line 201 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_hessians_collocation_eo |
Stores the shape hessians of the shape function space associated to the quadrature (collocation), given by FE_DGQ<1>(Quadrature<1>). For faster evaluation only the even-odd format is necessary.
Definition at line 208 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_data_on_face[2] |
Collects all data of 1D shape values evaluated at the point 0 and 1 (the vertices) in one data structure. Sorting is first the values, then gradients, then second derivatives.
Definition at line 215 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::values_within_subface[2] |
Stores one-dimensional values of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 221 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::gradients_within_subface[2] |
Stores one-dimensional gradients of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 227 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::hessians_within_subface[2] |
Stores one-dimensional gradients of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 233 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 242 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::fe_degree |
Stores the degree of the element.
Definition at line 247 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_q_points_1d |
Stores the number of quadrature points per dimension.
Definition at line 252 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 258 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 264 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 269 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 274 of file shape_info.h.
bool internal::MatrixFreeFunctions::ShapeInfo< Number >::nodal_at_cell_boundaries |
Indicates whether the basis functions are nodal in 0 and 1, i.e., the end points of the unit cell.
Definition at line 280 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 317 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 357 of file shape_info.h.