Reference documentation for deal.II version 8.5.1
Public Member Functions | Public Attributes | List of all members
internal::MatrixFreeFunctions::ShapeInfo< Number > Struct Template Reference

#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_gausslobatto ()
 

Public Attributes

ElementType element_type
 
AlignedVector< VectorizedArray< Number > > shape_values
 
AlignedVector< VectorizedArray< Number > > shape_gradients
 
AlignedVector< VectorizedArray< Number > > shape_hessians
 
AlignedVector< VectorizedArray< Number > > shape_val_evenodd
 
AlignedVector< VectorizedArray< Number > > shape_gra_evenodd
 
AlignedVector< VectorizedArray< Number > > shape_hes_evenodd
 
::Table< 2, unsigned int > face_indices
 
std::vector< Number > face_value [2]
 
std::vector< Number > face_gradient [2]
 
std::vector< Number > subface_value [2]
 
std::vector< Number > shape_values_number
 
std::vector< Number > shape_gradient_number
 
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_cell
 
unsigned int n_q_points_face
 
unsigned int dofs_per_face
 

Detailed Description

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

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.

Author
Katharina Kormann and Martin Kronbichler, 2010, 2011

Definition at line 59 of file shape_info.h.

Constructor & Destructor Documentation

◆ ShapeInfo() [1/2]

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

Empty constructor. Does nothing.

◆ ShapeInfo() [2/2]

template<typename Number >
template<int dim>
internal::MatrixFreeFunctions::ShapeInfo< Number >::ShapeInfo ( const Quadrature< 1 > &  quad,
const FiniteElement< dim > &  fe,
const unsigned int  base_element = 0 
)
inline

Constructor that initializes the data fields using the reinit method.

Definition at line 244 of file shape_info.h.

Member Function Documentation

◆ reinit()

template<typename Number>
template<int dim>
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.

◆ memory_consumption()

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

Return the memory consumption of this class in bytes.

◆ check_1d_shapes_symmetric()

template<typename Number>
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_???_evenodd fields.

◆ check_1d_shapes_gausslobatto()

template<typename Number>
bool internal::MatrixFreeFunctions::ShapeInfo< Number >::check_1d_shapes_gausslobatto ( )

Check whether symmetric 1D basis functions are such that the shape values form a diagonal matrix, which allows to use specialized algorithms that save some operations.

Member Data Documentation

◆ element_type

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

◆ shape_values

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

◆ shape_gradients

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

◆ shape_hessians

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

◆ shape_val_evenodd

template<typename Number>
AlignedVector<VectorizedArray<Number> > internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_val_evenodd

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

◆ shape_gra_evenodd

template<typename Number>
AlignedVector<VectorizedArray<Number> > internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_gra_evenodd

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

◆ shape_hes_evenodd

template<typename Number>
AlignedVector<VectorizedArray<Number> > internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_hes_evenodd

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

◆ face_indices

template<typename Number>
::Table<2,unsigned int> internal::MatrixFreeFunctions::ShapeInfo< Number >::face_indices

Stores the indices from cell DoFs to face DoFs. The rows go through the 2*dim faces, and the columns the DoFs on the faces.

Definition at line 151 of file shape_info.h.

◆ face_value

template<typename Number>
std::vector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::face_value[2]

Stores one-dimensional values of shape functions evaluated in zero and one, i.e., on the one-dimensional faces. Not vectorized.

Definition at line 157 of file shape_info.h.

◆ face_gradient

template<typename Number>
std::vector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::face_gradient[2]

Stores one-dimensional gradients of shape functions evaluated in zero and one, i.e., on the one-dimensional faces. Not vectorized.

Definition at line 163 of file shape_info.h.

◆ subface_value

template<typename Number>
std::vector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::subface_value[2]

Stores one-dimensional values of shape functions on subface. Since there are two subfaces, store two variants. Not vectorized.

Definition at line 169 of file shape_info.h.

◆ shape_values_number

template<typename Number>
std::vector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_values_number

Non-vectorized version of shape values. Needed when evaluating face info.

Definition at line 175 of file shape_info.h.

◆ shape_gradient_number

template<typename Number>
std::vector<Number> internal::MatrixFreeFunctions::ShapeInfo< Number >::shape_gradient_number

Non-vectorized version of shape gradients. Needed when evaluating face info.

Definition at line 181 of file shape_info.h.

◆ lexicographic_numbering

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

◆ fe_degree

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

Stores the degree of the element.

Definition at line 195 of file shape_info.h.

◆ n_q_points_1d

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

Stores the number of quadrature points per dimension.

Definition at line 200 of file shape_info.h.

◆ n_q_points

template<typename Number>
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_q_points

Stores the number of quadrature points in dim dimensions for a cell.

Definition at line 206 of file shape_info.h.

◆ dofs_per_cell

template<typename Number>
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::dofs_per_cell

Stores the number of DoFs per cell in dim dimensions.

Definition at line 211 of file shape_info.h.

◆ n_q_points_face

template<typename Number>
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::n_q_points_face

Stores the number of quadrature points per face in dim dimensions.

Definition at line 216 of file shape_info.h.

◆ dofs_per_face

template<typename Number>
unsigned int internal::MatrixFreeFunctions::ShapeInfo< Number >::dofs_per_face

Stores the number of DoFs per face in dim dimensions.

Definition at line 221 of file shape_info.h.


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