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

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

#include <lagrangeshapefunctions.hh>

Detailed Description

template<class ctype, int dimension, class Basis, class Scalar = double>
class Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >

Scalar Lagrange shape functions on the unit simplex.

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). Each node is referenced by its tuple index: a dim-tuple of integers with values between 0 and order (including), such that the total sum is at most order.

Template Parameters
ctypethe scalar type for coordinates
dimensionthe spatial dimension
Basisa class defining the actual interpolation nodes and shape functions, usually either EquidistantLagrange or LobattoGridLagrange
Scalarthe value type (a real number type)

Definition at line 382 of file lagrangeshapefunctions.hh.

Inheritance diagram for Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >:
Kaskade::ShapeFunction< ctype, dimension, double >

Public Types

typedef ctype CoordType
 
typedef Scalar ResultType
 

Public Member Functions

 LagrangeSimplexShapeFunction ()
 Default constructor. More...
 
 LagrangeSimplexShapeFunction (int local, std::shared_ptr< Basis const > basis_)
 Creates a shape function associated to the node with tuple index xi_. The entries in xi_ have to be nonnegative, and their sum must not exceed order. More...
 
 LagrangeSimplexShapeFunction (LagrangeSimplexShapeFunction const &other)=default
 Copy constructor. More...
 
virtual Dune::FieldVector< Scalar, 1 > evaluateFunction (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the shape function at point x. More...
 
virtual Dune::FieldMatrix< Scalar, 1, dimevaluateDerivative (Dune::FieldVector< CoordType, dim > const &x) const
 Evaluates the derivative of the shape function. More...
 
virtual Tensor< Scalar, 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...
 

Static Public Attributes

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

Member Typedef Documentation

◆ CoordType

template<class ctype , int dimension, class Basis , class Scalar = double>
typedef ctype Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::CoordType

Definition at line 388 of file lagrangeshapefunctions.hh.

◆ ResultType

template<class ctype , int dimension, class Basis , class Scalar = double>
typedef Scalar Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::ResultType

Definition at line 389 of file lagrangeshapefunctions.hh.

Constructor & Destructor Documentation

◆ LagrangeSimplexShapeFunction() [1/3]

template<class ctype , int dimension, class Basis , class Scalar = double>
Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::LagrangeSimplexShapeFunction ( )
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 398 of file lagrangeshapefunctions.hh.

◆ LagrangeSimplexShapeFunction() [2/3]

template<class ctype , int dimension, class Basis , class Scalar = double>
Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::LagrangeSimplexShapeFunction ( int  local,
std::shared_ptr< Basis const >  basis_ 
)
inlineexplicit

Creates a shape function associated to the node with tuple index xi_. The entries in xi_ have to be nonnegative, and their sum must not exceed order.

Definition at line 407 of file lagrangeshapefunctions.hh.

◆ LagrangeSimplexShapeFunction() [3/3]

template<class ctype , int dimension, class Basis , class Scalar = double>
Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::LagrangeSimplexShapeFunction ( LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar > const &  other)
default

Copy constructor.

Member Function Documentation

◆ evaluate2ndDerivative()

template<class ctype , int dimension, class Basis , class Scalar = double>
virtual Tensor< Scalar, 1, dim, dim > Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::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) \)

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, double >.

Definition at line 486 of file lagrangeshapefunctions.hh.

◆ evaluateDerivative()

template<class ctype , int dimension, class Basis , class Scalar = double>
virtual Dune::FieldMatrix< Scalar, 1, dim > Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::evaluateDerivative ( Dune::FieldVector< CoordType, dim > const &  x) const
inlinevirtual

Evaluates the derivative of the shape function.

The result is r[i][j] = \( \partial \phi_i / \partial x_j \)

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, double >.

Definition at line 469 of file lagrangeshapefunctions.hh.

◆ evaluateFunction()

template<class ctype , int dimension, class Basis , class Scalar = double>
virtual Dune::FieldVector< Scalar, 1 > Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::evaluateFunction ( Dune::FieldVector< CoordType, 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, double >.

Definition at line 455 of file lagrangeshapefunctions.hh.

◆ location()

template<class ctype , int dimension, class Basis , class Scalar = double>
virtual std::tuple< int, int, int, int > Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::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, double >.

Definition at line 501 of file lagrangeshapefunctions.hh.

Member Data Documentation

◆ comps

template<class ctype , int dimension, class Basis , class Scalar = double>
int const Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::comps = 1
static

Definition at line 386 of file lagrangeshapefunctions.hh.

◆ dim

template<class ctype , int dimension, class Basis , class Scalar = double>
int const Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::dim = dimension
static

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