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

Scalar Lagrange shape functions on the unit cube. More...

#include <lagrangeshapefunctions.hh>

Detailed Description

template<class ctype, int dimension, class Scalar, int O>
class Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >

Scalar Lagrange shape functions on the unit cube.

These are polynomial shape functions of fixed given order associated to certain interpolation nodes. A shape function is 1 at its associated node and 0 at all the other nodes.

The nodes form a cartesian grid with (order+1) nodes on each axis. The node placement is isotropic (the same node placement on each axis). In the current implementation, quasi Chebyshev nodes are used (due to their better interpolation properties compared to equidistant nodes). Each node is referenced by its index: a dim-tuple of integers with values between 0 and order (including), such that the maximum entry is at most order.

Template Parameters
ctypethe scalar type for coordinates
dimensionthe spatial dimension
Scalarthe value type (a real number type)
Othe polynomial order

Definition at line 613 of file lagrangeshapefunctions.hh.

Inheritance diagram for Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >:
Kaskade::ShapeFunction< ctype, dimension, Scalar >

Public Types

typedef ctype CoordType
 
typedef Scalar ResultType
 

Public Member Functions

 LagrangeCubeShapeFunction ()
 
 LagrangeCubeShapeFunction (LagrangeCubeShapeFunction const &other)
 
 LagrangeCubeShapeFunction (Dune::FieldVector< int, dimension > const &xi_)
 Creates a shape function associated to the node with indices xi_. More...
 
