Reference documentation for deal.II version 9.1.0-pre
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
FE_NedelecSZ< dim, spacedim > Class Template Reference

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

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

Classes

class  InternalData
 

Public Member Functions

 FE_NedelecSZ (const unsigned int degree)
 
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
virtual std::string get_name () const override
 
virtual std::unique_ptr< FiniteElement< dim, dim > > clone () 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
 
UpdateFlags update_once (const UpdateFlags flags) const
 
UpdateFlags update_each (const UpdateFlags flags) const
 
- 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 std::size_t memory_consumption () 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 bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
 
virtual const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
virtual const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
bool prolongation_is_implemented () const
 
bool isotropic_prolongation_is_implemented () const
 
bool restriction_is_implemented () const
 
bool isotropic_restriction_is_implemented () const
 
bool restriction_is_additive (const unsigned int index) const
 
const FullMatrix< double > & constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
bool constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
virtual bool hp_constraints_are_implemented () const
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual FiniteElementDomination::Domination compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const
 
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
 
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
 
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 (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from FiniteElementData< dim >
 FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
 
unsigned int n_dofs_per_vertex () const
 
unsigned int n_dofs_per_line () const
 
unsigned int n_dofs_per_quad () const
 
unsigned int n_dofs_per_hex () const
 
unsigned int n_dofs_per_face () const
 
unsigned int n_dofs_per_cell () const
 
template<int structdim>
unsigned int n_dofs_per_object () const
 
unsigned int n_components () const
 
unsigned int n_blocks () const
 
const BlockIndicesblock_indices () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Protected Member Functions

virtual std::unique_ptr< typename::FiniteElement< dim, spacedim >::InternalDataBase > get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) 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 &fedata,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &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 &fedata,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &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 &fedata,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &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< 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
 

Protected Attributes

MappingType mapping_type
 
- 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
 

Private Member Functions

void create_polynomials (const unsigned int degree)
 
unsigned int compute_num_dofs (const unsigned int degree) const
 
