Reference documentation for deal.II version 9.3.3
|
#include <deal.II/matrix_free/shape_info.h>
Public Member Functions | |
UnivariateShapeData () | |
std::size_t | memory_consumption () const |
This struct stores the shape functions, their gradients and Hessians evaluated for a one-dimensional section of a tensor product finite element and tensor product quadrature formula in reference coordinates. This data structure also includes the evaluation of quantities at the cell boundary and on the sub-interval \((0, 0.5)\) and \((0.5, 1)\) for face integrals.
Definition at line 113 of file shape_info.h.
internal::MatrixFreeFunctions::UnivariateShapeData< Number >::UnivariateShapeData | ( | ) |
Empty constructor. Sets default configuration.
std::size_t internal::MatrixFreeFunctions::UnivariateShapeData< Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
ElementType internal::MatrixFreeFunctions::UnivariateShapeData< 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 131 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_values |
Stores the shape values of the 1D finite element evaluated on all 1D quadrature points. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 139 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_gradients |
Stores the shape gradients of the 1D finite element evaluated on all 1D quadrature points. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 147 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_hessians |
Stores the shape Hessians of the 1D finite element evaluated on all 1D quadrature points. 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::UnivariateShapeData< Number >::shape_gradients_collocation |
Stores the shape gradients of the shape function space associated to the quadrature (collocation), given by FE_DGQ<1>(Quadrature<1>).
Definition at line 161 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_hessians_collocation |
Stores the shape hessians of the shape function space associated to the quadrature (collocation), given by FE_DGQ<1>(Quadrature<1>).
Definition at line 167 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< 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 174 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< 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 181 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< 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 188 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< 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>). This array provides an alternative representation of the shape_gradients_collocation field in the even-odd format.
Definition at line 196 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< 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>). This array provides an alternative representation of the shape_hessians_collocation field in the even-odd format.
Definition at line 204 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::inverse_shape_values |
Stores the inverse transformation from the data at quadrature points to the basis defined by the shape_values fields. The data at quadrature points is interpreted either implicitly by its polynomial interpolation, or explicitly in terms of separate polynomials such as with the _collocation
fields. The size of the array equals the layout of the shape_values
array, and it is combined with the shape values array such that this matrix is the pseudo inverse of shape_values. In case the number of 1D quadrature points equals the size of the basis, this array is exactly the inverse of the shape_values array. The length of this array is n_dofs_1d * n_q_points_1d
and quadrature points are the index running fastest.
Definition at line 220 of file shape_info.h.
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::inverse_shape_values_eo |
Stores the even-odd variant of the inverse_shape_values
field.
Definition at line 225 of file shape_info.h.
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_data_on_face |
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 232 of file shape_info.h.
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::quadrature_data_on_face |
Collects all data of 1D nodal shape values (defined by the Lagrange polynomials in the points of the quadrature rule) evaluated at the point 0 and 1 (the vertices) in one data structure.
This data structure can be used to interpolate from the cell to the face quadrature points.
Definition at line 244 of file shape_info.h.
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::values_within_subface |
Stores one-dimensional values of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 250 of file shape_info.h.
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::gradients_within_subface |
Stores one-dimensional gradients of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 256 of file shape_info.h.
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::hessians_within_subface |
Stores one-dimensional gradients of shape functions on subface. Since there are two subfaces, store two variants.
Definition at line 262 of file shape_info.h.
Quadrature<1> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::quadrature |
We store a copy of the one-dimensional quadrature formula used for initialization.
Definition at line 268 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::UnivariateShapeData< Number >::fe_degree |
Stores the degree of the element.
Definition at line 273 of file shape_info.h.
unsigned int internal::MatrixFreeFunctions::UnivariateShapeData< Number >::n_q_points_1d |
Stores the number of quadrature points per dimension.
Definition at line 278 of file shape_info.h.
bool internal::MatrixFreeFunctions::UnivariateShapeData< 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 284 of file shape_info.h.
Table<3, Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_values_face |
Stores the shape values of the finite element evaluated on all quadrature points for all faces and orientations (no tensor-product structure exploited).
Definition at line 291 of file shape_info.h.
Table<4, Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_gradients_face |
Stores the shape gradients of the finite element evaluated on all quadrature points for all faces, orientations, and directions (no tensor-product structure exploited).
Definition at line 298 of file shape_info.h.