ResultType evaluateFunction (int, Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the shape function at point x. More...
 
virtual Dune::FieldVector< Scalar, 1 > evaluateFunction (Dune::FieldVector< ctype, dim > const &x) const
 Evaluates the shape function (all components at once). More...
 
ResultType evaluateDerivative (int, int dir, Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the partial derivative of the shape function in coordinate direction dir.
More...
 
virtual Dune::FieldMatrix< ResultType, 1, dimevaluateDerivative (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the derivative of the shape function (all components and all directions at once). More...
 
virtual Tensor< ResultType, 1, dim, dimevaluate2ndDerivative (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the second derivative of the shape function (all components and all directions at once). More...
 
int localindex (int) const
 Returns the local index of the shape function. More...
 
int codim () const
 Returns the codimension of the subentity on which the associated node is located. More...
 
int entity () const
 Returns the subentity number of the smallest codimension subentity on which the associated node is located. More...
 
template<class Cell >
int entityIndex (Cell const &) const
 Returns the local index on the subentity. More...
 
Dune::FieldVector< CoordType, dimposition () const
 Returns the element-local position of the node associated to 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 = 1
 
static int const order = O
 

Member Typedef Documentation

◆ CoordType

template<class ctype , int dimension, class Scalar , int O>
typedef ctype Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::CoordType

Definition at line 620 of file lagrangeshapefunctions.hh.

◆ ResultType

template<class ctype , int dimension, class Scalar , int O>
typedef Scalar Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::ResultType

Definition at line 621 of file lagrangeshapefunctions.hh.

Constructor & Destructor Documentation

◆ LagrangeCubeShapeFunction() [1/3]

template<class ctype , int dimension, class Scalar , int O>
Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::LagrangeCubeShapeFunction ( )
inline

Default constructor. Constructs the shape function located at the origin of the reference element. Not particularly useful, but it's often convenient to have a default constructible class.

Definition at line 628 of file lagrangeshapefunctions.hh.

◆ LagrangeCubeShapeFunction() [2/3]

template<class ctype , int dimension, class Scalar , int O>
Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::LagrangeCubeShapeFunction ( LagrangeCubeShapeFunction< ctype, dimension, Scalar, O > const &  other)
inline

Definition at line 632 of file lagrangeshapefunctions.hh.

◆ LagrangeCubeShapeFunction() [3/3]

template<class ctype , int dimension, class Scalar , int O>
Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::LagrangeCubeShapeFunction ( Dune::FieldVector< int, dimension > const &  xi_)
explicit

Creates a shape function associated to the node with indices xi_.

The entries in xi_ have to be nonnegative, and their maximum must not exceed order.

Member Function Documentation

◆ codim()

template<class ctype , int dimension, class Scalar , int O>
int Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::codim ( ) const
inline

Returns the codimension of the subentity on which the associated node is located.

Definition at line 702 of file lagrangeshapefunctions.hh.

◆ entity()

template<class ctype , int dimension, class Scalar , int O>
int Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::entity ( ) const
inline

Returns the subentity number of the smallest codimension subentity on which the associated node is located.

Definition at line 708 of file lagrangeshapefunctions.hh.

◆ entityIndex()

template<class ctype , int dimension, class Scalar , int O>
template<class Cell >
int Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::entityIndex ( Cell const &  ) const
inline

Returns the local index on the subentity.

For higher order shape functions, globally unique numberings on intersections have to be used, therefore this depends on the actual element and an index set that covers this element.

Definition at line 718 of file lagrangeshapefunctions.hh.

◆ evaluate2ndDerivative()

template<class ctype , int dimension, class Scalar , int O>
virtual Tensor< ResultType, 1, dim, dim > Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::evaluate2ndDerivative ( Dune::FieldVector< CoordType, dim > const &  xi) const
inlinevirtual

Evaluates the second derivative of the shape function (all components and all directions at once).

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

Definition at line 685 of file lagrangeshapefunctions.hh.

◆ evaluateDerivative() [1/2]

template<class ctype , int dimension, class Scalar , int O>
virtual Dune::FieldMatrix< ResultType, 1, dim > Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::evaluateDerivative ( Dune::FieldVector< CoordType, dim > const &  xi) const
inlinevirtual

Evaluates the derivative of the shape function (all components and all directions at once).

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

Definition at line 676 of file lagrangeshapefunctions.hh.

◆ evaluateDerivative() [2/2]

template<class ctype , int dimension, class Scalar , int O>
ResultType Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::evaluateDerivative ( int  ,
int  dir,
Dune::FieldVector< CoordType, dim > const &  x 
) const
inline

Evaluates the partial derivative of the shape function in coordinate direction dir.

0 <= dir < dim has to hold. The point x is not restricted to be inside the unit cube, but the meaning of evaluating a shape function outside of the unit cube is questionable.

Definition at line 670 of file lagrangeshapefunctions.hh.

◆ evaluateFunction() [1/2]

template<class ctype , int dimension, class Scalar , int O>
virtual Dune::FieldVector< Scalar, 1 > Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::evaluateFunction ( Dune::FieldVector< ctype, dim > const &  xi) const
inlinevirtual

Evaluates the shape function (all components at once).

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

Definition at line 657 of file lagrangeshapefunctions.hh.

◆ evaluateFunction() [2/2]

template<class ctype , int dimension, class Scalar , int O>
ResultType Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::evaluateFunction ( int  ,
Dune::FieldVector< CoordType, dim > const &  x 
) const
inline

Evaluates the shape function at point x.

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

Definition at line 652 of file lagrangeshapefunctions.hh.

◆ localindex()

template<class ctype , int dimension, class Scalar , int O>
int Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::localindex ( int  ) const
inline

Returns the local index of the shape function.

Definition at line 697 of file lagrangeshapefunctions.hh.

◆ location()

template<class ctype , int dimension, class Scalar , int O>
virtual std::tuple< int, int, int, int > Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::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, Scalar >.

Definition at line 744 of file lagrangeshapefunctions.hh.

◆ position()

template<class ctype , int dimension, class Scalar , int O>
Dune::FieldVector< CoordType, dim > Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::position ( ) const
inline

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

Definition at line 736 of file lagrangeshapefunctions.hh.

Member Data Documentation

◆ comps

template<class ctype , int dimension, class Scalar , int O>
int const Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::comps = 1
static

Definition at line 617 of file lagrangeshapefunctions.hh.

◆ dim

template<class ctype , int dimension, class Scalar , int O>
int const Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::dim = dimension
static

◆ order

template<class ctype , int dimension, class Scalar , int O>
int const Kaskade::LagrangeCubeShapeFunction< ctype, dimension, Scalar, O >::order = O
static

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