KASKADE 7 development version
|
A shape function set of \( H(\text{div}) \) conforming functions. More...
#include <nedelecshapefunctions.hh>
A shape function set of \( H(\text{div}) \) conforming functions.
For each Lagrangian interpolation point \( \xi_i \) of each face, there is one shape function with non-vanishing normal component of size 1 exactly at this point. Those shape functions are associated to the corresponding face. The remaining basis functions that complete the polynomial space all have vanishing normal components and are associated to the cell.
Only the normal components at the face interpolation points are specified, the tangential components are well-defined but unspecified (implementation detail).
Definition at line 446 of file nedelecshapefunctions.hh.
Public Types | |
typedef HdivSimplexShapeFunction< ctype, dim, T > | value_type |
typedef double | Scalar |
Scalar field type. More... | |
typedef DynamicMatrix< Dune::FieldMatrix< double, 1, 1 > > | Matrix |
A matrix type mapping one-component coefficient vectors. More... | |
typedef std::vector< Dune::FieldVector< ctype, dimension > > | InterpolationNodes |
A container type for holding interpolation points in the reference elements. More... | |
typedef DynamicMatrix< Dune::FieldMatrix< double, comp, 1 > > | SfValueArray |
Public Member Functions | |
HdivSimplexShapeFunctionSet (int order) | |
Constructs a \( H(\text{div}) \) conforming shape function set of given polynomial order. More... | |
virtual value_type const & | operator[] (int i) const |
Random access to a shape function. More... | |
virtual Dune::GeometryType | type () const |
virtual void | interpolate (typename ShapeFunctionSet< ctype, dim, T, dim >::SfValueArray const &A, DynamicMatrix< Dune::FieldMatrix< T, 1, 1 > > &IA) const |
Left-multiplies the provided matrix with the interpolation matrix of the shape function set. More... | |
virtual int | size () const |
Number of shape functions in the set. More... | |
int | order () const |
Maximal polynomial order of shape functions. More... | |
auto | referenceElement () const |
void | initHierarchicalProjection (ShapeFunctionSet< ctype, dimension, double, comp > const *sfl) |
Initialize the hierarchical projection matrix based on the given lower order shape function set. More... | |
InterpolationNodes const & | interpolationNodes () const |
Interpolation points. More... | |
virtual void | interpolate (SfValueArray const &A, Matrix &IA) const=0 |
Left-multiplies the provided matrix with the interpolation matrix of the shape function set. More... | |
void | evaluate (InterpolationNodes const &iNodes, SfValueArray &phi) const |
Evaluate shape function set at a set of points. In notation of the LocalToGlobalMapperConcept, this gives the matrix \( \Phi \): the entry Phi[i][j] is the value of shape function j evaluated at iNodes[i]. More... | |
Matrix const & | hierarchicProjection () const |
Returns a square matrix that projects shape function coefficients to a subspace spanned by shape functions of lower order. This is intended to be used to implement embedded error estimators. The actual definition of this subspace depends on the shape function set specified in the call to initHierarchicalProjection(). More... | |
virtual void | removeShapeFunction (size_t index) |
Static Public Attributes | |
static int const | comps |
Number of components of the shape function values. More... | |
Protected Attributes | |
InterpolationNodes | iNodes |
Matrix | projection |
int | order_ |
int | size_ |
|
inherited |
A container type for holding interpolation points in the reference elements.
Definition at line 164 of file pshapefunctions.hh.
|
inherited |
A matrix type mapping one-component coefficient vectors.
Definition at line 119 of file pshapefunctions.hh.
|
inherited |
Scalar field type.
Definition at line 112 of file pshapefunctions.hh.
|
inherited |
A twodimensional array type for holding shape function values evaluated at a set of nodes.
Definition at line 168 of file pshapefunctions.hh.
typedef HdivSimplexShapeFunction<ctype,dim,T> Kaskade::HdivSimplexShapeFunctionSet< ctype, dimension, T >::value_type |
Definition at line 452 of file nedelecshapefunctions.hh.
Kaskade::HdivSimplexShapeFunctionSet< ctype, dimension, T >::HdivSimplexShapeFunctionSet | ( | int | order | ) |
Constructs a \( H(\text{div}) \) conforming shape function set of given polynomial order.
order | the polynomial order of the shape functions (at least 1) |
|
inlineinherited |
Evaluate shape function set at a set of points. In notation of the LocalToGlobalMapperConcept, this gives the matrix \( \Phi \): the entry Phi[i][j] is the value of shape function j evaluated at iNodes[i].
iNodes | the points at which the shape functions are to be evaluated. |
phi | the array that is filled with shape function values. The array will be resized if needed. |
Definition at line 253 of file pshapefunctions.hh.
|
inlineinherited |
Returns a square matrix that projects shape function coefficients to a subspace spanned by shape functions of lower order. This is intended to be used to implement embedded error estimators. The actual definition of this subspace depends on the shape function set specified in the call to initHierarchicalProjection().
Definition at line 276 of file pshapefunctions.hh.
|
inlineinherited |
Initialize the hierarchical projection matrix based on the given lower order shape function set.
Definition at line 174 of file pshapefunctions.hh.
|
pure virtualinherited |
Left-multiplies the provided matrix with the interpolation matrix of the shape function set.
Each column of A is interpreted as values of some function evaluated at this shape function set's interpolation points (see below). The columns of the output array IA then contain the shape function coefficients such that the corresponding linearcombination of shape functions "interpolates" that function in the interpolation points. What "interpolation" means is up to the actual implementation.
IA is automatically resized if needed.
Storage order: A[i][j] contains the value of function j at interpolation node i.
|
inlinevirtual |
Left-multiplies the provided matrix with the interpolation matrix of the shape function set.
Each column of A is interpreted as values of some function evaluated at this shape function set's interpolation points (edge midpoints). The columns of the output matrix IA then contain the shape function coefficients such that the corresponding linear combination of shape functions "interpolates" that function in the interpolation points. "Interpolation" means equal tangential components.
Definition at line 475 of file nedelecshapefunctions.hh.
|
inlineinherited |
Interpolation points.
Provides interpolation points such that the shape function coefficients can be computed from function values at interpolation nodes by multiplication by this matrix. The interpolation points are guaranteed to be inside the reference element associated to this shape function set.
Definition at line 218 of file pshapefunctions.hh.
|
inlinevirtual |
Random access to a shape function.
Implements Kaskade::ShapeFunctionSet< ctype, dimension, double, dimension >.
Definition at line 460 of file nedelecshapefunctions.hh.
|
inlineinherited |
Maximal polynomial order of shape functions.
Definition at line 149 of file pshapefunctions.hh.
|
inlineinherited |
Returns a reference to the reference element on which this shape function set is defined.
Definition at line 158 of file pshapefunctions.hh.
|
inlinevirtualinherited |
Definition at line 284 of file pshapefunctions.hh.
|
inlinevirtualinherited |
Number of shape functions in the set.
Definition at line 144 of file pshapefunctions.hh.
|
inlinevirtual |
Definition at line 462 of file nedelecshapefunctions.hh.
|
staticinherited |
Number of components of the shape function values.
Definition at line 132 of file pshapefunctions.hh.
|
protectedinherited |
Definition at line 287 of file pshapefunctions.hh.
|
protectedinherited |
Definition at line 289 of file pshapefunctions.hh.
|
protectedinherited |
Definition at line 288 of file pshapefunctions.hh.
|
protectedinherited |
Definition at line 290 of file pshapefunctions.hh.