dune-pdelab 2.7-git
Loading...
Searching...
No Matches
memberadaptor.hh
Go to the documentation of this file.
1// -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=8 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_FUNCTION_MEMBERADAPTOR_HH
5#define DUNE_PDELAB_FUNCTION_MEMBERADAPTOR_HH
6
7#include <cstddef>
8
9#include <dune/common/fvector.hh>
10
12
13namespace Dune {
14 namespace PDELab {
15
17
31 template<class Member, class Class,
32 class GV, class RF, std::size_t dimR = 1>
34 public GridFunctionBase<
35 GridFunctionTraits< GV, RF, dimR, FieldVector<RF, dimR> >,
36 MemberFunctionToGridFunctionAdaptor<Member, Class, GV, RF, dimR>
37 >
38 {
39 public:
42
43 private:
44 typedef GridFunctionBase<
45 Traits,
47 > Base;
48
49 const Class &obj;
50 Member member;
51 GV gv;
52
53 public:
55
63 MemberFunctionToGridFunctionAdaptor(const Class &obj_, Member member_,
64 const GV& gv_) :
65 obj(obj_), member(member_), gv(gv_)
66 { }
67
68 void evaluate(const typename Traits::ElementType &e,
69 const typename Traits::DomainType &x,
70 typename Traits::RangeType &y) const
71 {
72 (obj.*member)(e, x, y);
73 }
74
76 const GV& getGridView() const { return gv; }
77 };
78
80
90 template<class RF, std::size_t dimRange,
91 class GV, class Class, class Member>
92 MemberFunctionToGridFunctionAdaptor<Member, Class, GV, RF, dimRange>
93 makeMemberFunctionToGridFunctionAdaptor(const Class &obj, Member member,
94 const GV &gv)
95 {
97 <Member, Class, GV, RF, dimRange>(obj, member, gv);
98 }
99
101
116 template<class Member, class Class,
117 class GV, class RF, std::size_t dimR = 1>
119 public GridFunctionBase<
120 GridFunctionTraits< GV, RF, dimR, FieldVector<RF, dimR> >,
121 TwoArgsMemberFunctionToGridFunctionAdaptor<Member, Class, GV, RF, dimR>
122 >
123 {
124 public:
127
128 private:
129 typedef GridFunctionBase<
130 Traits,
132 > Base;
133
134 const Class &obj;
135 Member member;
136 GV gv;
137
138 public:
140
149 Member member_,
150 const GV& gv_) :
151 obj(obj_), member(member_), gv(gv_)
152 { }
153
154 void evaluate(const typename Traits::ElementType &e,
155 const typename Traits::DomainType &x,
156 typename Traits::RangeType &y) const
157 {
158 y = (obj.*member)(e, x);
159 }
160
162 const GV& getGridView() const { return gv; }
163 };
164
166
176 template<class RF, std::size_t dimRange,
177 class GV, class Class, class Member>
178 TwoArgsMemberFunctionToGridFunctionAdaptor<Member, Class, GV, RF,
179 dimRange>
181 Member member,
182 const GV &gv)
183 {
185 <Member, Class, GV, RF, dimRange>(obj, member, gv);
186 }
187
189
204 template<class Member, class Class,
205 class GV, class RF, std::size_t dimR = 1>
208 BoundaryGridFunctionTraits< GV, RF, dimR, FieldVector<RF, dimR> >,
209 TwoArgsMemberFunctionToBoundaryGridFunctionAdaptor<Member, Class, GV,
210 RF, dimR>
211 >
212 {
213 public:
217
218 private:
220 Traits,
222 RF, dimR>
223 > Base;
224
225 const Class &obj;
226 Member member;
227 GV gv;
228
229 public:
231
240 Member member_,
241 const GV& gv_) :
242 obj(obj_), member(member_), gv(gv_)
243 { }
244
245 template<typename Intersection>
246 void evaluate(const Intersection &is,
247 const typename Traits::DomainType &x,
248 typename Traits::RangeType &y) const
249 {
250 y = (obj.*member)(is, x);
251 }
252
254 const GV& getGridView() const { return gv; }
255 };
256
259
269 template<class RF, std::size_t dimRange,
270 class GV, class Class, class Member>
271 TwoArgsMemberFunctionToBoundaryGridFunctionAdaptor<Member, Class, GV, RF,
272 dimRange>
274 Member member,
275 const GV &gv)
276 {
278 <Member, Class, GV, RF, dimRange>(obj, member, gv);
279 }
280
281 } // namspace PDELab
282} // namspace Dune
283
284#endif // DUNE_PDELAB_FUNCTION_MEMBERADAPTOR_HH
For backward compatibility – Do not use this!
Definition adaptivity.hh:28
TwoArgsMemberFunctionToGridFunctionAdaptor< Member, Class, GV, RF, dimRange > makeTwoArgsMemberFunctionToGridFunctionAdaptor(const Class &obj, Member member, const GV &gv)
easy construction of a TwoArgsMemberFunctionToGridFunctionAdaptor
Definition memberadaptor.hh:180
TwoArgsMemberFunctionToBoundaryGridFunctionAdaptor< Member, Class, GV, RF, dimRange > make2ArgsMemberFunctionToBoundaryGridFunctionAdaptor(const Class &obj, Member member, const GV &gv)
easy construction of a TwoArgsMemberFunctionToBoundaryGridFunctionAdaptor
Definition memberadaptor.hh:273
MemberFunctionToGridFunctionAdaptor< Member, Class, GV, RF, dimRange > makeMemberFunctionToGridFunctionAdaptor(const Class &obj, Member member, const GV &gv)
easy construction of a MemberFunctionToGridFunctionAdaptor
Definition memberadaptor.hh:93
Dune::FieldVector< GV::Grid::ctype, GV::dimension > DomainType
domain type in dim-size coordinates
Definition function.hh:50
GV::Traits::template Codim< 0 >::Entity ElementType
codim 0 entity
Definition function.hh:119
traits class holding the function signature, same as in local function
Definition function.hh:183
traits class holding function signature, same as in local function
Definition function.hh:237
leaf of a function tree
Definition function.hh:302
leaf of a function tree
Definition function.hh:329
GridFunction implemented by a member function of some class.
Definition memberadaptor.hh:38
MemberFunctionToGridFunctionAdaptor(const Class &obj_, Member member_, const GV &gv_)
Construct an adaptor object.
Definition memberadaptor.hh:63
GridFunctionTraits< GV, RF, dimR, FieldVector< RF, dimR > > Traits
export traits class
Definition memberadaptor.hh:41
const GV & getGridView() const
get reference to the internal gridview.
Definition memberadaptor.hh:76
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Definition memberadaptor.hh:68
GridFunction implemented by a member function of some class.
Definition memberadaptor.hh:123
GridFunctionTraits< GV, RF, dimR, FieldVector< RF, dimR > > Traits
export traits class
Definition memberadaptor.hh:126
TwoArgsMemberFunctionToGridFunctionAdaptor(const Class &obj_, Member member_, const GV &gv_)
Construct an adaptor object.
Definition memberadaptor.hh:148
const GV & getGridView() const
get reference to the internal gridview.
Definition memberadaptor.hh:162
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Definition memberadaptor.hh:154
BoundaryGridFunction implemented by a member function of some class.
Definition memberadaptor.hh:212
const GV & getGridView() const
get reference to the internal gridview.
Definition memberadaptor.hh:254
void evaluate(const Intersection &is, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Definition memberadaptor.hh:246
BoundaryGridFunctionTraits< GV, RF, dimR, FieldVector< RF, dimR > > Traits
export traits class
Definition memberadaptor.hh:216
TwoArgsMemberFunctionToBoundaryGridFunctionAdaptor(const Class &obj_, Member member_, const GV &gv_)
Construct an adaptor object.
Definition memberadaptor.hh:239