Reference documentation for deal.II version GIT relicensing-426-g7976cfd195 2024-04-18 21:10:01+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 | Public Attributes | List of all members
internal::MatrixFreeFunctions::UnivariateShapeData< Number > Struct Template Reference

#include <deal.II/matrix_free/shape_info.h>

Inheritance diagram for internal::MatrixFreeFunctions::UnivariateShapeData< Number >:
Inheritance graph
[legend]

Public Member Functions

 UnivariateShapeData ()
 
std::size_t memory_consumption () const
 
template<int dim, int spacedim>
void evaluate_shape_functions (const FiniteElement< dim, spacedim > &fe, const Quadrature< 1 > &quad, const std::vector< unsigned int > &lexicographic, const unsigned int direction)
 
template<int dim, int spacedim>
void evaluate_collocation_space (const FiniteElement< dim, spacedim > &fe, const Quadrature< 1 > &quad, const std::vector< unsigned int > &lexicographic, const unsigned int direction)
 
bool check_and_set_shapes_symmetric ()
 
bool check_shapes_collocation () const
 

Public Attributes

ElementType element_type
 
AlignedVector< Number > shape_values
 
AlignedVector< Number > shape_gradients
 
AlignedVector< Number > shape_hessians
 
AlignedVector< Number > shape_gradients_collocation
 
AlignedVector< Number > shape_hessians_collocation
 
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 > inverse_shape_values
 
AlignedVector< Number > inverse_shape_values_eo
 
std::array< AlignedVector< Number >, 2 > shape_data_on_face
 
std::array< AlignedVector< Number >, 2 > quadrature_data_on_face
 
std::array< AlignedVector< Number >, 2 > values_within_subface
 
std::array< AlignedVector< Number >, 2 > gradients_within_subface
 
std::array< AlignedVector< Number >, 2 > hessians_within_subface
 
std::array< AlignedVector< Number >, 2 > subface_interpolation_matrices
 
std::array< AlignedVector< typename ::internal::VectorizedArrayTrait< Number >::value_type >, 2 > subface_interpolation_matrices_scalar
 
Quadrature< 1 > quadrature
 
unsigned int fe_degree
 
unsigned int n_q_points_1d
 
bool nodal_at_cell_boundaries
 
Table< 3, Number > shape_values_face
 
Table< 4, Number > shape_gradients_face
 

Detailed Description

template<typename Number>
struct internal::MatrixFreeFunctions::UnivariateShapeData< Number >

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 133 of file shape_info.h.

Constructor & Destructor Documentation

◆ UnivariateShapeData()

template<typename Number >
internal::MatrixFreeFunctions::UnivariateShapeData< Number >::UnivariateShapeData ( )

Empty constructor. Sets default configuration.

Member Function Documentation

◆ memory_consumption()

template<typename Number >
std::size_t internal::MatrixFreeFunctions::UnivariateShapeData< Number >::memory_consumption ( ) const

Return the memory consumption of this class in bytes.

◆ evaluate_shape_functions()

template<typename Number >
template<int dim, int spacedim>
void internal::MatrixFreeFunctions::UnivariateShapeData< Number >::evaluate_shape_functions ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< 1 > &  quad,
const std::vector< unsigned int > &  lexicographic,
const unsigned int  direction 
)

Evaluate the finite element shape functions at the points of the given quadrature formula, filling the fields shape_[values,gradients,hessians] and related information.

The two last arguments 'lexicographic' and 'direction' are used to describe the unknowns along a single dimension, and the respective direction of derivatives.

◆ evaluate_collocation_space()

template<typename Number >
template<int dim, int spacedim>
void internal::MatrixFreeFunctions::UnivariateShapeData< Number >::evaluate_collocation_space ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< 1 > &  quad,
const std::vector< unsigned int > &  lexicographic,
const unsigned int  direction 
)

Evaluate the auxiliary polynomial space associated with the Lagrange polynomials in points of the given quadrature formula, filling the fields shape_[gradients,hessians]_collocation and related information.

◆ check_and_set_shapes_symmetric()

template<typename Number >
bool internal::MatrixFreeFunctions::UnivariateShapeData< Number >::check_and_set_shapes_symmetric ( )

