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

Vectorial Nedelec shape functions on the unit simplex. More...

#include <nedelecshapefunctions.hh>

Detailed Description

template<class ctype, int dimension, class T>
class Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >

Vectorial Nedelec shape functions on the unit simplex.

These are linear shape functions associated to edges. In barycentric coordinates \( \lambda \), with \( \lambda_i = x_i \) for \( i=0,\dots, \mathrm{dim}-1 \) and \( \lambda_\mathrm{dim} = 1-\sum_{i=0}^{\mathrm{dim}-1} x_i \) on the unit simplex, the shape function associated to the edge \( e = (p_i,p_j) \) is

\[ \phi_e = \pm (\lambda_i \nabla \lambda_j - \lambda_j \nabla \lambda i). \]

Here, \( p_k \) is the vertex with \( \lambda_k = 1 \). The sign of the shape function depends on the numbering of vertices relative to the edge implemented by the grid manager.

Template Parameters
ctypethe coordinate type (a real number type)
Tthe value type (a real number type)
Dthe dimension of the unit simplex over which the shape functions are defined

Definition at line 45 of file nedelecshapefunctions.hh.

Inheritance diagram for Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >:
Kaskade::ShapeFunction< ctype, dimension, T, dimension >

Public Types

typedef ctype CoordType
 
typedef T ResultType
 

Public Member Functions

 NedelecSimplexShapeFunction ()
 
 NedelecSimplexShapeFunction (int e)
 
virtual Dune::FieldVector< T, dimevaluateFunction (Dune::FieldVector< ctype, dim > const &x) const
 Evaluates the shape function at point x. More...
 
virtual Dune::FieldMatrix< ResultType, dim, dimevaluateDerivative (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the derivative of the shape function. More...
 
virtual Tensor< ResultType, 1, dim, dimevaluate2ndDerivative (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the second derivative of the shape function. More...
 
virtual std::tuple< int, int, int, int > location () const
 Returns a tuple (nominalOrder,codim,entity,index) giving detailed information about the location of the shape function. More...
 
Dune::FieldVector< CoordType, dimposition () const
 Returns the element-local position of the node associated to the shape function. More...
 
Dune::FieldVector< CoordType, dimedge () const
 Returns the edge vector of the edge associated to this shape function. More...
 

Static Public Attributes

static int const dim = dimension
 
static int const comps = dim
 
static int const order = 1
 

Member Typedef Documentation

◆ CoordType

template<class ctype , int dimension, class T >
typedef ctype Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::CoordType

Definition at line 52 of file nedelecshapefunctions.hh.

◆ ResultType

template<class ctype , int dimension, class T >
typedef T Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::ResultType

Definition at line 53 of file nedelecshapefunctions.hh.

Constructor & Destructor Documentation

◆ NedelecSimplexShapeFunction() [1/2]

template<class ctype , int dimension, class T >
Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::NedelecSimplexShapeFunction ( )
inline

Default constructor. Constructs the shape function associated with the edge 0. Not particularly useful, but it's often convenient to have a default constructible class.

Definition at line 60 of file nedelecshapefunctions.hh.

◆ NedelecSimplexShapeFunction() [2/2]

template<class ctype , int dimension, class T >
Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::NedelecSimplexShapeFunction ( int  e)
inlineexplicit

Creates a shape function associated to the edge e.

Parameters
esubentity index of the edge (codim=dim-1 subentity index)

Definition at line 70 of file nedelecshapefunctions.hh.

Member Function Documentation

◆ edge()

template<class ctype , int dimension, class T >
Dune::FieldVector< CoordType, dim > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::edge ( ) const
inline

Returns the edge vector of the edge associated to this shape function.

Definition at line 175 of file nedelecshapefunctions.hh.

◆ evaluate2ndDerivative()

template<class ctype , int dimension, class T >
virtual Tensor< ResultType, 1, dim, dim > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::evaluate2ndDerivative ( Dune::FieldVector< CoordType, dim > const &  x) const
inlinevirtual

Evaluates the second derivative of the shape function.

The result is r[i][j][k] = \( \partial^2 \phi_i / (\partial x_j \partial x_k) = 0. \)

The point x is not restricted to be inside the unit simplex, but the meaning of evaluating a shape function outside of the unit simplex is questionable.

Reimplemented from Kaskade::ShapeFunction< ctype, dimension, T, dimension >.

Definition at line 156 of file nedelecshapefunctions.hh.

◆ evaluateDerivative()

template<class ctype , int dimension, class T >
virtual Dune::FieldMatrix< ResultType, dim, dim > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::evaluateDerivative ( Dune::FieldVector< CoordType, dim > const &  x) const
inlinevirtual

Evaluates the derivative of the shape function.

The point x is not restricted to be inside the unit simplex, but the meaning of evaluating a shape function outside of the unit simplex is questionable.

Implements Kaskade::ShapeFunction< ctype, dimension, T, dimension >.

Definition at line 115 of file nedelecshapefunctions.hh.

◆ evaluateFunction()

template<class ctype , int dimension, class T >
virtual Dune::FieldVector< T, dim > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::evaluateFunction ( Dune::FieldVector< ctype, dim > const &  x) const
inlinevirtual

Evaluates the shape function at point x.

The point x is not restricted to be inside the unit simplex, but the meaning of evaluating a shape function outside of the unit simplex is questionable.

Implements Kaskade::ShapeFunction< ctype, dimension, T, dimension >.

Definition at line 102 of file nedelecshapefunctions.hh.

◆ location()

template<class ctype , int dimension, class T >
virtual std::tuple< int, int, int, int > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::location ( ) const
inlinevirtual

Returns a tuple (nominalOrder,codim,entity,index) giving detailed information about the location of the shape function.

Each shape function is associated to a certain subentity of the element.

nominalOrder is a nonnegative ordering parameter that is usually the polynomial order of the shape function, but need not coincide.

codim is the codimension of the subentity to which the shape function is associated, entity is the number of the subentity, and index is the number of the shape function among those that are associated to the same subentity.

Implements Kaskade::ShapeFunction< ctype, dimension, T, dimension >.

Definition at line 162 of file nedelecshapefunctions.hh.

◆ position()

template<class ctype , int dimension, class T >
Dune::FieldVector< CoordType, dim > Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::position ( ) const
inline

Returns the element-local position of the node associated to the shape function.

Definition at line 170 of file nedelecshapefunctions.hh.

Member Data Documentation

◆ comps

template<class ctype , int dimension, class T >
int const Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::comps = dim
static

Definition at line 49 of file nedelecshapefunctions.hh.

◆ dim

template<class ctype , int dimension, class T >
int const Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::dim = dimension
static

◆ order

template<class ctype , int dimension, class T >
int const Kaskade::NedelecSimplexShapeFunction< ctype, dimension, T >::order = 1
static

Definition at line 50 of file nedelecshapefunctions.hh.


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