4 #ifndef DUNE_PDELAB_FUNCTIONWRAPPERS_HH
5 #define DUNE_PDELAB_FUNCTIONWRAPPERS_HH
9 #include <dune/typetree/nodetags.hh>
31 template<
typename Engine,
typename F0,
32 typename F1 = TypeTree::EmptyNode,
33 typename F2 = TypeTree::EmptyNode,
34 typename F3 = TypeTree::EmptyNode,
35 typename F4 = TypeTree::EmptyNode,
36 typename F5 = TypeTree::EmptyNode,
37 typename F6 = TypeTree::EmptyNode,
38 typename F7 = TypeTree::EmptyNode,
39 typename F8 = TypeTree::EmptyNode,
40 typename F9 = TypeTree::EmptyNode>
44 PointwiseGridFunctionAdapter<
46 F0, F1, F2, F3, F4, F5, F6, F7, F8, F9> >
75 static void inc(
unsigned &sum,
const F*
f)
83 static void inc(
unsigned &sum,
const TypeTree::EmptyNode*
f)
89 static unsigned count() {
117 static void evaluate(
const F &
f,
unsigned &index,
118 const typename Traits::ElementType&
e,
119 const typename Traits::DomainType& x,
120 std::vector<typename Traits::RangeType>& y) {
121 f.evaluate(e, x, y[index]);
136 static void evaluate(
const TypeTree::EmptyNode &f,
unsigned &index,
137 const typename Traits::ElementType& e,
138 const typename Traits::DomainType& x,
139 std::vector<typename Traits::RangeType>& y)
159 const F1& f1_ = TypeTree::EmptyNode(),
160 const F2& f2_ = TypeTree::EmptyNode(),
161 const F3& f3_ = TypeTree::EmptyNode(),
162 const F4& f4_ = TypeTree::EmptyNode(),
163 const F5& f5_ = TypeTree::EmptyNode(),
164 const F6& f6_ = TypeTree::EmptyNode(),
165 const F7& f7_ = TypeTree::EmptyNode(),
166 const F8& f8_ = TypeTree::EmptyNode(),
167 const F9& f9_ = TypeTree::EmptyNode())
168 : engine(engine_), f0(f0_)
169 , f1(f1_), f2(f2_), f3(f3_)
170 , f4(f4_), f5(f5_), f6(f6_)
171 , f7(f7_), f8(f8_), f9(f9_)
176 inline void evaluate (
const typename Traits::ElementType& e,
177 const typename Traits::DomainType& x,
178 typename Traits::RangeType& y)
const
180 std::vector<typename Traits::RangeType> in(size);
182 evaluate(f0, index, e, x, in);
183 evaluate(f1, index, e, x, in);
184 evaluate(f2, index, e, x, in);
185 evaluate(f3, index, e, x, in);
186 evaluate(f4, index, e, x, in);
187 evaluate(f5, index, e, x, in);
188 evaluate(f6, index, e, x, in);
189 evaluate(f7, index, e, x, in);
190 evaluate(f8, index, e, x, in);
191 evaluate(f9, index, e, x, in);
192 engine.evaluate(y, in);
197 return f0.getGridView();
203 template<
typename Engine,
typename F0>
204 PointwiseGridFunctionAdapter<Engine, F0>
214 template<
typename Engine,
typename F0,
typename F1>
215 PointwiseGridFunctionAdapter<Engine, F0, F1>
226 template<
typename Engine,
typename F0,
typename F1,
typename F2>
227 PointwiseGridFunctionAdapter<Engine, F0, F1, F2>
229 const F1& f1,
const F2& f2)
238 template<
typename Engine,
typename F0,
typename F1,
typename F2,
240 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3>
242 const F1& f1,
const F2& f2,
const F3& f3)
246 (engine,f0,f1,f2,f3);
251 template<
typename Engine,
typename F0,
typename F1,
typename F2,
252 typename F3,
typename F4>
253 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4>
255 const F1& f1,
const F2& f2,
const F3& f3,
259 <Engine,F0,F1,F2,F3,F4>
260 (engine,f0,f1,f2,f3,f4);
265 template<
typename Engine,
typename F0,
typename F1,
typename F2,
266 typename F3,
typename F4,
typename F5>
267 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4, F5>
269 const F1& f1,
const F2& f2,
const F3& f3,
270 const F4& f4,
const F5& f5)
273 <Engine,F0,F1,F2,F3,F4,F5>
274 (engine,f0,f1,f2,f3,f4,f5);
279 template<
typename Engine,
typename F0,
typename F1,
typename F2,
280 typename F3,
typename F4,
typename F5,
typename F6>
281 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4, F5, F6>
283 const F1& f1,
const F2& f2,
const F3& f3,
284 const F4& f4,
const F5& f5,
const F6& f6)
287 <Engine,F0,F1,F2,F3,F4,F5,F6>
288 (engine,f0,f1,f2,f3,f4,f5,f6);
293 template<
typename Engine,
typename F0,
typename F1,
typename F2,
294 typename F3,
typename F4,
typename F5,
typename F6,
typename F7>
295 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4, F5, F6, F7>
297 const F1& f1,
const F2& f2,
const F3& f3,
298 const F4& f4,
const F5& f5,
const F6& f6,
302 <Engine,F0,F1,F2,F3,F4,F5,F6,F7>
303 (engine,f0,f1,f2,f3,f4,f5,f6,f7);
308 template<
typename Engine,
typename F0,
typename F1,
typename F2,
309 typename F3,
typename F4,
typename F5,
typename F6,
typename F7,
311 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4, F5, F6, F7, F8>
313 const F1& f1,
const F2& f2,
const F3& f3,
314 const F4& f4,
const F5& f5,
const F6& f6,
315 const F7& f7,
const F8& f8)
318 <Engine,F0,F1,F2,F3,F4,F5,F6,F7,F8>
319 (engine,f0,f1,f2,f3,f4,f5,f6,f7,f8);
324 template<
typename Engine,
typename F0,
typename F1,
typename F2,
325 typename F3,
typename F4,
typename F5,
typename F6,
typename F7,
326 typename F8,
typename F9>
327 PointwiseGridFunctionAdapter<Engine, F0, F1, F2, F3, F4, F5, F6, F7, F8,
330 const F1& f1,
const F2& f2,
const F3& f3,
331 const F4& f4,
const F5& f5,
const F6& f6,
332 const F7& f7,
const F8& f8,
const F9& f9)
335 <Engine,F0,F1,F2,F3,F4,F5,F6,F7,F8,F9>
336 (engine,f0,f1,f2,f3,f4,f5,f6,f7,f8,f9);
357 template<
typename Domain,
typename Range>
359 const std::vector<Domain>& in)
const;
384 template<
typename Domain,
typename Range>
386 const std::vector<Domain>& in)
const {
387 assert(in.size() == 1);
397 PointwiseScaleAdapterEngine<S>
413 template<
typename Domain,
typename Range>
415 const std::vector<Domain>& in)
const {
417 for(
unsigned i = 0; i < in.size(); ++i)
427 #endif // DUNE_PDELAB_FUNCTIONWRAPPERS_HH
Interface of a pointwise adapter engine.
Definition: functionwrappers.hh:345
void evaluate(Range &out, const std::vector< Domain > &in) const
calculate the adapted value from a set of input values
Definition: functionwrappers.hh:385
F0::Traits Traits
Definition: functionwrappers.hh:49
void evaluate(Range &out, const std::vector< Domain > &in) const
calculate the adapted value from a set of input values
PointwiseGridFunctionAdapter(const Engine &engine_, const F0 &f0_,...)
construct a PointwiseGridFunctionAdapter
Definition: functionwrappers.hh:155
Sum all terms.
Definition: functionwrappers.hh:408
PointwiseGridFunctionAdapter< Engine, F0 > makePointwiseGridFunctionAdapter(const Engine &engine, const F0 &f0)
Definition: functionwrappers.hh:205
PointwiseScaleAdapterEngine(const S scale_)
create a PointwiseScaleAdapterEngine
Definition: functionwrappers.hh:379
PointwiseScaleAdapterEngine< S > makePointwiseScaleAdapterEngine(const S scale)
syntactic sugar to create a PointwiseScaleAdapterEngine
Definition: functionwrappers.hh:398
Scale the output value.
Definition: functionwrappers.hh:369
a GridFunction maps x in DomainType to y in RangeType
Definition: function.hh:187
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Definition: functionwrappers.hh:176
Definition: functionwrappers.hh:41
const F & f
Definition: common/constraints.hh:145
const Traits::GridViewType & getGridView() const
Definition: functionwrappers.hh:195
void evaluate(Range &out, const std::vector< Domain > &in) const
calculate the adapted value from a set of input values
Definition: functionwrappers.hh:414
const E & e
Definition: interpolate.hh:172