Check whether we have symmetries in the shape values. In that case, also fill the shape_???_eo fields.

◆ check_shapes_collocation()

template<typename Number >
bool internal::MatrixFreeFunctions::UnivariateShapeData< Number >::check_shapes_collocation ( ) const

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.

Member Data Documentation

◆ element_type

template<typename Number >
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 196 of file shape_info.h.

◆ shape_values

template<typename Number >
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_values

Stores the shape values of the 1d finite element evaluated at 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 204 of file shape_info.h.

◆ shape_gradients

template<typename Number >
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_gradients

Stores the shape gradients of the 1d finite element evaluated at 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 212 of file shape_info.h.

◆ shape_hessians

template<typename Number >
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_hessians

Stores the shape Hessians of the 1d finite element evaluated at 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 220 of file shape_info.h.

◆ shape_gradients_collocation

template<typename Number >
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 226 of file shape_info.h.

◆ shape_hessians_collocation

template<typename Number >
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 232 of file shape_info.h.

◆ shape_values_eo

template<typename Number >
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 239 of file shape_info.h.

◆ shape_gradients_eo

template<typename Number >
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 246 of file shape_info.h.

◆ shape_hessians_eo

template<typename Number >
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 253 of file shape_info.h.

◆ shape_gradients_collocation_eo

template<typename Number >
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 261 of file shape_info.h.

◆ shape_hessians_collocation_eo

template<typename Number >
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 269 of file shape_info.h.

◆ inverse_shape_values

template<typename Number >
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 285 of file shape_info.h.

◆ inverse_shape_values_eo

template<typename Number >
AlignedVector<Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::inverse_shape_values_eo

Stores the even-odd variant of the inverse_shape_values field.

Definition at line 290 of file shape_info.h.

◆ shape_data_on_face

template<typename Number >
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. The sorting of data is to start with the values, then gradients, then second derivatives.

Definition at line 297 of file shape_info.h.

◆ quadrature_data_on_face

template<typename Number >
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. The sorting of data is to start with the values, then gradients, then second derivatives.

Definition at line 308 of file shape_info.h.

◆ values_within_subface

template<typename Number >
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 314 of file shape_info.h.

◆ gradients_within_subface

template<typename Number >
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 320 of file shape_info.h.

◆ hessians_within_subface

template<typename Number >
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 326 of file shape_info.h.

◆ subface_interpolation_matrices

template<typename Number >
std::array<AlignedVector<Number>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::subface_interpolation_matrices

A 1d subface interpolation matrices to the first and second quadrant. This data structure is only set up for FE_Q for dim > 1.

Definition at line 332 of file shape_info.h.

◆ subface_interpolation_matrices_scalar

template<typename Number >
std::array<AlignedVector<typename ::internal::VectorizedArrayTrait< Number>::value_type>, 2> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::subface_interpolation_matrices_scalar

Same as above but stored in a scalar format independent of the type of Number

Definition at line 341 of file shape_info.h.

◆ quadrature

template<typename Number >
Quadrature<1> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::quadrature

We store a copy of the one-dimensional quadrature formula used for initialization.

Definition at line 347 of file shape_info.h.

◆ fe_degree

template<typename Number >
unsigned int internal::MatrixFreeFunctions::UnivariateShapeData< Number >::fe_degree

Stores the degree of the element.

Definition at line 352 of file shape_info.h.

◆ n_q_points_1d

template<typename Number >
unsigned int internal::MatrixFreeFunctions::UnivariateShapeData< Number >::n_q_points_1d

Stores the number of quadrature points per dimension.

Definition at line 357 of file shape_info.h.

◆ nodal_at_cell_boundaries

template<typename Number >
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 363 of file shape_info.h.

◆ shape_values_face

template<typename Number >
Table<3, Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_values_face

Stores the shape values of the finite element evaluated at all quadrature points for all faces and orientations (no tensor-product structure exploited).

Definition at line 370 of file shape_info.h.

◆ shape_gradients_face

template<typename Number >
Table<4, Number> internal::MatrixFreeFunctions::UnivariateShapeData< Number >::shape_gradients_face

Stores the shape gradients of the finite element evaluated at all quadrature points for all faces, orientations, and directions (no tensor-product structure exploited).

Definition at line 377 of file shape_info.h.


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