KASKADE 7 development version
|
Scalar Lagrange shape functions on the unit simplex. More...
#include <lagrangeshapefunctions.hh>
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.
ctype | the scalar type for coordinates |
dimension | the spatial dimension |
Basis | a class defining the actual interpolation nodes and shape functions, usually either EquidistantLagrange or LobattoGridLagrange |
Scalar | the value type (a real number type) |
Definition at line 382 of file lagrangeshapefunctions.hh.
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, dim > | evaluateDerivative (Dune::FieldVector< CoordType, dim > const &x) const |
Evaluates the derivative of the shape function. More... | |
virtual Tensor< Scalar, 1, dim, dim > | evaluate2ndDerivative (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 |
typedef ctype Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::CoordType |
Definition at line 388 of file lagrangeshapefunctions.hh.
typedef Scalar Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::ResultType |
Definition at line 389 of file lagrangeshapefunctions.hh.
|
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.
|
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.
|
default |
Copy constructor.
|
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.
|
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.
|
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.
|
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.
|
static |
Definition at line 386 of file lagrangeshapefunctions.hh.
|
static |
Definition at line 385 of file lagrangeshapefunctions.hh.
Referenced by Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::evaluate2ndDerivative(), Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::evaluateDerivative(), and Kaskade::LagrangeSimplexShapeFunction< ctype, dimension, Basis, Scalar >::LagrangeSimplexShapeFunction().