void fill_edge_values (const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
 
void fill_face_values (const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
 

Static Private Member Functions

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

Private Attributes

std::vector< Polynomials::Polynomial< double > > IntegratedLegendrePolynomials
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
- Static Public Member Functions inherited from FiniteElement< dim, dim >
static::ExceptionBase & ExcShapeFunctionNotPrimitive (int arg1)
 
static::ExceptionBase & ExcFENotPrimitive ()
 
static::ExceptionBase & ExcUnitShapeValuesDoNotExist ()
 
static::ExceptionBase & ExcFEHasNoSupportPoints ()
 
static::ExceptionBase & ExcEmbeddingVoid ()
 
static::ExceptionBase & ExcProjectionVoid ()
 
static::ExceptionBase & ExcWrongInterfaceMatrixSize (int arg1, int arg2)
 
static::ExceptionBase & ExcInterpolationNotImplemented ()
 
- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (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
 
- Static Protected Member Functions inherited from FiniteElement< dim, dim >
static std::vector< unsigned int > compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 

Detailed Description

template<int dim, int spacedim = dim>
class FE_NedelecSZ< dim, spacedim >

This class represents an implementation of the Hcurl-conforming Nédélec element described in the PhD thesis of S. Zaglmayr, High Order Finite Element Methods for Electromagnetic Field Computation, Johannes Kepler Universität Linz, 2006.

This element overcomes the sign conflict issues present in traditional Nédélec elements which arise from the edge and face parameterizations used in the basis functions. Therefore, this element should provide consistent results for general quadrilateral and hexahedral elements.

The way this element addresses the sign conflict problem is to assign local edges and faces a globally defined orientation for all cells. The local edge orientation is always chosen such that the first vertex defining the edge is always chosen to be that which has the highest global vertex numbering, with the second being that which has the lowest global vertex numbering.

Similarly, the face orientation is always chosen such that the first vertex is chosen to be that with the highest global vertex numbering of the four vertices making up the face. The third vertex is then chosen to be that which is geometrically opposite the first vertex, and the second and fourth vertices are decided such that the second has a higher global vertex numbering than the fourth.

Note that this element does not support non-conforming meshes at this time.

Further details on this element, including some benchmarking, can be in the paper R. Kynch, P. Ledger: Resolving the sign conflict problem for hp-hexahedral Nédélec elements with application to eddy current problems, Computers & Structures 181, 41-54, 2017 (see https://doi.org/10.1016/j.compstruc.2016.05.021).

Author
Ross Kynch
Date
2016, 2017, 2018

Definition at line 71 of file fe_nedelec_sz.h.

Constructor & Destructor Documentation

template<int dim, int spacedim>
FE_NedelecSZ< dim, spacedim >::FE_NedelecSZ ( const unsigned int  degree)

Constructor for an element of given degree.

Definition at line 24 of file fe_nedelec_sz.cc.

Member Function Documentation

template<int dim, int spacedim>
UpdateFlags FE_NedelecSZ< dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
overridevirtual

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

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

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

See also
UpdateFlags

Implements FiniteElement< dim, dim >.

Definition at line 2157 of file fe_nedelec_sz.cc.

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

Return a string that uniquely identifies a finite element. The general convention is that this is the class name, followed by the dimension in angle brackets, and the polynomial degree and whatever else is necessary in parentheses. For example, FE_Q<2>(3) is the value returned for a cubic element in 2d.

Systems of elements have their own naming convention, see the FESystem class.

Implements FiniteElement< dim, dim >.

Definition at line 2202 of file fe_nedelec_sz.cc.

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

Definition at line 2219 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
double FE_NedelecSZ< dim, spacedim >::shape_value ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

This element is vector-valued so this function will throw an exception.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 52 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
double FE_NedelecSZ< dim, spacedim >::shape_value_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Not implemented.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 64 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
Tensor< 1, dim > FE_NedelecSZ< dim, spacedim >::shape_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

This element is vector-valued so this function will throw an exception.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 78 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
Tensor< 1, dim > FE_NedelecSZ< dim, spacedim >::shape_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Not implemented.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 90 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
Tensor< 2, dim > FE_NedelecSZ< dim, spacedim >::shape_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

This element is vector-valued so this function will throw an exception.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 103 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
Tensor< 2, dim > FE_NedelecSZ< dim, spacedim >::shape_grad_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const
overridevirtual

Not implemented.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 115 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
UpdateFlags FE_NedelecSZ< dim, spacedim >::update_once ( const UpdateFlags  flags) const

Given flags, determines the values which must be computed only for the reference cell. Make sure, that mapping_type is set by the derived class, such that this function can operate correctly.

Definition at line 2165 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
UpdateFlags FE_NedelecSZ< dim, spacedim >::update_each ( const UpdateFlags  flags) const

Given flags, determines the values which must be computed in each cell cell. Make sure, that mapping_type is set by the derived class, such that this function can operate correctly.

Definition at line 2178 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
std::unique_ptr< typename::FiniteElement< dim, spacedim >::InternalDataBase > FE_NedelecSZ< dim, spacedim >::get_data ( const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const Quadrature< dim > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtual

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

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

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

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

Implements FiniteElement< dim, dim >.

Definition at line 127 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::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 &  fedata,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &  data 
) const
overrideprotectedvirtual

Compute information about the shape functions on the cell denoted by the first argument. Note that this function must recompute the cell-dependent degrees of freedom, and so is not thread-safe at this time.

Definition at line 1901 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::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 &  fedata,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &  data 
) const
overrideprotectedvirtual

Compute information about the shape functions on the cell and face denoted by the first two arguments. Note that this function must recompute the cell-dependent degrees of freedom, and so is not thread-safe at this time.

Definition at line 2013 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::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 &  fedata,
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &  data 
) const
overrideprotectedvirtual

Not implemented.

Definition at line 2139 of file fe_nedelec_sz.cc.

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

Internal function to return a vector of "dofs per object" where the components of the returned vector refer to: 0 = vertex 1 = edge 2 = face (which is a cell in 2D) 3 = cell

Definition at line 2226 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::create_polynomials ( const unsigned int  degree)
private

Internal function to populate the internal array of integrated Legendre polynomials.

Definition at line 2269 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
unsigned int FE_NedelecSZ< dim, spacedim >::compute_num_dofs ( const unsigned int  degree) const
private

Returns the number of DoFs in the basis set.

Definition at line 2247 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::fill_edge_values ( const typename Triangulation< dim, dim >::cell_iterator &  cell,
const Quadrature< dim > &  quadrature,
const InternalData fedata 
) const
private

Populates cell-dependent edge-based shape functions on the given InternalData object.

Definition at line 1102 of file fe_nedelec_sz.cc.

template<int dim, int spacedim>
void FE_NedelecSZ< dim, spacedim >::fill_face_values ( const typename Triangulation< dim, dim >::cell_iterator &  cell,
const Quadrature< dim > &  quadrature,
const InternalData fedata 
) const
private

Populates the cell-dependent face-based shape functions on the given InternalData object.

Definition at line 1527 of file fe_nedelec_sz.cc.

Member Data Documentation

template<int dim, int spacedim = dim>
MappingType FE_NedelecSZ< dim, spacedim >::mapping_type
protected

The mapping type to be used to map shape functions from the reference cell to the mesh cell.

Definition at line 157 of file fe_nedelec_sz.h.

template<int dim, int spacedim = dim>
std::vector<Polynomials::Polynomial<double> > FE_NedelecSZ< dim, spacedim >::IntegratedLegendrePolynomials
private

Internal storage for all required integrated Legendre polynomials.

Definition at line 417 of file fe_nedelec_sz.h.


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