Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/fe_nedelec.h>
Public Member Functions | |
FE_Nedelec (const unsigned int order) | |
virtual std::string | get_name () const |
virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const |
virtual bool | hp_constraints_are_implemented () const |
virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_vertex_dof_identities (const FiniteElement< dim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_line_dof_identities (const FiniteElement< dim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_quad_dof_identities (const FiniteElement< dim > &fe_other) const |
virtual void | get_face_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const |
virtual void | get_subface_interpolation_matrix (const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix) 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 |
virtual void | convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const |
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () const |
virtual std::size_t | memory_consumption () const |
virtual std::unique_ptr< FiniteElement< dim, dim > > | clone () const |
Public Member Functions inherited from FE_PolyTensor< PolynomialsNedelec< dim >, dim > | |
FE_PolyTensor (const unsigned int degree, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
virtual double | shape_value (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 1, dim > | shape_grad (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 2, dim > | shape_grad_grad (const unsigned int i, const Point< dim > &p) 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 ()=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 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 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 ComponentMask & | get_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_index > | system_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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () 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 BlockIndices & | block_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 () |
Static Private Member Functions | |
static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree, bool dg=false) |
Private Attributes | |
Table< 2, double > | boundary_weights |
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 ::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 |
Protected Member Functions inherited from FiniteElement< dim, dim > | |
void | reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false) |
TableIndices< 2 > | interface_constraints_size () const |
virtual std::unique_ptr< InternalDataBase > | get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
virtual std::unique_ptr< InternalDataBase > | get_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const |
virtual std::unique_ptr< InternalDataBase > | get_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const |
virtual void | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
Static Protected Member Functions inherited from FiniteElement< dim, dim > | |
static std::vector< unsigned int > | compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components) |
Protected Attributes inherited from FE_PolyTensor< PolynomialsNedelec< dim >, dim > | |
MappingType | mapping_type |
PolynomialsNedelec< dim > | poly_space |
FullMatrix< double > | inverse_node_matrix |
Threads::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< 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< ComponentMask > | nonzero_components |
const std::vector< unsigned int > | n_nonzero_components_table |
const bool | cached_primitivity |
Implementation of Nédélec elements, conforming with the space Hcurl. These elements generate vector fields with tangential components continuous between mesh cells.
We follow the convention that the degree of Nédélec elements denotes the polynomial degree of the largest complete polynomial subspace contained in the Nédélec space. This leads to the consistently numbered sequence of spaces
\[ Q_{k+1} \stackrel{\text{grad}}{\rightarrow} \text{Nedelec}_k \stackrel{\text{curl}}{\rightarrow} \text{RaviartThomas}_k \stackrel{\text{div}}{\rightarrow} DGQ_{k} \]
Consequently, approximation order of the Nédélec space equals the value degree given to the constructor. In this scheme, the lowest order element would be created by the call FE_Nedelec<dim>(0). 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
).
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.
The node values for an element of degree k on the reference cell are:
dim
-1 dimensional FE_Nedelec polynomials of degree k-1. 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).
Definition at line 111 of file fe_nedelec.h.
FE_Nedelec< dim >::FE_Nedelec | ( | const unsigned int | order | ) |
Constructor for the Nédélec 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).
Definition at line 67 of file fe_nedelec.cc.
|
virtual |
Return a string that uniquely identifies a finite element. This class returns FE_Nedelec<dim>(degree)
, with dim
and degree
replaced by appropriate values.
Implements FiniteElement< dim, dim >.
Definition at line 198 of file fe_nedelec.cc.
|
virtual |
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 2027 of file fe_nedelec.cc.
|
virtual |
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.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 2355 of file fe_nedelec.cc.
|
virtual |
Return whether this element dominates the one, which is given as argument.
Definition at line 2311 of file fe_nedelec.cc.
|
virtual |
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 2362 of file fe_nedelec.cc.
|
virtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
Definition at line 2372 of file fe_nedelec.cc.
|
virtual |
Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.
Definition at line 2413 of file fe_nedelec.cc.
|
virtual |
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 2468 of file fe_nedelec.cc.
|
virtual |
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 2573 of file fe_nedelec.cc.
|
virtual |
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 3042 of file fe_nedelec.cc.
|
virtual |
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 2991 of file fe_nedelec.cc.
|
virtual |
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
.
[in] | support_point_values | An 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_values | An array of size dofs_per_cell that contains the node functionals of the element applied to the given function. |
Reimplemented from FiniteElement< dim, dim >.
Definition at line 3100 of file fe_nedelec.cc.
|
virtual |
Return a list of constant modes of the element.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 4018 of file fe_nedelec.cc.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 4034 of file fe_nedelec.cc.
|
virtual |
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 216 of file fe_nedelec.cc.
|
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 1991 of file fe_nedelec.cc.
|
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.
|
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 525 of file fe_nedelec.cc.
Allow access from other dimensions.
Definition at line 322 of file fe_nedelec.h.
|
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 312 of file fe_nedelec.h.