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

Vectorial \( H(\text{div}) \) conforming shape functions on the unit simplex. More...

#include <nedelecshapefunctions.hh>

Detailed Description

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

Vectorial \( H(\text{div}) \) conforming shape functions on the unit simplex.

These form a basis of \( \mathbb{P}_p^d \) of the space of vectorial polynomials of degree at most \( p, p\ge 1 \). A shape function is represented as a linear combination of Lagrangian shape functions.

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

Definition at line 333 of file nedelecshapefunctions.hh.

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

Public Types

typedef ctype CoordType
 

Public Member Functions

 HdivSimplexShapeFunction ()=default
 Default constructor. Constructs a vanishing "shape function". Not particularly useful, but it's often convenient to have a default constructible class. More...
 
template<class Coefficients >
 HdivSimplexShapeFunction (Coefficients const &coeff_, std::tuple< int, int, int, int > loc_)
 Creates a shape function as a linear combination of Lagrangian shape functions. More...
 
virtual std::unique_ptr< ShapeFunction< ctype, dim, T, dim > > clone () const
 
virtual Dune::FieldVector< T, dimevaluateFunction (Dune::FieldVector< ctype, dim > const &xi) const
 Evaluates the shape function at point xi. More...
 
virtual Dune::FieldMatrix< T, dim, dimevaluateDerivative (Dune::FieldVector< CoordType, dim > const &xi) const
 Evaluates the derivative of the shape function. More...
 
virtual Tensor< T, dim, dim, dimevaluate2ndDerivative (Dune::FieldVector< CoordType, dim > const &xi) 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...
 

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::HdivSimplexShapeFunction< ctype, dimension, T >::CoordType

Definition at line 340 of file nedelecshapefunctions.hh.

Constructor & Destructor Documentation

◆ HdivSimplexShapeFunction() [1/2]

template<class ctype , int dimension, class T >
Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::HdivSimplexShapeFunction ( )
default

Default constructor. Constructs a vanishing "shape function". Not particularly useful, but it's often convenient to have a default constructible class.

Todo:
makes only sense if it is assignable - any real need for default constructibility?

Referenced by Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::clone().

◆ HdivSimplexShapeFunction() [2/2]

template<class ctype , int dimension, class T >
template<class Coefficients >
Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::HdivSimplexShapeFunction ( Coefficients const &  coeff_,
std::tuple< int, int, int, int >  loc_ 
)
inlineexplicit

Creates a shape function as a linear combination of Lagrangian shape functions.

Template Parameters
Coefficientsan STL container with value type Dune::FieldVector<T,comps>.
Parameters
coeffcontains coefficients of linear combinations of Lagrangian shape functions. Note that the coefficients are vectorial, as the Lagrangian shape functions are scalar, but we create vectorial shape functions.
locthe topological location information (nominalOrder,codim,entity,index)

Definition at line 357 of file nedelecshapefunctions.hh.

Member Function Documentation

◆ clone()

template<class ctype , int dimension, class T >
virtual std::unique_ptr< ShapeFunction< ctype, dim, T, dim > > Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::clone ( ) const
inlinevirtual

Definition at line 363 of file nedelecshapefunctions.hh.

◆ evaluate2ndDerivative()

template<class ctype , int dimension, class T >
virtual Tensor< T, dim, dim, dim > Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::evaluate2ndDerivative ( Dune::FieldVector< CoordType, dim > const &  xi) const
inlinevirtual

Evaluates the second derivative of the shape function.

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

The point xi 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 409 of file nedelecshapefunctions.hh.

◆ evaluateDerivative()

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

Evaluates the derivative of the shape function.

The point xi 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 389 of file nedelecshapefunctions.hh.

◆ evaluateFunction()

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

Evaluates the shape function at point xi.

The point xi 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 375 of file nedelecshapefunctions.hh.

Referenced by Kaskade::HdivSimplexShapeFunction< ctype, dimension, T >::evaluateFunction().

◆ location()

template<class ctype , int dimension, class T >
virtual std::tuple< int, int, int, int > Kaskade::HdivSimplexShapeFunction< 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 421 of file nedelecshapefunctions.hh.

Member Data Documentation

◆ comps

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

Definition at line 337 of file nedelecshapefunctions.hh.

◆ dim

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

◆ order

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

Definition at line 338 of file nedelecshapefunctions.hh.


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