#include <dune/fem/space/discontinuousgalerkin/space.hh>
|
| static constexpr unsigned int | size () |
| |
|
| int | order () const |
| | return order of shape functions
|
| |
| std::size_t constexpr | size () const |
| | return number of shape functions
|
| |
| void | evaluateEach (const Point &x, Functor functor) const |
| | evalute each shape function
|
| |
| void | jacobianEach (const Point &x, Functor functor) const |
| | evalute jacobian of each shape function
|
| |
| void | hessianEach (const Point &x, Functor functor) const |
| | evalute hessian of each shape function
|
| |
| static std::size_t constexpr | size (int order) |
| | please doc me
|
| |
◆ BaseType
template<class
FunctionSpace , class GridPart , int polOrder, class Storage >
◆ DomainType
◆ FunctionSpaceType
◆ HessianRangeType
◆ JacobianRangeType
◆ RangeType
◆ ScalarShapeFunctionSet()
template<class
FunctionSpace , class GridPart , int polOrder, class Storage >
◆ evaluateEach()
evalute each shape function
- Parameters
-
| [in] | x | coordinate or quadrature point |
| [in] | functor | functor call for evaluating each shape function |
The functor has to be a copyable object satisfying the following interface:
struct Functor
{
template< class Value >
void operator() ( const int shapeFunction, const Value &value );
};
◆ hessianEach()
evalute hessian of each shape function
- Parameters
-
| [in] | x | coordinate or quadrature point |
| [in] | functor | functor call for evaluating the hessian of each shape function |
The functor has to be a copyable object satisfying the following interface:
struct Functor
{
template< class Hessian >
void operator() ( const int shapeFunction, const Hessian &hessian );
};
◆ jacobianEach()
evalute jacobian of each shape function
- Parameters
-
| [in] | x | coordinate or quadrature point |
| [in] | functor | functor call for evaluating the jacobian of each shape function |
The functor has to be a copyable object satisfying the following interface:
struct Functor
{
template< class Jacobian >
void operator() ( const int shapeFunction, const Jacobian &jacobian );
};
◆ order()
return order of shape functions
◆ size() [1/3]
template<class
FunctionSpace , class GridPart , int polOrder, class Storage >
◆ size() [2/3]
return number of shape functions
◆ size() [3/3]
|
|
inlinestaticconstexprinherited |
◆ dimension
◆ numberShapeFunctions
template<class
FunctionSpace , class GridPart , int polOrder, class Storage >
Initial value:=
OrthonormalShapeFunctions< ScalarShapeFunctionSpaceType::dimDomain >::size(polOrder)
The documentation for this struct was generated from the following file: