DUNE-FEM (unstable)

proxy.hh
Go to the documentation of this file.
1#ifndef DUNE_FEM_SHAPEFUNCTIONSET_PROXY_HH
2#define DUNE_FEM_SHAPEFUNCTIONSET_PROXY_HH
3
4// C++ includes
5#include <cassert>
6#include <cstddef>
7
8#include <dune/fem/common/utility.hh>
9
17namespace Dune
18{
19
20 namespace Fem
21 {
22
23 // ShapeFunctionSetProxy
24 // ---------------------
25
26 /*
27 * \brief A proxy object converting a pointer to a shape function set to a object
28 *
29 * \tparam ShapeFunctionSet An implementation of Dune::Fem::ShapeFunctionSet
30 *
31 * \note This class has an implicit constructor from a pointer to a shape function set.
32 */
33 template< class ShapeFunctionSet >
34 class ShapeFunctionSetProxy
35 {
36 typedef ShapeFunctionSetProxy< ShapeFunctionSet > ThisType;
37
38 public:
39 typedef ShapeFunctionSet ImplementationType;
40 static const int pointSetId = detail::SelectPointSetId< ShapeFunctionSet >::value;
41
42 // if ScalarShapeFunctionSetType has a member variable codegenShapeFunctionSet then this is forwarded here
43 // otherwise this value defaults to false
44 static constexpr bool codegenShapeFunctionSet = detail::IsCodegenShapeFunctionSet< ImplementationType >::value;
45
47
48 typedef typename FunctionSpaceType::DomainType DomainType;
49 typedef typename FunctionSpaceType::RangeType RangeType;
50 typedef typename FunctionSpaceType::JacobianRangeType JacobianRangeType;
51 typedef typename FunctionSpaceType::HessianRangeType HessianRangeType;
52
53 const ImplementationType &impl () const
54 {
55 assert( shapeFunctionSet_ );
56 return *shapeFunctionSet_;
57 }
58
59 ShapeFunctionSetProxy ()
60 : shapeFunctionSet_( nullptr )
61 {}
62
63 ShapeFunctionSetProxy ( const ShapeFunctionSet *shapeFunctionSet )
64 : shapeFunctionSet_( shapeFunctionSet )
65 {}
66
67 int order () const { return impl().order(); }
68
69 std::size_t size () const { return impl().size(); }
70
71 template< class Point, class Functor >
72 void evaluateEach ( const Point &x, Functor functor ) const
73 {
74 impl().evaluateEach( x, functor );
75 }
76
77 template< class Point, class Functor >
78 void jacobianEach ( const Point &x, Functor functor ) const
79 {
80 impl().jacobianEach( x, functor );
81 }
82
83 template< class Point, class Functor >
84 void hessianEach ( const Point &x, Functor functor ) const
85 {
86 impl().hessianEach( x, functor );
87 }
88
89 private:
90 const ShapeFunctionSet *shapeFunctionSet_;
91 };
92
93 } // namespace Fem
94
95} // namespace Dune
96
97#endif // #ifndef DUNE_FEM_SHAPEFUNCTIONSET_PROXY_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
FunctionSpace FunctionSpaceType
function space type
Definition: shapefunctionset.hh:36
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)