KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | Protected Attributes | List of all members
Kaskade::ShapeFunctionSet< ctype, dimension, T, comp > Class Template Referenceabstract

A set of shape functions. More...

#include <pshapefunctions.hh>

Detailed Description

template<class ctype, int dimension, class T, int comp = 1>
class Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >

A set of shape functions.

Template Parameters
ctypescalar type for coordinates
dimensionspatial dimension
Tscalar type
compnumber of components

Definition at line 108 of file pshapefunctions.hh.

Public Types

typedef T Scalar
 Scalar field type. More...
 
typedef ShapeFunction< ctype, dimension, T, comp > value_type
 type of one shape function More...
 
typedef DynamicMatrix< Dune::FieldMatrix< T, 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< T, comp, 1 > > SfValueArray
 

Public Member Functions

 ShapeFunctionSet (Dune::GeometryType gt)
 Constructor. More...
 
 ShapeFunctionSet (ShapeFunctionSet const &other)=default
 
virtual ~ShapeFunctionSet ()
 
virtual value_type const & operator[] (int i) const =0
 Random access to a shape function. More...
 
virtual int size () const
 Number of shape functions in the set. More...
 
int order () const
 Maximal polynomial order of shape functions. More...
 
Dune::GeometryType type () const
 Type of geometry on which the shape functions are defined. More...
 
auto referenceElement () const
 
void initHierarchicalProjection (ShapeFunctionSet< ctype, dimension, T, 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 = comp
 Number of components of the shape function values. More...
 

Protected Attributes

InterpolationNodes iNodes
 
Matrix projection
 
int order_
 
int size_
 

Member Typedef Documentation

◆ InterpolationNodes

template<class ctype , int dimension, class T , int comp = 1>
typedef std::vector<Dune::FieldVector<ctype,dimension> > Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::InterpolationNodes

A container type for holding interpolation points in the reference elements.

Definition at line 164 of file pshapefunctions.hh.

◆ Matrix

template<class ctype , int dimension, class T , int comp = 1>
typedef DynamicMatrix<Dune::FieldMatrix<T,1,1> > Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::Matrix

A matrix type mapping one-component coefficient vectors.

Definition at line 119 of file pshapefunctions.hh.

◆ Scalar

template<class ctype , int dimension, class T , int comp = 1>
typedef T Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::Scalar

Scalar field type.

Definition at line 112 of file pshapefunctions.hh.

◆ SfValueArray

template<class ctype , int dimension, class T , int comp = 1>
typedef DynamicMatrix<Dune::FieldMatrix<T,comp,1> > Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::SfValueArray

A twodimensional array type for holding shape function values evaluated at a set of nodes.

Definition at line 168 of file pshapefunctions.hh.

◆ value_type

template<class ctype , int dimension, class T , int comp = 1>
typedef ShapeFunction<ctype,dimension,T,comp> Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::value_type

type of one shape function

Definition at line 114 of file pshapefunctions.hh.

Constructor & Destructor Documentation

◆ ShapeFunctionSet() [1/2]

template<class ctype , int dimension, class T , int comp = 1>
Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::ShapeFunctionSet ( Dune::GeometryType  gt)
inline

Constructor.

Definition at line 122 of file pshapefunctions.hh.

◆ ShapeFunctionSet() [2/2]

template<class ctype , int dimension, class T , int comp = 1>
Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::ShapeFunctionSet ( ShapeFunctionSet< ctype, dimension, T, comp > const &  other)
default

◆ ~ShapeFunctionSet()

template<class ctype , int dimension, class T , int comp = 1>
virtual Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::~ShapeFunctionSet ( )
inlinevirtual

Definition at line 135 of file pshapefunctions.hh.

Member Function Documentation

◆ evaluate()

template<class ctype , int dimension, class T , int comp = 1>
void Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::evaluate ( InterpolationNodes const &  iNodes,
SfValueArray phi 
) const
inline

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].

Parameters
iNodesthe points at which the shape functions are to be evaluated.
phithe array that is filled with shape function values. The array will be resized if needed.

Definition at line 253 of file pshapefunctions.hh.

Referenced by Kaskade::computeSimplexSfPermutation(), Kaskade::evaluateGlobalShapeFunctions(), and Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::initHierarchicalProjection().

◆ hierarchicProjection()

template<class ctype , int dimension, class T , int comp = 1>
Matrix const & Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::hierarchicProjection ( ) const
inline

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.

◆ initHierarchicalProjection()

template<class ctype , int dimension, class T , int comp = 1>
void Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::initHierarchicalProjection ( ShapeFunctionSet< ctype, dimension, T, comp > const *  sfl)
inline

◆ interpolate()

template<class ctype , int dimension, class T , int comp = 1>
virtual void Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::interpolate ( SfValueArray const &  A,
Matrix IA 
) const
pure virtual

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.

Referenced by Kaskade::approximateGlobalValues(), Kaskade::computeSimplexSfPermutation(), and Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::initHierarchicalProjection().

◆ interpolationNodes()

template<class ctype , int dimension, class T , int comp = 1>
InterpolationNodes const & Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::interpolationNodes ( ) const
inline

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.

Referenced by Kaskade::approximateGlobalValues(), Kaskade::computeSimplexSfPermutation(), and Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::initHierarchicalProjection().

◆ operator[]()

template<class ctype , int dimension, class T , int comp = 1>
virtual value_type const & Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::operator[] ( int  i) const
pure virtual

◆ order()

template<class ctype , int dimension, class T , int comp = 1>
int Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::order ( ) const
inline

Maximal polynomial order of shape functions.

Definition at line 149 of file pshapefunctions.hh.

◆ referenceElement()

template<class ctype , int dimension, class T , int comp = 1>
auto Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::referenceElement ( ) const
inline

Returns a reference to the reference element on which this shape function set is defined.

Definition at line 158 of file pshapefunctions.hh.

◆ removeShapeFunction()

template<class ctype , int dimension, class T , int comp = 1>
virtual void Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::removeShapeFunction ( size_t  index)
inlinevirtual

◆ size()

template<class ctype , int dimension, class T , int comp = 1>
virtual int Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::size ( ) const
inlinevirtual

◆ type()

template<class ctype , int dimension, class T , int comp = 1>
Dune::GeometryType Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::type ( ) const
inline

Type of geometry on which the shape functions are defined.

Definition at line 152 of file pshapefunctions.hh.

Member Data Documentation

◆ comps

template<class ctype , int dimension, class T , int comp = 1>
int const Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::comps = comp
static

Number of components of the shape function values.

Definition at line 132 of file pshapefunctions.hh.

◆ iNodes

template<class ctype , int dimension, class T , int comp = 1>
InterpolationNodes Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::iNodes
protected

◆ order_

template<class ctype , int dimension, class T , int comp = 1>
int Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::order_
protected

◆ projection

template<class ctype , int dimension, class T , int comp = 1>
Matrix Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::projection
protected

◆ size_

template<class ctype , int dimension, class T , int comp = 1>
int Kaskade::ShapeFunctionSet< ctype, dimension, T, comp >::size_
protected

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