16#ifndef dealii_matrix_free_shape_info_h
17#define dealii_matrix_free_shape_info_h
33template <
int dim,
int spacedim>
39 namespace MatrixFreeFunctions
132 template <
typename Number>
155 template <
int dim,
int spacedim>
159 const std::vector<unsigned int> &lexicographic,
160 const unsigned int direction);
168 template <
int dim,
int spacedim>
172 const std::vector<unsigned int> &lexicographic,
173 const unsigned int direction);
338 std::array<
AlignedVector<typename ::internal::VectorizedArrayTrait<
339 Number>::value_type>,
392 template <
typename Number>
410 template <
int dim,
int spacedim,
int dim_q>
413 const unsigned int base_element = 0);
423 template <
int dim,
int spacedim,
int dim_q>
427 const unsigned int base_element = 0);
445 template <
int dim,
int spacedim>
464 const unsigned int component = 0)
const;
485 std::vector<UnivariateShapeData<Number>>
data;
631 template <
typename Number>
634 const unsigned int component)
const
640 return *(data_access(dimension, component));
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
@ tensor_symmetric_no_collocation
@ tensor_symmetric_collocation
@ tensor_symmetric_plus_dg0
@ tensor_symmetric_hermite
unsigned int n_q_points_face
unsigned int n_dimensions
unsigned int dofs_per_component_on_cell
std::size_t memory_consumption() const
static Table< 2, unsigned int > compute_orientation_table(const unsigned int n_points_per_dim)
ShapeInfo(const Quadrature< dim_q > &quad, const FiniteElement< dim, spacedim > &fe, const unsigned int base_element=0)
::Table< 2, unsigned int > face_to_cell_index_nodal
void reinit(const Quadrature< dim_q > &quad, const FiniteElement< dim, spacedim > &fe_dim, const unsigned int base_element=0)
static bool is_supported(const FiniteElement< dim, spacedim > &fe)
::Table< 2, unsigned int > face_orientations_dofs
std::vector< UnivariateShapeData< Number > > data
::Table< 2, unsigned int > face_orientations_quad
const UnivariateShapeData< Number > & get_shape_data(const unsigned int dimension=0, const unsigned int component=0) const
unsigned int n_components
std::vector< unsigned int > lexicographic_numbering
::Table< 2, unsigned int > face_to_cell_index_hermite
unsigned int dofs_per_component_on_face
std::vector< unsigned int > n_q_points_faces
::Table< 2, UnivariateShapeData< Number > * > data_access
AlignedVector< Number > inverse_shape_values
void evaluate_collocation_space(const FiniteElement< dim, spacedim > &fe, const Quadrature< 1 > &quad, const std::vector< unsigned int > &lexicographic, const unsigned int direction)
AlignedVector< Number > shape_values
AlignedVector< Number > shape_hessians_collocation
std::array< AlignedVector< Number >, 2 > shape_data_on_face
AlignedVector< Number > shape_values_eo
AlignedVector< Number > inverse_shape_values_eo
std::array< AlignedVector< typename ::internal::VectorizedArrayTrait< Number >::value_type >, 2 > subface_interpolation_matrices_scalar
AlignedVector< Number > shape_hessians_eo
bool check_and_set_shapes_symmetric()
std::size_t memory_consumption() const
Table< 3, Number > shape_gradients_face
AlignedVector< Number > shape_gradients_collocation_eo
unsigned int n_q_points_1d
AlignedVector< Number > shape_gradients_eo
std::array< AlignedVector< Number >, 2 > subface_interpolation_matrices
bool check_shapes_collocation() const
void evaluate_shape_functions(const FiniteElement< dim, spacedim > &fe, const Quadrature< 1 > &quad, const std::vector< unsigned int > &lexicographic, const unsigned int direction)
AlignedVector< Number > shape_hessians
AlignedVector< Number > shape_gradients
std::array< AlignedVector< Number >, 2 > hessians_within_subface
std::array< AlignedVector< Number >, 2 > values_within_subface
std::array< AlignedVector< Number >, 2 > quadrature_data_on_face
AlignedVector< Number > shape_gradients_collocation
Quadrature< 1 > quadrature
bool nodal_at_cell_boundaries
Table< 3, Number > shape_values_face
AlignedVector< Number > shape_hessians_collocation_eo
std::array< AlignedVector< Number >, 2 > gradients_within_subface