DUNE-FEM (unstable)

shapefunctionset.hh
Go to the documentation of this file.
1#ifndef DUNE_FEM_SHAPEFUNCTIONSET_SHAPEFUNCTIONSET_HH
2#define DUNE_FEM_SHAPEFUNCTIONSET_SHAPEFUNCTIONSET_HH
3
4//- C++ includes
5#include <cstddef>
6
13namespace Dune
14{
15
16 namespace Fem
17 {
18
19 // ShapeFunctionSet
20 // ----------------
21
31 template< class FunctionSpace >
33 {
34 public:
37
46
48 int order () const;
49
51 std::size_t size () const;
52
69 template< class Point, class Functor >
70 void evaluateEach ( const Point &x, Functor functor ) const;
71
88 template< class Point, class Functor >
89 void jacobianEach ( const Point &x, Functor functor ) const;
90
107 template< class Point, class Functor >
108 void hessianEach ( const Point &x, Functor functor ) const;
109 };
110
111 } // namespace Fem
112
113} // namespace Dune
114
115#endif // #ifndef DUNE_FEM_SHAPEFUNCTIONSET_SHAPEFUNCTIONSET_HH
Definition: explicitfieldvector.hh:75
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
FunctionSpaceTraits::LinearMappingType JacobianRangeType
Intrinsic type used for the jacobian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:75
FunctionSpaceTraits::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:67
A vector valued function space.
Definition: functionspace.hh:60
Interface class for shape function sets.
Definition: shapefunctionset.hh:33
FunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian range type
Definition: shapefunctionset.hh:43
FunctionSpace FunctionSpaceType
function space type
Definition: shapefunctionset.hh:36
FunctionSpaceType::HessianRangeType HessianRangeType
hessian range type
Definition: shapefunctionset.hh:45
FunctionSpaceType::RangeType RangeType
range type
Definition: shapefunctionset.hh:41
void hessianEach(const Point &x, Functor functor) const
evalute hessian of each shape function
void evaluateEach(const Point &x, Functor functor) const
evalute each shape function
std::size_t size() const
return number of shape functions
int order() const
return order of shape functions
FunctionSpaceType::DomainType DomainType
domain type
Definition: shapefunctionset.hh:39
void jacobianEach(const Point &x, Functor functor) const
evalute jacobian of each shape function
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)