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

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

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

Public Member Functions

 FE_TraceQ (unsigned int p)
 
virtual std::string get_name () const override
 
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone () const override
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const override
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const override
 
virtual bool hp_constraints_are_implemented () const override
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
 
- Public Member Functions inherited from FE_PolyFace< TensorProductPolynomials< dim - 1 >, dim, spacedim >
 FE_PolyFace (const TensorProductPolynomials< dim - 1 > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
 
unsigned int get_degree () const
 
- Public Member Functions inherited from FiniteElement< dim, spacedim >
 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 std::size_t memory_consumption () const
 
virtual double shape_value (const unsigned int i, const Point< dim > &p) const
 
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 3, dim > shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 3, dim > shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 4, dim > shape_4th_derivative (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 4, dim > shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) 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 void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, 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
 
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
 
- 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
 

Static Private Member Functions

static std::vector< unsigned int > get_dpo_vector (const unsigned int deg)
 

Private Attributes

FE_Q< dim, spacedim > fe_q
 

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, spacedim >
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, spacedim >
static const unsigned int space_dimension = spacedim
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Protected Member Functions inherited from FiniteElement< dim, spacedim >
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 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, spacedim >
static std::vector< unsigned int > compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FE_PolyFace< TensorProductPolynomials< dim - 1 >, dim, spacedim >
TensorProductPolynomials< dim - 1 > poly_space
 
- Protected Attributes inherited from FiniteElement< dim, spacedim >
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, int spacedim = dim>
class FE_TraceQ< dim, spacedim >

A finite element, which is the trace of FE_Q elements, that is a tensor product of polynomials on the faces, undefined in the interior of the cells and continuous. The basis functions on the faces are formed by a tensor product of 1D Lagrange polynomials with equidistant points up to degree 2 and Gauss-Lobatto points starting from degree 3.

This finite element is the trace space of FE_Q on the faces.

Note
Since these are only finite elements on faces, only FEFaceValues and FESubfaceValues will be able to extract reasonable values from any face polynomial. In order to make the use of FESystem simpler, FEValues objects will not fail using this finite element space, but all shape function values extracted will equal to zero.

Definition at line 46 of file fe_trace.h.

Constructor & Destructor Documentation

◆ FE_TraceQ()

template<int dim, int spacedim>
FE_TraceQ< dim, spacedim >::FE_TraceQ ( unsigned int  p)

Constructor for tensor product polynomials of degree p. The shape functions created using this constructor correspond to Legendre polynomials in each coordinate direction.

Definition at line 38 of file fe_trace.cc.

Member Function Documentation

◆ get_name()

template<int dim, int spacedim>
std::string FE_TraceQ< dim, spacedim >::get_name ( ) const
overridevirtual

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

Implements FiniteElement< dim, spacedim >.

Definition at line 84 of file fe_trace.cc.

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< FiniteElement< dim, spacedim > > FE_TraceQ< dim, spacedim >::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, spacedim >.

Definition at line 75 of file fe_trace.cc.

◆ convert_generalized_support_point_values_to_dof_values()

template<int dim, int spacedim>
void FE_TraceQ< dim, spacedim >::convert_generalized_support_point_values_to_dof_values ( const std::vector< Vector< double >> &  support_point_values,
std::vector< double > &  nodal_values 
) const
overridevirtual

Implementation of the corresponding function in the FiniteElement class. Since the current element is interpolatory, the nodal values are exactly the support point values. Furthermore, since the current element is scalar, the support point values need to be vectors of length 1.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 134 of file fe_trace.cc.

◆ has_support_on_face()

template<int dim, int spacedim>
bool FE_TraceQ< dim, spacedim >::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, spacedim >.

Definition at line 101 of file fe_trace.cc.

◆ get_constant_modes()

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

Return a list of constant modes of the element. For this element, it simply returns one row with all entries set to true.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 122 of file fe_trace.cc.

◆ hp_constraints_are_implemented()

template<int dim, int spacedim>
bool FE_TraceQ< dim, spacedim >::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".

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 172 of file fe_trace.cc.

◆ get_face_interpolation_matrix()

template<int dim, int spacedim>
void FE_TraceQ< dim, spacedim >::get_face_interpolation_matrix ( const FiniteElement< dim, spacedim > &  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. This element only provides interpolation matrices for elements of the same type and FE_Nothing. For all other elements, an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 219 of file fe_trace.cc.

◆ get_subface_interpolation_matrix()

template<int dim, int spacedim>
void FE_TraceQ< dim, spacedim >::get_subface_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
const unsigned int  subface,
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. This element only provides interpolation matrices for elements of the same type and FE_Nothing. For all other elements, an exception of type FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 232 of file fe_trace.cc.

◆ compare_for_domination()

template<int dim, int spacedim>
FiniteElementDomination::Domination FE_TraceQ< dim, spacedim >::compare_for_domination ( const FiniteElement< dim, spacedim > &  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.

Reimplemented from FiniteElement< dim, spacedim >.

Definition at line 180 of file fe_trace.cc.

◆ get_dpo_vector()

template<int dim, int spacedim>
std::vector< unsigned int > FE_TraceQ< dim, spacedim >::get_dpo_vector ( const unsigned int  deg)
staticprivate

Return vector with dofs per vertex, line, quad, hex.

Definition at line 154 of file fe_trace.cc.

Member Data Documentation

◆ fe_q

template<int dim, int spacedim = dim>
FE_Q<dim, spacedim> FE_TraceQ< dim, spacedim >::fe_q
private

Store a copy of FE_Q for delegating the hp-constraints functionality.

Definition at line 138 of file fe_trace.h.


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