Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
FE_Q_Hierarchical< dim > Class Template Reference

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

Inheritance diagram for FE_Q_Hierarchical< dim >:
[legend]

Public Member Functions

 FE_Q_Hierarchical (const unsigned int p)
 
virtual std::string get_name () const override
 
virtual std::unique_ptr< FiniteElement< dim, dim > > clone () const override
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const override
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
 
virtual std::size_t memory_consumption () const override
 
std::vector< unsigned int > get_embedding_dofs (const unsigned int sub_degree) const
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const override
 
template<>
bool has_support_on_face (const unsigned int, const unsigned int) const
 
template<>
bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
 
Functions to support hp
virtual bool hp_constraints_are_implemented () const override
 
virtual void get_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const override
 
virtual const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities (const FiniteElement< dim > &fe_other) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities (const FiniteElement< dim > &fe_other) const override
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities (const FiniteElement< dim > &fe_other) const override
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
 
- Public Member Functions inherited from FE_Poly< TensorProductPolynomials< dim >, dim >
 FE_Poly (const TensorProductPolynomials< dim > &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
 
std::vector< unsigned int > get_poly_space_numbering () const
 
std::vector< unsigned int > get_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
 
- Public Member Functions inherited from FiniteElement< dim, dim >
 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 int > operator^ (const unsigned int multiplicity) const
 
const FiniteElement< dim, spacedim > & operator[] (const unsigned int fe_index) const
 
virtual bool operator== (const FiniteElement< dim, spacedim > &fe) const
 
bool operator!= (const FiniteElement< dim, spacedim > &) const
 
virtual const FullMatrix< double > & get_restriction_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 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
 
std::pair< unsigned int, unsigned int > system_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 int > face_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 int > system_to_base_index (const unsigned int index) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_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 int > component_to_base_index (const unsigned int component) const
 
std::pair< unsigned int, unsigned int > block_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
 
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
 
const std::vector< Point< dim - 1 > > & get_generalized_face_support_points () const
 
bool has_generalized_face_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
 
- 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
 

Private Member Functions

void build_dofs_cell (std::vector< FullMatrix< double >> &dofs_cell, std::vector< FullMatrix< double >> &dofs_subcell) const
 
void initialize_constraints (const std::vector< FullMatrix< double >> &dofs_subcell)
 
void initialize_embedding_and_restriction (const std::vector< FullMatrix< double >> &dofs_cell, const std::vector< FullMatrix< double >> &dofs_subcell)
 
void initialize_generalized_support_points ()
 
void initialize_generalized_face_support_points ()
 

Static Private Member Functions

static std::vector< unsigned int > get_dpo_vector (const unsigned int degree)
 
static std::vector< unsigned int > hierarchic_to_fe_q_hierarchical_numbering (const FiniteElementData< dim > &fe)
 
static std::vector< unsigned int > face_fe_q_hierarchical_to_hierarchic_numbering (const unsigned int degree)
 

Private Attributes

const std::vector< unsigned int > face_renumber
 

Friends

template<int dim1>
class FE_Q_Hierarchical
 

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< dim, dim >
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< dim, dim >
static const unsigned int space_dimension
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Protected Member Functions inherited from FE_Poly< TensorProductPolynomials< dim >, dim >
void correct_third_derivatives (internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const unsigned int n_q_points, const unsigned int dof) const
 
- Protected Member Functions inherited from FiniteElement< dim, dim >
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< InternalDataBase > get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0
 
virtual std::unique_ptr< InternalDataBase > get_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< InternalDataBase > get_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
 
- Static Protected Member Functions inherited from FiniteElement< dim, dim >
static std::vector< unsigned int > compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FE_Poly< TensorProductPolynomials< dim >, dim >
TensorProductPolynomials< dim > poly_space
 
- Protected Attributes inherited from FiniteElement< dim, dim >
std::vector< std::vector< FullMatrix< double > > > restriction
 
std::vector< std::vector< FullMatrix< double > > > prolongation
 
FullMatrix< double > interface_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, int > adjust_quad_dof_index_for_face_orientation_table
 
std::vector< int > adjust_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< bool > restriction_is_additive_flags
 
const std::vector< ComponentMasknonzero_components
 
const std::vector< unsigned int > n_nonzero_components_table
 
const bool cached_primitivity
 

Detailed Description

template<int dim>
class FE_Q_Hierarchical< dim >

Implementation of hierarchical Qp shape functions that yield the finite element space of continuous, piecewise polynomials of degree p. This class is realized using tensor product polynomials based on a hierarchical basis Polynomials::Hierarchical on the interval [0,1] which is suitable for building an hp tensor product finite element if we assume that each element has a single degree.

The constructor of this class takes the degree p of this finite element.

This class is not implemented for the codimension one case (spacedim != dim).

Implementation

The constructor creates a TensorProductPolynomials object that includes the tensor product of Hierarchical polynomials of degree p. This TensorProductPolynomials object provides all values and derivatives of the shape functions.

Numbering of the degrees of freedom (DoFs)

The original ordering of the shape functions represented by the TensorProductPolynomials is a tensor product numbering. However, the shape functions on a cell are renumbered beginning with the shape functions whose support points are at the vertices, then on the line, on the quads, and finally (for 3d) on the hexes. To be explicit, these numberings are listed in the following:

Q1 elements

The \(Q_1^H\) element is of polynomial degree one and, consequently, is exactly the same as the \(Q_1\) element in class FE_Q. In particular, the shape function are defined in the exact same way:

In 2d, these shape functions look as follows:

Q1H_shape0000.png

Q1H_shape0001.png

\(Q_1^H\) element, shape function 0

\(Q_1^H\) element, shape function 1

Q1H_shape0002.png

Q1H_shape0003.png

\(Q_1^H\) element, shape function 2

\(Q_1^H\) element, shape function 3

Q2 elements

In 2d, these shape functions look as follows (the black plane corresponds to zero; negative shape function values may not be visible):

Q2H_shape0000.png

Q2H_shape0001.png

\(Q_2^H\) element, shape function 0

\(Q_2^H\) element, shape function 1

Q2H_shape0002.png

Q2H_shape0003.png

\(Q_2^H\) element, shape function 2

\(Q_2^H\) element, shape function 3

Q2H_shape0004.png

Q2H_shape0005.png

\(Q_2^H\) element, shape function 4

\(Q_2^H\) element, shape function 5

Q2H_shape0006.png

Q2H_shape0007.png

\(Q_2^H\) element, shape function 6

\(Q_2^H\) element, shape function 7

Q2H_shape0008.png

\(Q_2^H\) element, shape function 8

Q3 elements

In 2d, these shape functions look as follows (the black plane corresponds to zero; negative shape function values may not be visible):

Q3H_shape0000.png

Q3H_shape0001.png

\(Q_3^H\) element, shape function 0

\(Q_3^H\) element, shape function 1

Q3H_shape0002.png

Q3H_shape0003.png

\(Q_3^H\) element, shape function 2

\(Q_3^H\) element, shape function 3

Q3H_shape0004.png

Q3H_shape0005.png

\(Q_3^H\) element, shape function 4

\(Q_3^H\) element, shape function 5

Q3H_shape0006.png

Q3H_shape0007.png

\(Q_3^H\) element, shape function 6

\(Q_3^H\) element, shape function 7

Q3H_shape0008.png

Q3H_shape0009.png

\(Q_3^H\) element, shape function 8

\(Q_3^H\) element, shape function 9

Q3H_shape0010.png

Q3H_shape0011.png

\(Q_3^H\) element, shape function 10

\(Q_3^H\) element, shape function 11

Q3H_shape0012.png

Q3H_shape0013.png

\(Q_3^H\) element, shape function 12

\(Q_3^H\) element, shape function 13

Q3H_shape0014.png

Q3H_shape0015.png

\(Q_3^H\) element, shape function 14

\(Q_3^H\) element, shape function 15

Q4 elements

In 2d, these shape functions look as follows (the black plane corresponds to zero; negative shape function values may not be visible):

Q4H_shape0000.png

Q4H_shape0001.png

\(Q_4^H\) element, shape function 0

\(Q_4^H\) element, shape function 1

Q4H_shape0002.png

Q4H_shape0003.png

\(Q_4^H\) element, shape function 2

\(Q_4^H\) element, shape function 3

Q4H_shape0004.png

Q4H_shape0005.png

\(Q_4^H\) element, shape function 4

\(Q_4^H\) element, shape function 5

Q4H_shape0006.png

Q4H_shape0007.png

\(Q_4^H\) element, shape function 6

\(Q_4^H\) element, shape function 7

Q4H_shape0008.png

Q4H_shape0009.png

\(Q_4^H\) element, shape function 8

\(Q_4^H\) element, shape function 9

Q4H_shape0010.png

Q4H_shape0011.png

\(Q_4^H\) element, shape function 10

\(Q_4^H\) element, shape function 11

Q4H_shape0012.png

Q4H_shape0013.png

\(Q_4^H\) element, shape function 12

\(Q_4^H\) element, shape function 13

Q4H_shape0014.png

Q4H_shape0015.png

\(Q_4^H\) element, shape function 14

\(Q_4^H\) element, shape function 15

Q4H_shape0016.png

Q4H_shape0017.png

\(Q_4^H\) element, shape function 16

\(Q_4^H\) element, shape function 17

Q4H_shape0018.png

Q4H_shape0019.png

\(Q_4^H\) element, shape function 18

\(Q_4^H\) element, shape function 19

Q4H_shape0020.png

Q4H_shape0021.png

\(Q_4^H\) element, shape function 20

\(Q_4^H\) element, shape function 21

Q4H_shape0022.png

Q4H_shape0023.png

\(Q_4^H\) element, shape function 22

\(Q_4^H\) element, shape function 23

Q4H_shape0024.png

\(Q_4^H\) element, shape function 24

Author
Brian Carnes, 2002, Ralf Hartmann 2004, 2005, Denis Davydov, 2015

Definition at line 543 of file fe_q_hierarchical.h.

Constructor & Destructor Documentation

◆ FE_Q_Hierarchical()

template<int dim>
FE_Q_Hierarchical< dim >::FE_Q_Hierarchical ( const unsigned int  p)

Constructor for tensor product polynomials of degree p.

Definition at line 61 of file fe_q_hierarchical.cc.

Member Function Documentation

◆ get_name()

template<int dim>
std::string FE_Q_Hierarchical< dim >::get_name ( ) const
overridevirtual

Return a string that uniquely identifies a finite element. This class returns FE_Q_Hierarchical<dim>(degree), with dim and degree replaced by appropriate values.

Implements FiniteElement< dim, dim >.

Definition at line 121 of file fe_q_hierarchical.cc.

◆ clone()

template<int dim>
std::unique_ptr< FiniteElement< dim, dim > > FE_Q_Hierarchical< dim >::clone ( ) const
overridevirtual

A sort of virtual copy constructor, this function returns a copy of the finite element object. Derived classes need to override the function here in this base class and return an object of the same type as the derived class.

Some places in the library, for example the constructors of FESystem as well as the hp::FECollection class, need to make copies of finite elements without knowing their exact type. They do so through this function.

Implements FiniteElement< dim, dim >.

Definition at line 140 of file fe_q_hierarchical.cc.

◆ has_support_on_face() [1/3]

template<int dim>
bool FE_Q_Hierarchical< dim >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const
overridevirtual

This function returns true, if the shape function shape_index has non-zero function values somewhere on the face face_index.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 2166 of file fe_q_hierarchical.cc.

◆ hp_constraints_are_implemented()

template<int dim>
bool FE_Q_Hierarchical< dim >::hp_constraints_are_implemented ( ) const
overridevirtual

Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible".

For the FE_Q_Hierarchical class the result is always true (independent of the degree of the element), as it implements the complete set of functions necessary for hp capability.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 220 of file fe_q_hierarchical.cc.

◆ get_interpolation_matrix()

template<int dim>
void FE_Q_Hierarchical< dim >::get_interpolation_matrix ( const FiniteElement< dim > &  source,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from the given finite element to the present one. Interpolation only between FE_Q_Hierarchical is supported.

Definition at line 149 of file fe_q_hierarchical.cc.

◆ get_prolongation_matrix()

template<int dim>
const FullMatrix< double > & FE_Q_Hierarchical< dim >::get_prolongation_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case = RefinementCase<dim>::isotropic_refinement 
) const
overridevirtual

Embedding matrix between grids. Only isotropic refinement is supported.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 200 of file fe_q_hierarchical.cc.

◆ hp_vertex_dof_identities()

template<int dim>
std::vector< std::pair< unsigned int, unsigned int > > FE_Q_Hierarchical< dim >::hp_vertex_dof_identities ( const FiniteElement< dim > &  fe_other) const
overridevirtual

If, on a vertex, several finite elements are active, the hp code first assigns the degrees of freedom of each of these FEs different global indices. It then calls this function to find out which of them should get identical values, and consequently can receive the same global DoF index. This function therefore returns a list of identities between DoFs of the present finite element object with the DoFs of fe_other, which is a reference to a finite element object representing one of the other finite elements active on this particular vertex. The function computes which of the degrees of freedom of the two finite element objects are equivalent, both numbered between zero and the corresponding value of dofs_per_vertex of the two finite elements. The first index of each pair denotes one of the vertex dofs of the present element, whereas the second is the corresponding index of the other finite element.

Definition at line 228 of file fe_q_hierarchical.cc.

◆ hp_line_dof_identities()

template<int dim>
std::vector< std::pair< unsigned int, unsigned int > > FE_Q_Hierarchical< dim >::hp_line_dof_identities ( const FiniteElement< dim > &  fe_other) const
overridevirtual

Same as above but for lines.

Definition at line 261 of file fe_q_hierarchical.cc.

◆ hp_quad_dof_identities()

template<int dim>
std::vector< std::pair< unsigned int, unsigned int > > FE_Q_Hierarchical< dim >::hp_quad_dof_identities ( const FiniteElement< dim > &  fe_other) const
overridevirtual

Same as above but for faces.

Definition at line 299 of file fe_q_hierarchical.cc.

◆ compare_for_domination()

template<int dim>
FiniteElementDomination::Domination FE_Q_Hierarchical< dim >::compare_for_domination ( const FiniteElement< dim > &  fe_other,
const unsigned int  codim = 0 
) const
finaloverridevirtual

Return whether this element dominates another one given as argument fe_other, whether it is the other way around, whether neither dominates, or if either could dominate. The codim parameter describes the codimension of the investigated subspace and specifies that it is subject to this comparison. For example, if codim==0 then this function compares which element dominates at the cell level. If codim==1, then the elements are compared at faces, i.e., the comparison happens between the function spaces of the two finite elements as restricted to a face. Larger values of codim work correspondingly.

For a definition of domination, see FiniteElementDomination::Domination and in particular the hp paper.

Definition at line 337 of file fe_q_hierarchical.cc.

◆ get_face_interpolation_matrix()

template<int dim>
void FE_Q_Hierarchical< dim >::get_face_interpolation_matrix ( const FiniteElement< dim > &  source,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from a face of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented.

Definition at line 915 of file fe_q_hierarchical.cc.

◆ get_subface_interpolation_matrix()

template<int dim>
void FE_Q_Hierarchical< dim >::get_subface_interpolation_matrix ( const FiniteElement< dim > &  source,
const unsigned int  subface,
FullMatrix< double > &  matrix 
) const
overridevirtual

Return the matrix interpolating from a face of one element to the subface of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type ExcInterpolationNotImplemented.

Definition at line 1001 of file fe_q_hierarchical.cc.

◆ memory_consumption()

template<int dim>
std::size_t FE_Q_Hierarchical< dim >::memory_consumption ( ) const
overridevirtual

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

This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 2411 of file fe_q_hierarchical.cc.

◆ get_embedding_dofs()

template<int dim>
std::vector< unsigned int > FE_Q_Hierarchical< dim >::get_embedding_dofs ( const unsigned int  sub_degree) const

For a finite element of degree sub_degree < degree, we return a vector which maps the numbering on an FE of degree sub_degree into the numbering on this element.

Definition at line 2254 of file fe_q_hierarchical.cc.

◆ get_constant_modes()

template<int dim>
std::pair< Table< 2, bool >, std::vector< unsigned int > > FE_Q_Hierarchical< dim >::get_constant_modes ( ) const
overridevirtual

Return a list of constant modes of the element. For this element, the list consists of true arguments for the first vertex shape functions and false for the remaining ones.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 2394 of file fe_q_hierarchical.cc.

◆ get_dpo_vector()

template<int dim>
std::vector< unsigned int > FE_Q_Hierarchical< dim >::get_dpo_vector ( const unsigned int  degree)
staticprivate

Only for internal use. Its full name is get_dofs_per_object_vector function and it creates the dofs_per_object vector that is needed within the constructor to be passed to the constructor of FiniteElementData.

Definition at line 1934 of file fe_q_hierarchical.cc.

◆ hierarchic_to_fe_q_hierarchical_numbering()

template<int dim>
std::vector< unsigned int > FE_Q_Hierarchical< dim >::hierarchic_to_fe_q_hierarchical_numbering ( const FiniteElementData< dim > &  fe)
staticprivate

The numbering of the degrees of freedom in continuous finite elements is hierarchic, i.e. in such a way that we first number the vertex dofs, in the order of the vertices as defined by the triangulation, then the line dofs in the order and respecting the direction of the lines, then the dofs on quads, etc.

The dofs associated with 1d hierarchical polynomials are ordered with the vertices first ( \(phi_0(x)=1-x\) and \(phi_1(x)=x\)) and then the line dofs (the higher degree polynomials). The 2d and 3d hierarchical polynomials originate from the 1d hierarchical polynomials by tensor product. In the following, the resulting numbering of dofs will be denoted by fe_q_hierarchical numbering.

This function constructs a table which fe_q_hierarchical index each degree of freedom in the hierarchic numbering would have.

This function is analogous to the FETools::hierarchic_to_lexicographic_numbering() function. However, in contrast to the fe_q_hierarchical numbering defined above, the lexicographic numbering originates from the tensor products of consecutive numbered dofs (like for LagrangeEquidistant).

It is assumed that the size of the output argument already matches the correct size, which is equal to the number of degrees of freedom in the finite element.

Definition at line 1946 of file fe_q_hierarchical.cc.

◆ face_fe_q_hierarchical_to_hierarchic_numbering()

template<int dim>
std::vector< unsigned int > FE_Q_Hierarchical< dim >::face_fe_q_hierarchical_to_hierarchic_numbering ( const unsigned int  degree)
staticprivate

This is an analogon to the previous function, but working on faces.

Definition at line 2120 of file fe_q_hierarchical.cc.

◆ build_dofs_cell()

template<int dim>
void FE_Q_Hierarchical< dim >::build_dofs_cell ( std::vector< FullMatrix< double >> &  dofs_cell,
std::vector< FullMatrix< double >> &  dofs_subcell 
) const
private

Initialize two auxiliary fields that will be used in setting up the various matrices in the constructor.

Definition at line 389 of file fe_q_hierarchical.cc.

◆ initialize_constraints()

template<int dim>
void FE_Q_Hierarchical< dim >::initialize_constraints ( const std::vector< FullMatrix< double >> &  dofs_subcell)
private

Initialize the hanging node constraints matrices. Called from the constructor.

Definition at line 503 of file fe_q_hierarchical.cc.

◆ initialize_embedding_and_restriction()

template<int dim>
void FE_Q_Hierarchical< dim >::initialize_embedding_and_restriction ( const std::vector< FullMatrix< double >> &  dofs_cell,
const std::vector< FullMatrix< double >> &  dofs_subcell 
)
private

Initialize the embedding matrices. Called from the constructor.

Definition at line 673 of file fe_q_hierarchical.cc.

◆ initialize_generalized_support_points()

template<int dim>
void FE_Q_Hierarchical< dim >::initialize_generalized_support_points ( )
private

Initialize the generalized_support_points field of the FiniteElement class. Called from the constructor.

Definition at line 808 of file fe_q_hierarchical.cc.

◆ initialize_generalized_face_support_points()

template<int dim>
void FE_Q_Hierarchical< dim >::initialize_generalized_face_support_points ( )
private

Initialize the generalized_face_support_points field of the FiniteElement class. Called from the constructor.

Definition at line 1884 of file fe_q_hierarchical.cc.

◆ has_support_on_face() [2/3]

template<>
bool FE_Q_Hierarchical< 1 >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const
virtual

This function returns true, if the shape function shape_index has non-zero function values somewhere on the face face_index. The function is typically used to determine whether some matrix elements resulting from face integrals can be assumed to be zero and may therefore be omitted from integration.

A default implementation is provided in this base class which always returns true. This is the safe way to go.

Reimplemented from FiniteElement< dim, dim >.

◆ has_support_on_face() [3/3]

template<>
bool FE_Q_Hierarchical< 1 >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const
virtual

This function returns true, if the shape function shape_index has non-zero function values somewhere on the face face_index. The function is typically used to determine whether some matrix elements resulting from face integrals can be assumed to be zero and may therefore be omitted from integration.

A default implementation is provided in this base class which always returns true. This is the safe way to go.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 2143 of file fe_q_hierarchical.cc.

Friends And Related Function Documentation

◆ FE_Q_Hierarchical

template<int dim>
template<int dim1>
friend class FE_Q_Hierarchical
friend

Allow access from other dimensions. We need this since we want to call the functions get_dpo_vector and lexicographic_to_hierarchic_numbering for the faces of the finite element of dimension dim+1.

Definition at line 795 of file fe_q_hierarchical.h.

Member Data Documentation

◆ face_renumber

template<int dim>
const std::vector<unsigned int> FE_Q_Hierarchical< dim >::face_renumber
private

Mapping from lexicographic to shape function numbering on first face.

Definition at line 786 of file fe_q_hierarchical.h.


The documentation for this class was generated from the following files: