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\}}\)
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
FE_Nedelec< dim > Class Template Reference

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

Inheritance diagram for FE_Nedelec< dim >:
[legend]

Public Member Functions

 FE_Nedelec (const unsigned int order)
 
virtual std::string get_name () const override
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override
 
virtual bool hp_constraints_are_implemented () const override
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
 
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 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 const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) 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 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 std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const override
 
virtual std::size_t memory_consumption () const override
 
virtual std::unique_ptr< FiniteElement< dim, dim > > clone () const override
 
void get_subface_interpolation_matrix (const FiniteElement< 1, 1 > &, const unsigned int, FullMatrix< double > &) const
 
- Public Member Functions inherited from FE_PolyTensor< dim >
 FE_PolyTensor (const TensorPolynomialsBase< dim > &polynomials, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
 FE_PolyTensor (const FE_PolyTensor &fe)
 
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
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
 
- 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 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 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
 
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 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
 

Private Member Functions

void initialize_support_points (const unsigned int order)
 
void initialize_restriction ()
 
void initialize_support_points (const unsigned int)
 
void initialize_support_points (const unsigned int order)
 
void initialize_support_points (const unsigned int order)
 
void initialize_restriction ()
 

Static Private Member Functions

static std::vector< unsigned intget_dpo_vector (const unsigned int degree, bool dg=false)
 

Private Attributes

Table< 2, doubleboundary_weights
 
Threads::Mutex mutex
 

Friends

template<int dim1>
class FE_Nedelec
 

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_PolyTensor< dim >
bool single_mapping_kind () const
 
MappingKind get_mapping_kind (const unsigned int i) const
 
virtual std::unique_ptr< typename FiniteElement< dim, dim >::InternalDataBaseget_data (const UpdateFlags update_flags, const Mapping< dim, dim > &, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &) const override
 
virtual void fill_fe_values (const typename Triangulation< dim, dim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &output_data) const override
 
virtual void fill_fe_face_values (const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &output_data) const override
 
virtual void fill_fe_subface_values (const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &output_data) const override
 
- 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< InternalDataBaseget_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< 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
 
- Static Protected Member Functions inherited from FiniteElement< dim, dim >
static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FE_PolyTensor< dim >
std::vector< MappingKindmapping_kind
 
const std::unique_ptr< const TensorPolynomialsBase< dim > > poly_space
 
FullMatrix< doubleinverse_node_matrix
 
std::mutex cache_mutex
 
Point< dim > cached_point
 
std::vector< Tensor< 1, dim > > cached_values
 
std::vector< Tensor< 2, dim > > cached_grads
 
std::vector< Tensor< 3, dim > > cached_grad_grads
 
- Protected Attributes inherited from FiniteElement< dim, dim >
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
 

Detailed Description

template<int dim>
class FE_Nedelec< dim >

Warning
Several aspects of the implementation are experimental. For the moment, it is safe to use the element on globally refined meshes with consistent orientation of faces. See the todo entries below for more detailed caveats.

Implementation of Nédélec elements. The Nédélec space is designed to solve problems in which the solution only lives in the space \(H^\text{curl}=\{ {\mathbf u} \in L_2: \text{curl}\, {\mathbf u} \in L_2\}\), rather than in the more commonly used space \(H^1=\{ u \in L_2: \nabla u \in L_2\}\). In other words, the solution must be a vector field whose curl is square integrable, but for which the gradient may not be square integrable. The typical application for this space (and these elements) is to the Maxwell equations and corresponding simplifications, such as the reduced version of the Maxwell equation that only involves the electric field \(\mathbf E\) which has to satisfy the equation \(\text{curl}\, \text{curl}\, {\mathbf E} = 0\) in the time independent case when no currents are present, or the equation \(\text{curl}\,\text{curl}\,{\mathbf A} = 4\pi{\mathbf j}\) that the magnetic vector potential \(\mathbf A\) has to satisfy in the time independent case.

The defining characteristic of functions in \(H^\text{curl}\) is that they are in general discontinuous – but that if you draw a line in 2d (or a surface in 3d), then the tangential component(s) of the vector field must be continuous across the line (or surface) even though the normal component may not be. As a consequence, the Nédélec element is constructed in such a way that (i) it is vector-valued, (ii) the shape functions are discontinuous, but (iii) the tangential component(s) of the vector field represented by each shape function are continuous across the faces of cells.

Other properties of the Nédélec element are that (i) it is not a primitive element; (ii) the shape functions are defined so that certain integrals over the faces are either zero or one, rather than the common case of certain point values being either zero or one.

We follow the commonly used – though confusing – definition of the "degree" of Nédélec elements. Specifically, the "degree" of the element denotes the polynomial degree of the largest complete polynomial subspace contained in the finite element space, even if the space may contain shape functions of higher polynomial degree. The lowest order element is consequently FE_Nedelec(0), i.e., the Raviart-Thomas element "of degree zero", even though the functions of this space are in general polynomials of degree one in each variable. This choice of "degree" implies that the approximation order of the function itself is degree+1, as with usual polynomial spaces. The numbering so chosen implies the sequence

\[ Q_{k+1} \stackrel{\text{grad}}{\rightarrow} \text{Nedelec}_k \stackrel{\text{curl}}{\rightarrow} \text{RaviartThomas}_k \stackrel{\text{div}}{\rightarrow} DGQ_{k} \]

Note that this follows the convention of Brezzi and Raviart, though not the one used in the original paper by Nédélec.

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

Todo:
Even if this element is implemented for two and three space dimensions, the definition of the node values relies on consistently oriented faces in 3D. Therefore, care should be taken on complicated meshes.

Interpolation

The interpolation operators associated with the Nédélec element are constructed such that interpolation and computing the curl are commuting operations on rectangular mesh cells. We require this from interpolating arbitrary functions as well as the restriction matrices.

Node values

The node values for an element of degree k on the reference cell are:

  1. On edges: the moments of the tangential component with respect to polynomials of degree k.
  2. On faces: the moments of the tangential components with respect to dim-1 dimensional FE_Nedelec polynomials of degree k-1.
  3. In cells: the moments with respect to gradients of polynomials in FE_Q of degree k.

Generalized support points

The node values above rely on integrals, which will be computed by quadrature rules themselves. The generalized support points are a set of points such that this quadrature can be performed with sufficient accuracy. The points needed are those of QGaussk+1 on each edge and QGaussk+2 on each face and in the interior of the cell (or none for N1).

Author
Markus Bürg
Date
2009, 2010, 2011

Definition at line 147 of file fe_nedelec.h.

Constructor & Destructor Documentation

◆ FE_Nedelec()

template<int dim>
FE_Nedelec< dim >::FE_Nedelec ( const unsigned int  order)

Constructor for the Nedelec element of given order. The maximal polynomial degree of the shape functions is order+1 (in each variable; the total polynomial degree may be higher). If order = 0, the element is linear and has degrees of freedom only on the edges. If order >=1 the element has degrees of freedom on the edges, faces and volume. For example the 3D version of FE_Nedelec has 12 degrees of freedom for order = 0 and 54 for degree = 1. It is important to have enough quadrature points in order to perform the quadrature with sufficient accuracy. For example QGauss<dim>(order + 2) can be used for the quadrature formula, where order is the order of FE_Nedelec.

Definition at line 71 of file fe_nedelec.cc.

Member Function Documentation

◆ get_name()

template<int dim>
std::string FE_Nedelec< dim >::get_name
overridevirtual

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

Definition at line 207 of file fe_nedelec.cc.

◆ has_support_on_face()

template<int dim>
bool FE_Nedelec< 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 2035 of file fe_nedelec.cc.

◆ hp_constraints_are_implemented()

template<int dim>
bool FE_Nedelec< dim >::hp_constraints_are_implemented
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_Nedelec 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.

Definition at line 2310 of file fe_nedelec.cc.

◆ compare_for_domination()

template<int dim>
FiniteElementDomination::Domination FE_Nedelec< 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 2274 of file fe_nedelec.cc.

◆ hp_vertex_dof_identities()

template<int dim>
std::vector< std::pair< unsigned int, unsigned int > > FE_Nedelec< 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 2317 of file fe_nedelec.cc.

◆ hp_line_dof_identities()

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

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.

Definition at line 2326 of file fe_nedelec.cc.

◆ hp_quad_dof_identities()

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

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.

Definition at line 2368 of file fe_nedelec.cc.

◆ get_face_interpolation_matrix()

template<int dim>
void FE_Nedelec< 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 2422 of file fe_nedelec.cc.

◆ get_subface_interpolation_matrix() [1/2]

template<int dim>
void FE_Nedelec< 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 2524 of file fe_nedelec.cc.

◆ get_restriction_matrix()

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

Projection from a fine grid space onto a coarse grid space. If this projection operator is associated with a matrix P, then the restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation or the sum of the cell matrices P_i, depending on the restriction_is_additive_flags. This distinguishes interpolation (concatenation) and projection with respect to scalar products (summation).

Row and column indices are related to coarse grid and fine grid spaces, respectively, consistent with the definition of the associated operator.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 3010 of file fe_nedelec.cc.

◆ get_prolongation_matrix()

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

Embedding matrix between grids.

The identity operator from a coarse grid space into a fine grid space is associated with a matrix P. The restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation, not the sum of the cell matrices P_i. That is, if the same non-zero entry j,k exists in two different child matrices P_i, the value should be the same in both matrices and it is copied into the matrix P only once.

Row and column indices are related to fine grid and coarse grid spaces, respectively, consistent with the definition of the associated operator.

These matrices are used by routines assembling the prolongation matrix for multi-level methods. Upon assembling the transfer matrix between cells using this matrix array, zero elements in the prolongation matrix are discarded and will not fill up the transfer matrix.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 2955 of file fe_nedelec.cc.

◆ convert_generalized_support_point_values_to_dof_values()

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

Given the values of a function \(f(\mathbf x)\) at the (generalized) support points of the reference cell, this function then computes what the nodal values of the element are, i.e., \(\Psi_i[f]\), where \(\Psi_i\) are the node functionals of the element (see also Node values or node functionals). The values \(\Psi_i[f]\) are then the expansion coefficients for the shape functions of the finite element function that interpolates the given function \(f(x)\), i.e., \( f_h(\mathbf x) = \sum_i \Psi_i[f] \varphi_i(\mathbf x) \) is the finite element interpolant of \(f\) with the current element. The operation described here is used, for example, in the FETools::compute_node_matrix() function.

In more detail, let us assume that the generalized support points (see this glossary entry ) of the current element are \(\hat{\mathbf x}_i\) and that the node functionals associated with the current element are \(\Psi_i[\cdot]\). Then, the fact that the element is based on generalized support points, implies that if we apply \(\Psi_i\) to a (possibly vector-valued) finite element function \(\varphi\), the result must have the form \(\Psi_i[\varphi] = f_i(\varphi(\hat{\mathbf x}_i))\) – in other words, the value of the node functional \(\Psi_i\) applied to \(\varphi\) only depends on the values of \(\varphi\) at \(\hat{\mathbf x}_i\) and not on values anywhere else, or integrals of \(\varphi\), or any other kind of information.

The exact form of \(f_i\) depends on the element. For example, for scalar Lagrange elements, we have that in fact \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)\). If you combine multiple scalar Lagrange elements via an FESystem object, then \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)_{c(i)}\) where \(c(i)\) is the result of the FiniteElement::system_to_component_index() function's return value's first component. In these two cases, \(f_i\) is therefore simply the identity (in the scalar case) or a function that selects a particular vector component of its argument. On the other hand, for Raviart-Thomas elements, one would have that \(f_i(\mathbf y) = \mathbf y \cdot \mathbf n_i\) where \(\mathbf n_i\) is the normal vector of the face at which the shape function is defined.

Given all of this, what this function does is the following: If you input a list of values of a function \(\varphi\) at all generalized support points (where each value is in fact a vector of values with as many components as the element has), then this function returns a vector of values obtained by applying the node functionals to these values. In other words, if you pass in \(\{\varphi(\hat{\mathbf x}_i)\}_{i=0}^{N-1}\) then you will get out a vector \(\{\Psi[\varphi]\}_{i=0}^{N-1}\) where \(N\) equals dofs_per_cell.

Parameters
[in]support_point_valuesAn array of size dofs_per_cell (which equals the number of points the get_generalized_support_points() function will return) where each element is a vector with as many entries as the element has vector components. This array should contain the values of a function at the generalized support points of the current element.
[out]nodal_valuesAn array of size dofs_per_cell that contains the node functionals of the element applied to the given function.
Note
It is safe to call this function for (transformed) values on the real cell only for elements with trivial MappingKind. For all other elements (for example for H(curl), or H(div) conforming elements) vector values have to be transformed to the reference cell first.
Given what the function is supposed to do, the function clearly can only work for elements that actually implement (generalized) support points. Elements that do not have generalized support points – e.g., elements whose nodal functionals evaluate integrals or moments of functions (such as FE_Q_Hierarchical) – can in general not make sense of the operation that is required for this function. They consequently may not implement it.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 3073 of file fe_nedelec.cc.

◆ get_constant_modes()

template<int dim>
std::pair< Table< 2, bool >, std::vector< unsigned int > > FE_Nedelec< dim >::get_constant_modes
overridevirtual

Return a list of constant modes of the element.

Definition at line 4019 of file fe_nedelec.cc.

◆ memory_consumption()

template<int dim>
std::size_t FE_Nedelec< dim >::memory_consumption
overridevirtual

Definition at line 4035 of file fe_nedelec.cc.

◆ clone()

template<int dim>
std::unique_ptr< FiniteElement< dim, dim > > FE_Nedelec< dim >::clone
overridevirtual

Definition at line 225 of file fe_nedelec.cc.

◆ get_dpo_vector()

template<int dim>
std::vector< unsigned int > FE_Nedelec< dim >::get_dpo_vector ( const unsigned int  degree,
bool  dg = false 
)
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.

If the optional argument dg is true, the vector returned will have all degrees of freedom assigned to the cell, none on the faces and edges.

Definition at line 1999 of file fe_nedelec.cc.

◆ initialize_support_points() [1/4]

template<int dim>
void FE_Nedelec< dim >::initialize_support_points ( const unsigned int  order)
private

Initialize the generalized_support_points field of the FiniteElement class and fill the tables with interpolation weights (boundary_weights and interior_weights). Called from the constructor.

◆ initialize_restriction() [1/2]

template<int dim>
void FE_Nedelec< dim >::initialize_restriction
private

Initialize the interpolation from functions on refined mesh cells onto the father cell. According to the philosophy of the Nédélec element, this restriction operator preserves the curl of a function weakly.

Definition at line 515 of file fe_nedelec.cc.

◆ initialize_support_points() [2/4]

void FE_Nedelec< 1 >::initialize_support_points ( const unsigned int  )
private

Definition at line 244 of file fe_nedelec.cc.

◆ initialize_support_points() [3/4]

void FE_Nedelec< 2 >::initialize_support_points ( const unsigned int  order)
private

Definition at line 253 of file fe_nedelec.cc.

◆ initialize_support_points() [4/4]

void FE_Nedelec< 3 >::initialize_support_points ( const unsigned int  order)
private

Definition at line 335 of file fe_nedelec.cc.

◆ initialize_restriction() [2/2]

void FE_Nedelec< 1 >::initialize_restriction ( )
private

Definition at line 502 of file fe_nedelec.cc.

◆ get_subface_interpolation_matrix() [2/2]

void FE_Nedelec< 1 >::get_subface_interpolation_matrix ( const FiniteElement< 1, 1 > &  ,
const unsigned int  ,
FullMatrix< double > &   
) const

Definition at line 2497 of file fe_nedelec.cc.

Friends And Related Function Documentation

◆ FE_Nedelec

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

Definition at line 373 of file fe_nedelec.h.

Member Data Documentation

◆ boundary_weights

template<int dim>
Table<2, double> FE_Nedelec< dim >::boundary_weights
private

These are the factors multiplied to a function in the generalized_face_support_points when computing the integration.

See the glossary entry on generalized support points for more information.

Definition at line 364 of file fe_nedelec.h.

◆ mutex

template<int dim>
Threads::Mutex FE_Nedelec< dim >::mutex
mutableprivate

Mutex for protecting initialization of restriction and embedding matrix.

Definition at line 369 of file fe_nedelec.h.


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