Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
FE_Poly< PolynomialType, dim, spacedim > Class Template Reference

#include <deal.II/fe/fe_poly.h>

Inheritance diagram for FE_Poly< PolynomialType, dim, spacedim >:
[legend]

Classes

class  InternalData
 

Public Member Functions

 FE_Poly (const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
unsigned int get_degree () const
 
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
std::vector< unsigned intget_poly_space_numbering () const
 
std::vector< unsigned intget_poly_space_numbering_inverse () const
 
virtual double shape_value (const unsigned int i, const Point< dim > &p) const override
 
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 3, dim > shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 3, dim > shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual Tensor< 4, dim > shape_4th_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 4, dim > shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const override
 
virtual std::size_t memory_consumption () const override
 
- Public Member Functions inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
 FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
 FiniteElement (FiniteElement< dim, spacedim > &&)=default
 
 FiniteElement (const FiniteElement< dim, spacedim > &)=default
 
virtual ~FiniteElement () override=default
 
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned intoperator^ (const unsigned int multiplicity) const
 
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone () const=0
 
virtual std::string get_name () const=0
 
const FiniteElement< dim, spacedim > & operator[] (const unsigned int fe_index) const
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
 
virtual const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
virtual const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
bool prolongation_is_implemented () const
 
bool isotropic_prolongation_is_implemented () const
 
bool restriction_is_implemented () const
 
bool isotropic_restriction_is_implemented () const
 
bool restriction_is_additive (const unsigned int index) const
 
const FullMatrix< double > & constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
bool constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
virtual bool hp_constraints_are_implemented () const
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual FiniteElementDomination::Domination compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const final
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
 
virtual bool operator== (const FiniteElement< dim, spacedim > &fe) const
 
bool operator!= (const FiniteElement< dim, spacedim > &) const
 
std::pair< unsigned int, unsigned intsystem_to_component_index (const unsigned int index) const
 
unsigned int component_to_system_index (const unsigned int component, const unsigned int index) const
 
std::pair< unsigned int, unsigned intface_system_to_component_index (const unsigned int index) const
 
unsigned int adjust_quad_dof_index_for_face_orientation (const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
 
virtual unsigned int face_to_cell_index (const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
 
unsigned int adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const
 
const ComponentMaskget_nonzero_components (const unsigned int i) const
 
unsigned int n_nonzero_components (const unsigned int i) const
 
bool is_primitive () const
 
bool is_primitive (const unsigned int i) const
 
unsigned int n_base_elements () const
 
virtual const FiniteElement< dim, spacedim > & base_element (const unsigned int index) const
 
unsigned int element_multiplicity (const unsigned int index) const
 
const FiniteElement< dim, spacedim > & get_sub_fe (const ComponentMask &mask) const
 
virtual const FiniteElement< dim, spacedim > & get_sub_fe (const unsigned int first_component, const unsigned int n_selected_components) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned intsystem_to_base_index (const unsigned int index) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned intface_system_to_base_index (const unsigned int index) const
 
types::global_dof_index first_block_of_base (const unsigned int b) const
 
std::pair< unsigned int, unsigned intcomponent_to_base_index (const unsigned int component) const
 
std::pair< unsigned int, unsigned intblock_to_base_index (const unsigned int block) const
 
std::pair< unsigned int, types::global_dof_indexsystem_to_block_index (const unsigned int component) const
 
unsigned int component_to_block_index (const unsigned int component) const
 
ComponentMask component_mask (const FEValuesExtractors::Scalar &scalar) const
 
ComponentMask component_mask (const FEValuesExtractors::Vector &vector) const
 
ComponentMask component_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
ComponentMask component_mask (const BlockMask &block_mask) const
 
BlockMask block_mask (const FEValuesExtractors::Scalar &scalar) const
 
BlockMask block_mask (const FEValuesExtractors::Vector &vector) const
 
BlockMask block_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
BlockMask block_mask (const ComponentMask &component_mask) const
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const
 
const std::vector< Point< dim > > & get_unit_support_points () const
 
bool has_support_points () const
 
virtual Point< dim > unit_support_point (const unsigned int index) const
 
const std::vector< Point< dim - 1 > > & get_unit_face_support_points () const
 
bool has_face_support_points () const
 
virtual Point< dim - 1 > unit_face_support_point (const unsigned int index) const
 
const std::vector< Point< dim > > & get_generalized_support_points () const
 
bool has_generalized_support_points () const
 
GeometryPrimitive get_associated_geometry_primitive (const unsigned int cell_dof_index) const
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const
 
virtual std::size_t memory_consumption () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from FiniteElementData< dim >
 FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
 
unsigned int n_dofs_per_vertex () const
 
unsigned int n_dofs_per_line () const
 
unsigned int n_dofs_per_quad () const
 
unsigned int n_dofs_per_hex () const
 
unsigned int n_dofs_per_face () const
 
unsigned int n_dofs_per_cell () const
 
template<int structdim>
unsigned int n_dofs_per_object () const
 
unsigned int n_components () const
 
unsigned int n_blocks () const
 
const BlockIndicesblock_indices () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Protected Member Functions

virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBaseget_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
virtual void fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
 
void correct_hessians (internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
 
void correct_third_derivatives (internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
 
void fill_fe_values (const Triangulation< 1, 2 >::cell_iterator &, const CellSimilarity::Similarity cell_similarity, const Quadrature< 1 > &quadrature, const Mapping< 1, 2 > &mapping, const Mapping< 1, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 1, 2 > &mapping_data, const FiniteElement< 1, 2 >::InternalDataBase &fe_internal,::internal::FEValuesImplementation::FiniteElementRelatedData< 1, 2 > &output_data) const
 
void fill_fe_values (const Triangulation< 2, 3 >::cell_iterator &, const CellSimilarity::Similarity cell_similarity, const Quadrature< 2 > &quadrature, const Mapping< 2, 3 > &mapping, const Mapping< 2, 3 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 3 > &mapping_data, const FiniteElement< 2, 3 >::InternalDataBase &fe_internal,::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 3 > &output_data) const
 
- Protected Member Functions inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
void reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
 
TableIndices< 2 > interface_constraints_size () const
 
virtual std::unique_ptr< InternalDataBaseget_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
virtual std::unique_ptr< InternalDataBaseget_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
virtual void fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0
 
virtual void fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0
 
virtual void fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0
 

Protected Attributes

PolynomialType poly_space
 
- Protected Attributes inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
std::vector< std::vector< FullMatrix< double > > > restriction
 
std::vector< std::vector< FullMatrix< double > > > prolongation
 
FullMatrix< doubleinterface_constraints
 
std::vector< Point< dim > > unit_support_points
 
std::vector< Point< dim - 1 > > unit_face_support_points
 
std::vector< Point< dim > > generalized_support_points
 
std::vector< Point< dim - 1 > > generalized_face_support_points
 
Table< 2, intadjust_quad_dof_index_for_face_orientation_table
 
std::vector< intadjust_line_dof_index_for_line_orientation_table
 
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
 
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
 
BlockIndices base_to_block_indices
 
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
 
const std::vector< boolrestriction_is_additive_flags
 
const std::vector< ComponentMasknonzero_components
 
const std::vector< unsigned intn_nonzero_components_table
 
const bool cached_primitivity
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
enum  Conformity {
  unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04,
  H1 = Hcurl | Hdiv, H2 = 0x0e
}
 
- Static Public Member Functions inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
static ::ExceptionBaseExcShapeFunctionNotPrimitive (int arg1)
 
static ::ExceptionBaseExcFENotPrimitive ()
 
static ::ExceptionBaseExcUnitShapeValuesDoNotExist ()
 
static ::ExceptionBaseExcFEHasNoSupportPoints ()
 
static ::ExceptionBaseExcEmbeddingVoid ()
 
static ::ExceptionBaseExcProjectionVoid ()
 
static ::ExceptionBaseExcWrongInterfaceMatrixSize (int arg1, int arg2)
 
static ::ExceptionBaseExcInterpolationNotImplemented ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Public Attributes inherited from FiniteElementData< dim >
const unsigned int dofs_per_vertex
 
const unsigned int dofs_per_line
 
const unsigned int dofs_per_quad
 
const unsigned int dofs_per_hex
 
const unsigned int first_line_index
 
const unsigned int first_quad_index
 
const unsigned int first_hex_index
 
const unsigned int first_face_line_index
 
const unsigned int first_face_quad_index
 
const unsigned int dofs_per_face
 
const unsigned int dofs_per_cell
 
const unsigned int components
 
const unsigned int degree
 
const Conformity conforming_space
 
const BlockIndices block_indices_data
 
- Static Public Attributes inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
static const unsigned int space_dimension
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Static Protected Member Functions inherited from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >
static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 

Detailed Description

template<class PolynomialType, int dim = PolynomialType::dimension, int spacedim = dim>
class FE_Poly< PolynomialType, dim, spacedim >

This class gives a unified framework for the implementation of FiniteElement classes based on scalar polynomial spaces like the TensorProductPolynomials or PolynomialSpace classes. This class has a corresponding class for tensor-valued finite elements in the FE_PolyTensor class.

Every class that has the following public member variables and functions can be used as template parameter PolynomialType.

static const unsigned int dimension;
void evaluate (const Point<dim> &unit_point,
std::vector<double> &values,
std::vector<Tensor<1,dim> > &grads,
std::vector<Tensor<2,dim> > &grad_grads,
std::vector<Tensor<3,dim> > &third_derivatives,
std::vector<Tensor<4,dim> > &fourth_derivatives) const;
double compute_value (const unsigned int i,
const Point<dim> &p) const;
template <int order>
Tensor<order,dim> compute_derivative (const unsigned int i,
const Point<dim> &p) const;

Example classes are TensorProductPolynomials, PolynomialSpace or PolynomialsP.

This class is not a fully implemented FiniteElement class. Instead there are several pure virtual functions declared in the FiniteElement and FiniteElement classes which cannot be implemented by this class but are left for implementation in derived classes.

Todo:
Since nearly all functions for spacedim != dim are specialized, this class needs cleaning up.
Author
Ralf Hartmann 2004, Guido Kanschat, 2009

Definition at line 76 of file fe_poly.h.

Constructor & Destructor Documentation

◆ FE_Poly()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
FE_Poly< PolynomialType, dim, spacedim >::FE_Poly ( const PolynomialType poly_space,
const FiniteElementData< dim > &  fe_data,
const std::vector< bool > &  restriction_is_additive_flags,
const std::vector< ComponentMask > &  nonzero_components 
)

Constructor.

Member Function Documentation

◆ get_degree()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
unsigned int FE_Poly< PolynomialType, dim, spacedim >::get_degree ( ) const

Return the polynomial degree of this finite element, i.e. the value passed to the constructor.

◆ requires_update_flags()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual UpdateFlags FE_Poly< PolynomialType, dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
overridevirtual

Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.

As an example, if update_flags contains update_gradients a finite element class will typically require the computation of the inverse of the Jacobian matrix in order to rotate the gradient of shape functions on the reference cell to the real cell. It would then return not just update_gradients, but also update_covariant_transformation, the flag that makes the mapping class produce the inverse of the Jacobian matrix.

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.

See also
UpdateFlags

Implements FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ get_poly_space_numbering()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
std::vector<unsigned int> FE_Poly< PolynomialType, dim, spacedim >::get_poly_space_numbering ( ) const

Return the numbering of the underlying polynomial space compared to lexicographic ordering of the basis functions. Returns PolynomialType::get_numbering().

◆ get_poly_space_numbering_inverse()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
std::vector<unsigned int> FE_Poly< PolynomialType, dim, spacedim >::get_poly_space_numbering_inverse ( ) const

Return the inverse numbering of the underlying polynomial space. Returns PolynomialType::get_numbering_inverse().

◆ shape_value()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual double FE_Poly< PolynomialType, dim, spacedim >::shape_value ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the value of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_value_component()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual double FE_Poly< PolynomialType, dim, spacedim >::shape_value_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the value of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is scalar, the returned value is the same as if the function without the _component suffix were called, provided that the specified component is zero.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_grad()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<1, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the gradient of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_grad_component()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<1, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the gradient of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is scalar, the returned value is the same as if the function without the _component suffix were called, provided that the specified component is zero.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_grad_grad()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<2, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of second derivatives of the ith shape function at point p on the unit cell. See the FiniteElement base class for more information about the semantics of this function.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_grad_grad_component()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<2, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_grad_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the second derivative of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is scalar, the returned value is the same as if the function without the _component suffix were called, provided that the specified component is zero.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_3rd_derivative()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<3, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_3rd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of third derivatives of the ith shape function at point p on the unit cell. See the FiniteElement base class for more information about the semantics of this function.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_3rd_derivative_component()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<3, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_3rd_derivative_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the third derivative of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is scalar, the returned value is the same as if the function without the _component suffix were called, provided that the specified component is zero.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_4th_derivative()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<4, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_4th_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Return the tensor of fourth derivatives of the ith shape function at point p on the unit cell. See the FiniteElement base class for more information about the semantics of this function.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ shape_4th_derivative_component()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual Tensor<4, dim> FE_Poly< PolynomialType, dim, spacedim >::shape_4th_derivative_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Return the fourth derivative of the componentth vector component of the ith shape function at the point p. See the FiniteElement base class for more information about the semantics of this function.

Since this element is scalar, the returned value is the same as if the function without the _component suffix were called, provided that the specified component is zero.

Reimplemented from FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

◆ memory_consumption()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual std::size_t FE_Poly< PolynomialType, dim, spacedim >::memory_consumption ( ) const
overridevirtual

Return an estimate (in bytes) for the memory consumption of this object.

◆ get_data()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual std::unique_ptr< typename FiniteElement<dim, spacedim>::InternalDataBase> FE_Poly< PolynomialType, dim, spacedim >::get_data ( const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
inlineoverrideprotectedvirtual

Create an internal data object and return a pointer to it of which the caller of this function then assumes ownership. This object will then be passed to the FiniteElement::fill_fe_values() every time the finite element shape functions and their derivatives are evaluated on a concrete cell. The object created here is therefore used by derived classes as a place for scratch objects that are used in evaluating shape functions, as well as to store information that can be pre-computed once and re-used on every cell (e.g., for evaluating the values and gradients of shape functions on the reference cell, for later re-use when transforming these values to a concrete cell).

This function is the first one called in the process of initializing a FEValues object for a given mapping and finite element object. The returned object will later be passed to FiniteElement::fill_fe_values() for a concrete cell, which will itself place its output into an object of type internal::FEValuesImplementation::FiniteElementRelatedData. Since there may be data that can already be computed in its final form on the reference cell, this function also receives a reference to the internal::FEValuesImplementation::FiniteElementRelatedData object as its last argument. This output argument is guaranteed to always be the same one when used with the InternalDataBase object returned by this function. In other words, the subdivision of scratch data and final data in the returned object and the output_data object is as follows: If data can be pre- computed on the reference cell in the exact form in which it will later be needed on a concrete cell, then this function should already emplace it in the output_data object. An example are the values of shape functions at quadrature points for the usual Lagrange elements which on a concrete cell are identical to the ones on the reference cell. On the other hand, if some data can be pre-computed to make computations on a concrete cell cheaper, then it should be put into the returned object for later re-use in a derive class's implementation of FiniteElement::fill_fe_values(). An example are the gradients of shape functions on the reference cell for Lagrange elements: to compute the gradients of the shape functions on a concrete cell, one has to multiply the gradients on the reference cell by the inverse of the Jacobian of the mapping; consequently, we cannot already compute the gradients on a concrete cell at the time the current function is called, but we can at least pre-compute the gradients on the reference cell, and store it in the object returned.

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module. See also the documentation of the InternalDataBase class.

Parameters
[in]update_flagsA set of UpdateFlags values that describe what kind of information the FEValues object requests the finite element to compute. This set of flags may also include information that the finite element can not compute, e.g., flags that pertain to data produced by the mapping. An implementation of this function needs to set up all data fields in the returned object that are necessary to produce the finite- element related data specified by these flags, and may already pre- compute part of this information as discussed above. Elements may want to store these update flags (or a subset of these flags) in InternalDataBase::update_each so they know at the time when FiniteElement::fill_fe_values() is called what they are supposed to compute
[in]mappingA reference to the mapping used for computing values and derivatives of shape functions.
[in]quadratureA reference to the object that describes where the shape functions should be evaluated.
[out]output_dataA reference to the object that FEValues will use in conjunction with the object returned here and where an implementation of FiniteElement::fill_fe_values() will place the requested information. This allows the current function to already pre-compute pieces of information that can be computed on the reference cell, as discussed above. FEValues guarantees that this output object and the object returned by the current function will always be used together.
Returns
A pointer to an object of a type derived from InternalDataBase and that derived classes can use to store scratch data that can be pre- computed, or for scratch arrays that then only need to be allocated once. The calling site assumes ownership of this object and will delete it when it is no longer necessary.

Implements FiniteElement< PolynomialType::dimension, PolynomialType::dimension >.

Definition at line 246 of file fe_poly.h.

◆ fill_fe_values() [1/3]

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual void FE_Poly< PolynomialType, dim, spacedim >::fill_fe_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< dim > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

◆ fill_fe_face_values()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual void FE_Poly< PolynomialType, dim, spacedim >::fill_fe_face_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const Quadrature< dim - 1 > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

◆ fill_fe_subface_values()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
virtual void FE_Poly< PolynomialType, dim, spacedim >::fill_fe_subface_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const Quadrature< dim - 1 > &  quadrature,
const Mapping< dim, spacedim > &  mapping,
const typename Mapping< dim, spacedim >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const typename FiniteElement< dim, spacedim >::InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

◆ correct_hessians()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
void FE_Poly< PolynomialType, dim, spacedim >::correct_hessians ( internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data,
const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const unsigned int  n_q_points 
) const
protected

Correct the shape Hessians by subtracting the terms corresponding to the Jacobian pushed forward gradient.

Before the correction, the Hessians would be given by

\[ D_{ijk} = \frac{d^2\phi_i}{d \hat x_J d \hat x_K} (J_{jJ})^{-1} (J_{kK})^{-1}, \]

where \(J_{iI}=\frac{d x_i}{d \hat x_I}\). After the correction, the correct Hessians would be given by

\[ \frac{d^2 \phi_i}{d x_j d x_k} = D_{ijk} - H_{mjk} \frac{d \phi_i}{d x_m}, \]

where \(H_{ijk}\) is the Jacobian pushed-forward derivative.

◆ correct_third_derivatives()

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
void FE_Poly< PolynomialType, dim, spacedim >::correct_third_derivatives ( internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data,
const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  mapping_data,
const unsigned int  n_q_points 
) const
protected

Correct the shape third derivatives by subtracting the terms corresponding to the Jacobian pushed forward gradient and second derivative.

Before the correction, the third derivatives would be given by

\[ D_{ijkl} = \frac{d^3\phi_i}{d \hat x_J d \hat x_K d \hat x_L} (J_{jJ})^{-1} (J_{kK})^{-1} (J_{lL})^{-1}, \]

where \(J_{iI}=\frac{d x_i}{d \hat x_I}\). After the correction, the correct third derivative would be given by

\[ \frac{d^3\phi_i}{d x_j d x_k d x_l} = D_{ijkl} - H_{mjl} \frac{d^2 \phi_i}{d x_k d x_m} - H_{mkl} \frac{d^2 \phi_i}{d x_j d x_m} - H_{mjk} \frac{d^2 \phi_i}{d x_l d x_m} - K_{mjkl} \frac{d \phi_i}{d x_m}, \]

where \(H_{ijk}\) is the Jacobian pushed-forward derivative and \(K_{ijkl}\) is the Jacobian pushed-forward second derivative.

◆ fill_fe_values() [2/3]

void FE_Poly< TensorProductPolynomials< 1 >, 1, 2 >::fill_fe_values ( const Triangulation< 1, 2 >::cell_iterator &  ,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< 1 > &  quadrature,
const Mapping< 1, 2 > &  mapping,
const Mapping< 1, 2 >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< 1, 2 > &  mapping_data,
const FiniteElement< 1, 2 >::InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< 1, 2 > &  output_data 
) const
protected

Definition at line 34 of file fe_poly.cc.

◆ fill_fe_values() [3/3]

void FE_Poly< TensorProductPolynomials< 2 >, 2, 3 >::fill_fe_values ( const Triangulation< 2, 3 >::cell_iterator &  ,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< 2 > &  quadrature,
const Mapping< 2, 3 > &  mapping,
const Mapping< 2, 3 >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< 2, 3 > &  mapping_data,
const FiniteElement< 2, 3 >::InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 3 > &  output_data 
) const
protected

Definition at line 102 of file fe_poly.cc.

Member Data Documentation

◆ poly_space

template<class PolynomialType , int dim = PolynomialType::dimension, int spacedim = dim>
PolynomialType FE_Poly< PolynomialType, dim, spacedim >::poly_space
protected

The polynomial space. Its type is given by the template parameter PolynomialType.

Definition at line 516 of file fe_poly.h.


The documentation for this class was generated from the following file:
FiniteElementData::dimension
static const unsigned int dimension
Definition: fe_base.h:220
Tensor< 1, dim >
internal::PolynomialsRannacherTurekImplementation::compute_derivative
Tensor< order, dim > compute_derivative(const unsigned int, const Point< dim > &)
Definition: polynomials_rannacher_turek.h:118
Point< dim >