DUNE-ACFEM (2.5.1)

gridfunctionexpression.hh
1#ifndef __DUNE_ACFEM_GRIDFUNCTIONEXPRESSION_HH__
2#define __DUNE_ACFEM_GRIDFUNCTIONEXPRESSION_HH__
3
4#include <dune/fem/version.hh>
5#include <dune/fem/function/common/discretefunction.hh>
6
7#include "../functions/gridfunctionexpressionbase.hh"
8#include "../functions/parameterfunction.hh"
9#include "../functions/functionoperations.hh"
10#include "../functions/functionexpression.hh"
11#include "../functions/localfunctionwrapper.hh"
12#include "../functions/gridfunctionwrapper.hh"
13#include "../functions/gridfunctionwrapperexpression.hh"
14#include "../functions/boundarysupportedfunction.hh"
15#include "../expressions/parameterinterface.hh"
16#include "../expressions/parameterexpression.hh"
17
18namespace Dune
19{
20
21 namespace ACFem
22 {
23
58 template<class GridFunction>
59 static inline
60 const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>&
61 asExprFunction(const GridFunction& arg)
62 {
63 typedef Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction> Interface;
64 return static_cast<const Interface&>(arg);
65 };
66
71 template<class UnOp, class Function>
72 static inline
73 const UnaryGridFunctionExpression<UnOp, Function>&
75 {
76 return arg;
77 };
78
83 template<class BinOp, class LeftFunction, class RightFunction>
84 static inline
85 const BinaryGridFunctionExpression<BinOp, LeftFunction, RightFunction>&
87 {
88 return arg;
89 };
90
110 template<class F1, class F2>
111 static inline
112 ZeroGridFunction<
113 typename FunctionMultiplicationResultTraits<typename F1::FunctionSpaceType,
114 typename F2::FunctionSpaceType>::FunctionSpaceType,
115 typename F1::GridPartType>
116 zeroProduct(const F1& f1, const F2& f2)
117 {
118 typedef
119 FunctionMultiplicationResultTraits<typename F1::FunctionSpaceType,
120 typename F2::FunctionSpaceType>
121 TraitsType;
122 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
123 typedef typename F1::GridPartType GridPartType;
124
125 return
127 }
128
130 template<class GridFunction>
131 typename GridFunction::LocalFunctionType localFunction(const GridFunction& f)
132 {
133 return typename GridFunction::LocalFunctionType(f);
134 }
135
136 template<class DiscreteSpace, class Dof = typename DiscreteSpace::RangeFieldType>
137 Fem::TemporaryLocalFunction<DiscreteSpace, Dof>
138 temporaryLocalFunction(const DiscreteSpace& space, const Dof& dummy = Dof())
139 {
140 return Fem::TemporaryLocalFunction<DiscreteSpace, Dof>(space);
141 }
142
155 template<class FunctionSpace, class GridPart>
156 FractionGridFunction<typename FunctionSpace::FunctionSpaceType::ScalarFunctionSpaceType, GridPart, 1L, 1UL>
157 oneFunction(const FunctionSpace& space, const GridPart& gridPart)
158 {
159 typedef
161 ResultType;
162
163 return ResultType(gridPart);
164 }
165
174 template<class F>
175 auto oneFunction(const F& f)
176 -> decltype(oneFunction(f.space(), f.gridPart()))
177 {
178 return oneFunction(f.space(), f.gridPart());
179 }
180
182 template<class F>
183 auto zeroFunction(const Fem::Function<typename F::FunctionSpaceType, F>& f_)
185 {
187 }
188
190 template<class F>
191 auto constantFunction(const Fem::Function<typename F::FunctionSpaceType, F>& f_, const typename F::FunctionSpaceType::RangeType& c)
193 {
195 }
196
197 /******************************************************************************
198 *
199 * Traits and helpers
200 *
201 */
202
203 template<class BinOp, class LeftFunctionType, class RightFunctionType>
204 struct BinaryGridFunctionExpressionBaseTraits;
205
206 template<class LeftFunction, class RightFunction>
207 struct BinaryGridFunctionExpressionBaseTraits<PlusOperation, LeftFunction, RightFunction>
208 : public BinaryFunctionExpressionTraits<PlusOperation, LeftFunction, RightFunction>
209 {
210 typedef LeftFunction LeftFunctionType;
211 typedef RightFunction RightFunctionType;
212 typedef typename LeftFunctionType::DiscreteFunctionSpaceType LeftDiscreteFunctionSpaceType;
213 typedef typename RightFunctionType::DiscreteFunctionSpaceType RightDiscreteFunctionSpaceType;
214
215 static_assert(std::is_same<typename LeftFunctionType::FunctionSpaceType,
216 typename RightFunctionType::FunctionSpaceType>::value,
217 "Cannot add functions over differing function spaces");
218 static_assert(std::is_same<typename LeftFunctionType::GridPartType,
219 typename RightFunctionType::GridPartType>::value,
220 "Cannot add functions over different grid-parts");
221 typedef
222 Fem::DiscreteFunctionSpaceAdapter<typename LeftFunctionType::FunctionSpaceType,
223 typename LeftFunctionType::GridPartType>
224 DiscreteFunctionSpaceType;
225
226 // p.w. constant tag
227 typedef
228 std::integral_constant<bool, std::is_base_of<IsPieceWiseConstant, LeftFunction>::value>
229 LeftPieceWiseConstantType;
230 typedef
231 std::integral_constant<bool, std::is_base_of<IsPieceWiseConstant, RightFunction>::value>
232 RightPieceWiseConstantType;
233 typedef
234 std::integral_constant<bool, LeftPieceWiseConstantType::value && RightPieceWiseConstantType::value>
235 PieceWiseConstantType;
236
237 typedef
238 typename std::conditional<PieceWiseConstantType::value,
239 IsPieceWiseConstant, IsNotPieceWiseConstant>::type
240 PieceWiseConstantTagType;
241
242 // max. polynomial order
243 constexpr static unsigned combinedOrder(unsigned order1, unsigned order2) {
244 return PieceWiseConstantType() ? 0 : (order1 > order2 ? order1 : order2);
245 }
246
247 // Note that std::max is not const ... grin.
248 static const unsigned int order = combinedOrder(LeftDiscreteFunctionSpaceType::polynomialOrder,
249 RightDiscreteFunctionSpaceType::polynomialOrder);
250
251 // BoundarySupportedFunction support
252 typedef BoundaryFunctionTraits<LeftFunction> LeftBoundaryTraits;
253 typedef BoundaryFunctionTraits<RightFunction> RightBoundaryTraits;
254 enum {
255 hasBoundarySupport = LeftBoundaryTraits::hasBoundarySupport || RightBoundaryTraits::hasBoundarySupport
256 };
257 typedef
258 typename std::conditional<hasBoundarySupport,
259 HasBoundarySupport, HasNoBoundarySupport>::type
260 BoundarySupportTagType;
261
262 typedef typename LeftBoundaryTraits::IndicatorType LeftIndicatorType;
263 typedef typename RightBoundaryTraits::IndicatorType RightIndicatorType;
264 typedef
265 typename std::remove_reference<
266 typename std::remove_const<decltype(std::declval<LeftIndicatorType>()
267 ||
268 std::declval<RightIndicatorType>())>::type
269 >::type
270 IndicatorType;
271
272 enum {
273 emptySupport = ExpressionTraits<IndicatorType>::isZero,
274 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
275 };
276
277 static const bool boundarySupported(bool left, bool right)
278 {
279 return globalSupport || (!emptySupport && (left || right));
280 }
281
282 static IndicatorType indicator(const LeftFunction& left, const RightFunction& right)
283 {
284 return LeftBoundaryTraits::indicator(left) || RightBoundaryTraits::indicator(right);
285 }
286 };
287
288 template<class LeftFunction, class RightFunction>
289 struct BinaryGridFunctionExpressionBaseTraits<MinusOperation, LeftFunction, RightFunction>
290 : public BinaryGridFunctionExpressionBaseTraits<PlusOperation, LeftFunction, RightFunction>
291 {};
292
293 template<class LeftFunction, class RightFunction>
294 struct BinaryGridFunctionExpressionBaseTraits<MultiplyOperation, LeftFunction, RightFunction>
295 {
296 typedef LeftFunction LeftFunctionType;
297 typedef RightFunction RightFunctionType;
298 typedef typename LeftFunctionType::FunctionSpaceType LeftFunctionSpaceType;
299 typedef typename RightFunctionType::FunctionSpaceType RightFunctionSpaceType;
300
301 typedef typename
302 FunctionMultiplicationResultTraits<LeftFunctionSpaceType, RightFunctionSpaceType>::FunctionSpaceType
303 FunctionSpaceType;
304
305 typedef typename LeftFunctionType::DiscreteFunctionSpaceType LeftDiscreteFunctionSpaceType;
306 typedef typename RightFunctionType::DiscreteFunctionSpaceType RightDiscreteFunctionSpaceType;
307
308 typedef typename LeftDiscreteFunctionSpaceType::GridPartType GridPartType;
309
310 typedef Fem::DiscreteFunctionSpaceAdapter<FunctionSpaceType, GridPartType> DiscreteFunctionSpaceType;
311
312 // p.w. constant tag
313 typedef
314 std::integral_constant<bool, std::is_base_of<IsPieceWiseConstant, LeftFunction>::value>
315 LeftPieceWiseConstantType;
316 typedef
317 std::integral_constant<bool, std::is_base_of<IsPieceWiseConstant, RightFunction>::value>
318 RightPieceWiseConstantType;
319 typedef
320 std::integral_constant<bool, LeftPieceWiseConstantType::value && RightPieceWiseConstantType::value>
321 PieceWiseConstantType;
322
323 typedef
324 typename std::conditional<PieceWiseConstantType::value,
325 IsPieceWiseConstant, IsNotPieceWiseConstant>::type
326 PieceWiseConstantTagType;
327
328 // Compute the new polynomial order "guess"
329 constexpr static unsigned combinedOrder(unsigned order1, unsigned order2)
330 {
331 return PieceWiseConstantType() ? 0 : order1 + order2;
332 }
333
334 // polynomialOrder is only an upper bound, and in the special case
335 // of DiscreteFunctionSpaceAdapter tied to 111 (why not 42?)
336 static const unsigned int order = combinedOrder(LeftDiscreteFunctionSpaceType::polynomialOrder,
337 RightDiscreteFunctionSpaceType::polynomialOrder);
338
339 // BoundarySupportedFunction support
340 typedef BoundaryFunctionTraits<LeftFunction> LeftBoundaryTraits;
341 typedef BoundaryFunctionTraits<RightFunction> RightBoundaryTraits;
342 enum {
343 hasBoundarySupport = LeftBoundaryTraits::hasBoundarySupport || RightBoundaryTraits::hasBoundarySupport
344 };
345 typedef
346 typename std::conditional<hasBoundarySupport,
347 HasBoundarySupport, HasNoBoundarySupport>::type
348 BoundarySupportTagType;
349 typedef typename LeftBoundaryTraits::IndicatorType LeftIndicatorType;
350 typedef typename RightBoundaryTraits::IndicatorType RightIndicatorType;
351 typedef
352 typename std::remove_reference<
353 typename std::remove_const<decltype(std::declval<LeftIndicatorType>()
354 &&
355 std::declval<RightIndicatorType>())>::type
356 >::type
357 IndicatorType;
358 enum {
359 emptySupport = ExpressionTraits<IndicatorType>::isZero,
360 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
361 };
362
363 static const bool boundarySupported(bool left, bool right)
364 {
365 return globalSupport || (!emptySupport && (left && right));
366 }
367
368 static IndicatorType indicator(const LeftFunction& left, const RightFunction& right)
369 {
370 return LeftBoundaryTraits::indicator(left) && RightBoundaryTraits::indicator(right);
371 }
372 };
373
375 template<class BinOp, class LeftFunction, class RightFunction>
377 : public BinaryGridFunctionExpressionBaseTraits<BinOp, LeftFunction, RightFunction>
378 {
379 typedef
380 BinaryGridFunctionExpressionBaseTraits<BinOp, LeftFunction, RightFunction>
381 BaseType;
382
383 typedef LeftFunction LeftFunctionType;
384 typedef RightFunction RightFunctionType;
385
386 typedef typename BaseType::FunctionSpaceType FunctionSpaceType;
387 typedef typename BaseType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
388
389 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
390 typedef typename FunctionSpaceType::DomainFieldType DomainFieldType;
391 typedef typename FunctionSpaceType::RangeType RangeType;
392 typedef typename FunctionSpaceType::DomainType DomainType;
393 typedef typename FunctionSpaceType::JacobianRangeType JacobianRangeType;
394 typedef typename FunctionSpaceType::HessianRangeType HessianRangeType;
395
396 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
397 typedef typename GridPartType::GridType GridType;
398 typedef typename GridPartType::template Codim<0>::EntityType EntityType;
400 typedef typename GridPartType::template Codim<0>::IteratorType IteratorType;
402 typedef typename GridPartType::IndexSetType IndexSetType;
403
405
406 using BaseType::order;
407 using BaseType::combinedOrder;
408 };
409
410 // BinaryGridFunctionExpression
411 // -------------------
412
416 template<class BinOp, class LeftFunction, class RightFunction>
418 : public GridFunctionExpression<typename
419 BinaryGridFunctionExpressionTraits<BinOp, LeftFunction, RightFunction>::FunctionSpaceType,
420 BinaryGridFunctionExpression<BinOp, LeftFunction, RightFunction> >,
421 public BinaryGridFunctionExpressionTraits<BinOp, LeftFunction, RightFunction>::BoundarySupportTagType,
422 public BinaryGridFunctionExpressionTraits<BinOp, LeftFunction, RightFunction>::PieceWiseConstantTagType
423 {
424 typedef LeftFunction LeftFunctionType;
425 typedef RightFunction RightFunctionType;
426 public:
429 private:
431 typedef Fem::Function<typename TraitsType::FunctionSpaceType, ThisType> BaseType;
432 typedef BinOp OperationTagType;
433 typedef BinaryFunctionExpressionOperation<BinOp> OperationType;
434
435 typedef typename TraitsType::LeftPieceWiseConstantType LeftPieceWiseConstantType;
436 typedef typename TraitsType::RightPieceWiseConstantType RightPieceWiseConstantType;
437 typedef typename TraitsType::PieceWiseConstantType PieceWiseConstantType;
438
439 public:
442
444 typedef typename TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
445
446 // type of discrete function space
447 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
448
450 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
451
453 typedef typename DiscreteFunctionSpaceType::GridType GridType;
454
455 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
456 static const int dimDomain = GridPartType::GridType::dimensionworld;
457 static const int dimLocal = GridPartType::GridType::dimension;
458
460 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
462 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
464 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
466 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
468 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
470 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
471
473 typedef typename TraitsType::EntityType EntityType;
474
476 typedef typename GridPartType::IntersectionType IntersectionType;
477
478 private:
479 class LocalFunction;
480 //class LocalFunctionStorage;
482 using ExpressionBaseType::expressionName_;
483
484 public:
487
489 BinaryGridFunctionExpression(const LeftFunctionType &f,
490 const RightFunctionType &g)
491 : space_(f.space().gridPart(), TraitsType::combinedOrder(f.space().order(), g.space().order())),
492 leftFunction_(f),
493 rightFunction_(g)
494 {}
495
496 // We DO allow copying
497 BinaryGridFunctionExpression(const ThisType &other)
498 : space_(other.space_),
499 leftFunction_(other.leftFunction_),
500 rightFunction_(other.rightFunction_)
501 {}
502
504 void evaluate(const DomainType &global, RangeType &result) const
505 {
506 OperationType::evaluate(leftFunction_(), rightFunction_(), global, result);
507 }
508
510 void jacobian(const DomainType &global, JacobianRangeType &result) const
511 {
512 if (PieceWiseConstantType()) {
513 result = 0.;
514 return;
515 }
516 OperationType::jacobian(leftFunction_(), rightFunction_(), global, result,
517 LeftPieceWiseConstantType(), RightPieceWiseConstantType());
518 }
519
521 void hessian(const DomainType &global, HessianRangeType &result) const
522 {
523 if (PieceWiseConstantType()) {
524 result = 0.;
525 return;
526 }
527 OperationType::hessian(leftFunction_(), rightFunction_(), global, result);
528 }
529
531 const LocalFunctionType localFunction(const EntityType &entity) const
532 {
533 return LocalFunctionType(entity, *this);
534 }
535
538 {
539 return LocalFunctionType(entity, *this);
540 }
541
543 const std::string &expressionName() const
544 {
545 expressionName_ = operationName(OperationTagType(), leftFunction().name(), rightFunction().name());
546 return expressionName_;
547 }
548
551 {
552 return space_;
553 }
554
555 const GridPartType &gridPart() const
556 {
557 return space().gridPart();
558 }
559
565 typedef typename TraitsType::IndicatorType IndicatorType;
566
567 enum {
571 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
572 };
573
576 const IntersectionType &intersection) const
577 {
578 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
579 return LocalFunctionType(entity, intersection, *this);
580 }
581
583 LocalFunctionType localFunction(const EntityType& entity, const IntersectionType &intersection)
584 {
585 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
586 return LocalFunctionType(entity, intersection, *this);
587 }
588
590 {
591 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
592 return TraitsType::indicator(leftFunction_(), rightFunction_());
593 }
594
596
597 public:
598 const LeftFunctionType& leftFunction() const { return leftFunction_(); }
599 const RightFunctionType& rightFunction() const { return rightFunction_(); }
600 protected:
602 ExpressionStorage<LeftFunctionType> leftFunction_;
603 ExpressionStorage<RightFunctionType> rightFunction_;
604 mutable std::string name_;
605 };
606
607 // BinaryGridFunctionExpression::LocalFunction
608 // ----------------------------------
609
611 template<class BinOp, class LeftFunctionType, class RightFunctionType>
612 class BinaryGridFunctionExpression<BinOp, LeftFunctionType, RightFunctionType>::LocalFunction
613 {
614 typedef LocalFunction ThisType;
615 typedef BinOp OperationTagType;
616 typedef BinaryFunctionExpressionOperation<BinOp> OperationType;
617 public:
620 private:
621 typedef typename LeftFunctionType::LocalFunctionType LeftLocalFunctionType;
622 typedef typename RightFunctionType::LocalFunctionType RightLocalFunctionType;
623 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
624 typedef typename EntityType::Geometry GeometryType;
625
626 public:
627
628 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
629 static const int dimDomain = GridPartType::GridType::dimensionworld;
630 static const int dimLocal = GridPartType::GridType::dimension;
631
633 typedef typename DiscreteFunctionSpaceType::FunctionSpaceType FunctionSpaceType;
634
636 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
638 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
640 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
642 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
644 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
646 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
647
650 : intersection_(0),
651 supported_(true),
652 df_(df),
653 leftLocalFunction_(df.leftFunction().localFunction(entity)),
654 rightLocalFunction_(df.rightFunction().localFunction(entity)),
655 entity_(&entity),
656 order_(df.space().order())
657 {}
658
660 : intersection_(0),
661 supported_(true),
662 df_(df),
663 leftLocalFunction_(df.leftFunction()),
664 rightLocalFunction_(df.rightFunction()),
665 entity_(0),
666 order_(df.space().order())
667 {}
668
670 template<class PointType>
671 void evaluate(const PointType &x, RangeType &ret) const
672 {
673 if (!supported()) {
674 ret = 0.;
675 return;
676 }
677 OperationType::evaluate(leftLocalFunction_, rightLocalFunction_, x, ret);
678 }
679
681 template<class PointType>
682 void jacobian(const PointType &x, JacobianRangeType &ret) const
683 {
684 if (DiscreteFunctionType::PieceWiseConstantType() || !supported()) {
685 ret = 0.;
686 return;
687 }
688 OperationType::jacobian(leftLocalFunction_, rightLocalFunction_, x, ret);
689 }
690
692 template<class PointType>
693 void hessian(const PointType &x, HessianRangeType &ret) const
694 {
695 if (DiscreteFunctionType::PieceWiseConstantType() || !supported()) {
696 ret = 0.;
697 return;
698 }
699 OperationType::hessian(leftLocalFunction_, rightLocalFunction_, x, ret);
700 }
701
703 template<class QuadratureType, class VectorType>
704 void evaluateQuadrature(const QuadratureType& quadrature, VectorType& values) const
705 {
706 assert(values.size() == quadrature.nop());
707 DiscreteFunctionType::evaluateQuadratureImp(*this, quadrature, values, values[0]);
708 }
709
711 int order() const { return order_; }
712
714 void init(const EntityType &entity)
715 {
716 leftLocalFunction_.init(entity);
717 rightLocalFunction_.init(entity);
718 entity_ = &entity;
719 supported_ = true;
720 }
721
723 const EntityType &entity() const
724 {
725 assert(entity_);
726 return *entity_;
727 }
728
733 enum {
738 };
739
742 const IntersectionType &intersection,
743 const DiscreteFunctionType &df)
744 : df_(df),
745 leftLocalFunction_(LeftLocalFunctionType(df.leftFunction())),
746 rightLocalFunction_(RightLocalFunctionType(df.rightFunction())),
747 order_(df.space().order())
748 {
749 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
750 init(entity, intersection);
751 }
752
754 void init(const EntityType& entity, const IntersectionType& intersection)
755 {
756 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
757
758 BoundaryFunctionTraits<LeftFunctionType>::init(leftLocalFunction_, entity, intersection);
759 BoundaryFunctionTraits<RightFunctionType>::init(rightLocalFunction_, entity, intersection);
760
761 supported_ = TraitsType::boundarySupported(BoundaryFunctionTraits<LeftFunctionType>::supported(leftLocalFunction_),
763
764 entity_ = &entity;
765 intersection_ = &intersection;
766 }
767
770 {
771 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
772 assert(intersection_);
773 return *intersection_;
774 }
775
776 bool supported() const
777 {
778 return globalSupport || (!emptySupport && supported_);
779 }
780
781 private:
782 const IntersectionType *intersection_;
783 bool supported_;
784
786
787 private:
788 const DiscreteFunctionType& df_;
789 LeftLocalFunctionType leftLocalFunction_;
790 RightLocalFunctionType rightLocalFunction_;
791 const EntityType *entity_;
792 int order_;
793 };
794
796 //
797 // Specialization for left- and right multiplication. Seems better
798 // than to construct first scalar function with constant values ...
799
801 template<class Factor, class Function>
803 : public BinaryFunctionExpressionTraits<SMultiplyOperation, Factor, Function>
804 {
805 typedef Factor FactorType;
806 typedef Function FunctionType;
807
808 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
809 typedef typename FunctionType::DiscreteFunctionSpaceType FunctionDiscreteFunctionSpaceType;
810
811 typedef typename FunctionDiscreteFunctionSpaceType::GridPartType GridPartType;
812
813 // ALWAYS use a DiscreteFunctionSpaceAdapter
814 typedef Fem::DiscreteFunctionSpaceAdapter<FunctionSpaceType, GridPartType> DiscreteFunctionSpaceType;
815
816 typedef typename GridPartType::GridType GridType;
817 typedef typename GridPartType::template Codim<0>::EntityType EntityType;
819 typedef typename GridPartType::template Codim<0>::IteratorType IteratorType;
821 typedef typename GridPartType::IndexSetType IndexSetType;
822
824
825 // p.w. constant tag
826 enum {
827 isPieceWiseConstant = std::is_base_of<IsPieceWiseConstant, Function>::value
828 };
829
830 // p.w. constant tag
831 typedef
832 typename std::conditional<isPieceWiseConstant,
834 PieceWiseConstantTagType;
835
836 static const unsigned int order =
837 isPieceWiseConstant ? 0 : DiscreteFunctionSpaceType::polynomialOrder;
838
839 // Boundary function stuff
840 typedef BoundaryFunctionTraits<Function> BoundaryTraits;
841 enum {
842 hasBoundarySupport = BoundaryTraits::hasBoundarySupport
843 };
844 typedef
845 typename std::conditional<hasBoundarySupport,
846 HasBoundarySupport, HasNoBoundarySupport>::type
847 BoundarySupportTagType;
848
849 typedef
850 typename BoundaryFunctionTraits<FunctionType>::IndicatorType
851 IndicatorType;
852
853 enum {
854 emptySupport = ExpressionTraits<IndicatorType>::isZero,
855 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
856 };
857 };
858
859 // BinaryGridFunctionExpression
860 // -------------------
861
863 template<class Factor, class Function>
865 : public GridFunctionExpression<typename Function::FunctionSpaceType,
866 BinaryGridFunctionExpression<SMultiplyOperation, Factor, Function> >,
867 public BinaryGridFunctionExpressionTraits<SMultiplyOperation, Factor, Function>::BoundarySupportTagType,
868 public BinaryGridFunctionExpressionTraits<SMultiplyOperation, Factor, Function>::PieceWiseConstantTagType
869 {
871 typedef Fem::Function<typename Function::FunctionSpaceType, ThisType> BaseType;
874
875 // Make sure the functions are discrete functions
876 static_assert((std::is_base_of<Fem::HasLocalFunction, Function>::value),
877 "FunctionType must be a discrete function type.");
878
879 public:
881 typedef Factor FactorType;
882
884 typedef Function ContainedFunctionType;
885
888
890 enum {
891 isPieceWiseConstant = TraitsType::isPieceWiseConstant
892 };
893
896
898 typedef typename TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
899
900 // type of discrete function space
901 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
902
904 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
905
907 typedef typename DiscreteFunctionSpaceType::GridType GridType;
908
910 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
912 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
914 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
916 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
918 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
920 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
921
923 typedef typename TraitsType::EntityType EntityType;
924 typedef typename GridPartType::IntersectionType IntersectionType;
925
926 private:
927 class LocalFunction;
929 using ExpressionBaseType::expressionName_;
930
931 public:
933 typedef LocalFunction LocalFunctionType;
934
935 // reference to function this local belongs to
937 : space_(f.space().gridPart(), f.space().order()),
938 scalar_(s),
939 function_(f)
940 {}
941
942 // We DO allow copying of expressions
943 BinaryGridFunctionExpression(const ThisType &other)
944 : space_(other.space_),
945 scalar_(other.scalar_),
946 function_(other.function_)
947 {}
948
950 void evaluate(const DomainType &global, RangeType &result) const
951 {
952 function().evaluate(global, result);
953 result *= parameterValue(scalar());
954 }
955
957 void jacobian(const DomainType &global, JacobianRangeType &result) const
958 {
959 if (isPieceWiseConstant) {
960 result = 0.;
961 return;
962 }
963 function().jacobian(global, result);
964 result *= parameterValue(scalar());
965 }
966
968 void hessian(const DomainType &global, HessianRangeType &result) const
969 {
970 if (isPieceWiseConstant) {
971 result = 0.;
972 return;
973 }
974 function().hessian(global, result);
975 for (int k = 0; k < HessianRangeType::dimension; ++k) {
976 result[k] *= parameterValue(scalar());
977 }
978 }
979
981 const LocalFunctionType localFunction(const EntityType &entity) const
982 {
983 return LocalFunctionType(entity, *this);
984 }
985
988 {
989 return LocalFunctionType(entity, *this);
990 }
991
993 const std::string &expressionName() const
994 {
995 std::string factorName = std::to_string(parameterValue(scalar()));
996 if (ParameterValue<FactorType>::isParameter) {
997 factorName = "P(" + factorName + ")";
998 }
999 expressionName_ = operationName(OperationTagType(), factorName, function().name());
1000 return expressionName_;
1001 }
1002
1004 const DiscreteFunctionSpaceType &space() const { return space_; }
1005
1006 const GridPartType &gridPart() const { return space().gridPart(); }
1007
1014
1015 enum {
1019 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
1021
1023 const IntersectionType &intersection) const
1024 {
1025 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1026 return LocalFunctionType(entity, intersection, *this);
1027 }
1028
1030 const IntersectionType &intersection)
1031 {
1032 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1033 return LocalFunctionType(entity, intersection, *this);
1034 }
1035
1037 {
1038 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1040 }
1041
1043
1044 public:
1045 const FactorType& scalar() const { return scalar_(); }
1046 const ContainedFunctionType& function() const { return function_(); }
1047 protected:
1049 ExpressionStorage<FactorType> scalar_;
1050 ExpressionStorage<ContainedFunctionType> function_;
1051 mutable std::string name_;
1052 };
1053
1055 template<class FactorType, class FunctionType>
1056 class BinaryGridFunctionExpression<SMultiplyOperation, FactorType, FunctionType>::LocalFunction
1057 {
1058 typedef LocalFunction ThisType;
1059 typedef SMultiplyOperation OperationTagType;
1060 typedef BinaryFunctionExpressionOperation<SMultiplyOperation> OperationType;
1061 public:
1062 typedef BinaryGridFunctionExpression<SMultiplyOperation, FactorType, FunctionType> DiscreteFunctionType;
1064 typedef typename DiscreteFunctionType::FunctionSpaceType FunctionSpaceType;
1065 private:
1066 typedef typename FunctionType::LocalFunctionType SubLocalFunctionType;
1067 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1068 typedef typename EntityType::Geometry GeometryType;
1069 public:
1070
1071 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
1072 static const int dimDomain = GridPartType::GridType::dimensionworld;
1073 static const int dimLocal = GridPartType::GridType::dimension;
1074
1076 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
1078 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
1080 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
1082 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
1084 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
1086 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
1087
1090 : intersection_(0),
1091 supported_(true),
1092 scalar_(df.scalar()),
1093 localFunction_(df.function().localFunction(entity)),
1094 entity_(&entity),
1095 order_(df.space().order())
1096 {}
1097
1099 : intersection_(0),
1100 supported_(true),
1101 scalar_(df.scalar()),
1102 localFunction_(df.function()),
1103 entity_(0),
1104 order_(df.space().order())
1105 {}
1106
1108 template<class PointType>
1109 void evaluate(const PointType &x, RangeType &ret) const
1110 {
1111 if (!supported()) {
1112 ret = 0.;
1113 return;
1114 }
1115 localFunction_.evaluate(x, ret);
1116 ret *= parameterValue(scalar_);
1117 }
1118
1120 template<class PointType>
1121 void jacobian(const PointType &x, JacobianRangeType &ret) const
1122 {
1123 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1124 ret = 0.;
1125 return;
1126 }
1127 localFunction_.jacobian(x, ret);
1128 ret *= parameterValue(scalar_);
1129 }
1130
1132 template<class PointType>
1133 void hessian(const PointType &x, HessianRangeType &ret) const
1134 {
1135 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1136 ret = 0.;
1137 return;
1138 }
1139 localFunction_.hessian(x, ret);
1140 for (int k = 0; k < HessianRangeType::dimension; ++k) {
1141 ret[k] *= parameterValue(scalar_);
1142 }
1143 }
1144
1146 template<class QuadratureType, class VectorType>
1147 void evaluateQuadrature(const QuadratureType& quadrature, VectorType& values) const
1148 {
1149 assert(values.size() == quadrature.nop());
1150 DiscreteFunctionType::evaluateQuadratureImp(*this, quadrature, values, values[0]);
1151 }
1152
1154 int order() const
1155 {
1156 return DiscreteFunctionType::isPieceWiseConstant ? 0 : order_;
1157 }
1158
1160 void init(const EntityType &entity)
1161 {
1162 // scalar_ is as a constant unaffected from the entity
1163 localFunction_.init(entity);
1164 entity_ = &entity;
1165 supported_ = true;
1166 }
1167
1168 const EntityType &entity() const
1169 {
1170 assert(entity_);
1171 return *entity_;
1172 }
1173
1178 enum {
1183 };
1184
1186 const IntersectionType &intersection, const DiscreteFunctionType &df)
1187 : scalar_(df.scalar()),
1188 localFunction_(LocalFunctiontype(df.function())),
1189 order_(df.space().order())
1190 {
1191 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1192 init(entity, intersection);
1193 }
1194
1196 void init(const EntityType& entity, const IntersectionType& intersection)
1197 {
1198 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1199
1200 BoundaryFunctionTraits<FunctionType>::init(localFunction_, entity, intersection);
1201 supported_ = BoundaryFunctionTraits<FunctionType>::supported(localFunction_);
1202
1203 entity_ = &entity;
1204 intersection_ = &intersection;
1205 }
1206
1209 {
1210 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1211 assert(intersection_);
1212 return *intersection_;
1213 }
1214
1215 bool supported() const
1216 {
1217 return globalSupport || (!emptySupport && supported_);
1218 }
1219
1221
1222 private:
1223 const IntersectionType *intersection_;
1224 bool supported_;
1225
1226 private:
1227 const FactorType& scalar_;
1228 SubLocalFunctionType localFunction_;
1229 const EntityType *entity_;
1230 int order_;
1231 };
1232
1234 //
1235 // Unary expressions
1236
1237 // The result space may depend on the operation. The default is to
1238 // have the result space equal to the function space (inversion,
1239 // negation etc.)
1240 template<class UnOp, class FunctionType>
1241 struct UnaryGridFunctionExpressionBaseTraits
1242 {
1243 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
1244 typedef typename FunctionType::DiscreteFunctionSpaceType FunctionDiscreteFunctionSpaceType;
1245 typedef typename FunctionDiscreteFunctionSpaceType::GridPartType GridPartType;
1246
1247 typedef Fem::DiscreteFunctionSpaceAdapter<FunctionSpaceType, GridPartType> DiscreteFunctionSpaceType;
1248
1249 static const unsigned int order =
1251 FunctionDiscreteFunctionSpaceType::polynomialOrder);
1252 };
1253
1254 template<class FunctionType>
1255 struct UnaryGridFunctionExpressionBaseTraits<SquareOperation, FunctionType>
1256 {
1257 typedef typename FunctionType::DiscreteFunctionSpaceType::ScalarFunctionSpaceType FunctionSpaceType;
1258 typedef typename FunctionType::DiscreteFunctionSpaceType FunctionDiscreteFunctionSpaceType;
1259 typedef typename FunctionDiscreteFunctionSpaceType::GridPartType GridPartType;
1260
1261 typedef Fem::DiscreteFunctionSpaceAdapter<FunctionSpaceType, GridPartType> DiscreteFunctionSpaceType;
1262
1263 static const unsigned int order = 2*FunctionDiscreteFunctionSpaceType::polynomialOrder;
1264 };
1265
1267 template <class UnOp, class Function>
1269 : public UnaryGridFunctionExpressionBaseTraits<UnOp, Function>
1270 {
1271 typedef UnaryGridFunctionExpressionBaseTraits<UnOp, Function> BaseType;
1272 typedef Function FunctionType;
1273 typedef typename BaseType::FunctionSpaceType FunctionSpaceType;
1274
1275 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
1276 typedef typename FunctionSpaceType::DomainFieldType DomainFieldType;
1277 typedef typename FunctionSpaceType::RangeType RangeType;
1278 typedef typename FunctionSpaceType::DomainType DomainType;
1279 typedef typename FunctionSpaceType::JacobianRangeType JacobianRangeType;
1280
1281 typedef typename BaseType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
1282
1283 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1284 typedef typename GridPartType::GridType GridType;
1285 typedef typename GridPartType::template Codim<0>::EntityType EntityType;
1287 typedef typename GridPartType::template Codim<0>::IteratorType IteratorType;
1289 typedef typename GridPartType::IndexSetType IndexSetType;
1290
1292
1293 // p.w. constant tag
1294 enum {
1295 isPieceWiseConstant = std::is_base_of<IsPieceWiseConstant, Function>::value
1296 };
1297
1298 typedef
1299 typename std::conditional<isPieceWiseConstant,
1301 PieceWiseConstantTagType;
1302
1303 static const unsigned int order = isPieceWiseConstant ? 0 : BaseType::order;
1304
1305 // Boundary function stuff
1306 typedef BoundaryFunctionTraits<Function> BoundaryTraits;
1307 enum {
1308 hasBoundarySupport = BoundaryTraits::hasBoundarySupport
1309 };
1310 typedef
1311 typename std::conditional<hasBoundarySupport,
1312 HasBoundarySupport, HasNoBoundarySupport>::type
1313 BoundarySupportTagType;
1314 };
1315
1316 // UnaryGridFunctionExpression
1317 // -------------------
1318
1322 template<class UnOp, class FunctionType>
1324 : public GridFunctionExpression<typename UnaryGridFunctionExpressionTraits<UnOp, FunctionType>::FunctionSpaceType,
1325 UnaryGridFunctionExpression<UnOp, FunctionType> >,
1326 public UnaryGridFunctionExpressionTraits<UnOp, FunctionType>::BoundarySupportTagType,
1327 public UnaryGridFunctionExpressionTraits<UnOp, FunctionType>::PieceWiseConstantTagType
1328 {
1329 public:
1332 private:
1334 typedef Fem::Function<typename TraitsType::FunctionSpaceType, ThisType> BaseType;
1335 typedef UnOp OperationTagType;
1336 typedef UnaryFunctionExpressionOperation<UnOp> OperationType;
1337
1338 // Make sure the functions are discrete functions
1339 static_assert((std::is_base_of<Fem::HasLocalFunction, FunctionType>::value),
1340 "FunctionType must have a LocalFunctionType type.");
1341
1342 enum {
1343 isPieceWiseConstant = TraitsType::isPieceWiseConstant
1344 };
1345
1346 public:
1349
1351 typedef typename TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
1352
1353 // type of discrete function space
1354 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
1355
1357 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1358
1360 typedef typename DiscreteFunctionSpaceType::GridType GridType;
1361
1363 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
1365 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
1367 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
1369 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
1371 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
1373 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
1374
1376 typedef typename TraitsType::EntityType EntityType;
1377 typedef typename GridPartType::IntersectionType IntersectionType;
1378
1379 private:
1380 class LocalFunction;
1382 using ExpressionBaseType::expressionName_;
1383
1384 public:
1386 typedef LocalFunction LocalFunctionType;
1387
1388 // reference to function this local belongs to
1389 UnaryGridFunctionExpression(const FunctionType &f)
1390 : space_(f.space().gridPart(), TraitsType::order),
1391 function_(f)
1392 {}
1393
1394 // reference to function this local belongs to
1395 UnaryGridFunctionExpression(const ThisType &other)
1396 : space_(other.space_),
1397 function_(other.function_)
1398 {}
1399
1401 void evaluate(const DomainType &global, RangeType &result) const
1402 {
1403 OperationType::evaluate(function(), global, result);
1404 }
1405
1407 void jacobian(const DomainType &global, JacobianRangeType &result) const
1408 {
1409 if (isPieceWiseConstant) {
1410 result = 0.;
1411 return;
1412 }
1413 OperationType::jacobian(function(), global, result);
1414 }
1415
1417 void hessian(const DomainType &global, HessianRangeType &result) const
1418 {
1419 if (isPieceWiseConstant) {
1420 result = 0.;
1421 return;
1422 }
1423 OperationType::hessian(function(), global, result);
1424 }
1425
1427 const LocalFunctionType localFunction(const EntityType &entity) const
1428 {
1429 return LocalFunctionType(entity, *this);
1430 }
1431
1434 {
1435 return LocalFunctionType(entity, *this);
1436 }
1437
1439 const std::string &expressionName() const
1440 {
1441 expressionName_ = operationName(OperationTagType(), function().name());
1442 return expressionName_;
1443 }
1444
1446 const DiscreteFunctionSpaceType &space() const { return space_; }
1447
1448 const GridPartType &gridPart() const { return space().gridPart(); }
1449
1456
1457 enum {
1461 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
1463
1466 const IntersectionType &intersection) const
1467 {
1468 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1469 return LocalFunctionType(entity, intersection, *this);
1470 }
1471
1474 const IntersectionType &intersection)
1475 {
1476 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1477 return LocalFunctionType(entity, intersection, *this);
1478 }
1479
1481 {
1482 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1484 }
1485
1487
1488 public:
1489 const FunctionType& function() const { return function_(); }
1490 private:
1492 ExpressionStorage<FunctionType> function_;
1493 mutable std::string name_;
1494 };
1495
1496 // UnaryGridFunctionExpression::LocalFunction
1497 // ----------------------------------
1498
1499 template<class UnOp, class FunctionType>
1500 class UnaryGridFunctionExpression<UnOp, FunctionType>::LocalFunction
1501 {
1502 typedef LocalFunction ThisType;
1503 typedef UnOp OperationTagType;
1504 typedef UnaryFunctionExpressionOperation<UnOp> OperationType;
1505 public:
1506 typedef UnaryGridFunctionExpression<UnOp, FunctionType> DiscreteFunctionType;
1508 typedef typename DiscreteFunctionType::FunctionSpaceType FunctionSpaceType;
1509 private:
1510 typedef typename FunctionType::LocalFunctionType LocalFunctionType;
1511 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1512 typedef typename EntityType::Geometry GeometryType;
1513 public:
1514
1515 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
1516 static const int dimDomain = GridPartType::GridType::dimensionworld;
1517 static const int dimLocal = GridPartType::GridType::dimension;
1518
1520 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
1522 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
1524 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
1526 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
1528 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
1530 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
1531
1533 LocalFunction(const EntityType &entity, const DiscreteFunctionType &df)
1534 : intersection_(0),
1535 supported_(true),
1536 localFunction_(df.function().localFunction(entity)),
1537 order_(df.space().order())
1538 {}
1539
1540 LocalFunction(const DiscreteFunctionType &df)
1541 : intersection_(0),
1542 supported_(true),
1543 localFunction_(df.function()),
1544 order_(df.space().order())
1545 {}
1546
1548 template<class PointType>
1549 void evaluate(const PointType &x, RangeType &ret) const
1550 {
1551 if (!supported()) {
1552 ret = 0.;
1553 return;
1554 }
1555 OperationType::evaluate(localFunction_, x, ret);
1556 }
1557
1559 template<class PointType>
1560 void jacobian(const PointType &x, JacobianRangeType &ret) const
1561 {
1562 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1563 ret = 0.;
1564 return;
1565 }
1566 OperationType::jacobian(localFunction_, x, ret);
1567 }
1568
1570 template<class PointType>
1571 void hessian(const PointType &x, HessianRangeType &ret) const
1572 {
1573 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1574 ret = 0.;
1575 return;
1576 }
1577 OperationType::hessian(localFunction_, x, ret);
1578 }
1579
1581 template<class QuadratureType, class VectorType>
1582 void evaluateQuadrature(const QuadratureType& quadrature, VectorType& values) const
1583 {
1584 assert(values.size() == quadrature.nop());
1585 DiscreteFunctionType::evaluateQuadratureImp(*this, quadrature, values, values[0]);
1586 }
1587
1589 int order() const { return order_; }
1590
1592 void init(const EntityType &entity)
1593 {
1594 localFunction_.init(entity);
1595 }
1596
1597 const EntityType &entity() const
1598 {
1599 localFunction_.entity();
1600 }
1601
1606 enum {
1611 };
1612
1613 LocalFunction(const EntityType& entity,
1614 const IntersectionType &intersection,
1615 const DiscreteFunctionType &df)
1616 : localFunction_(LocalFunctiontype(df.function())),
1617 order_(df.space().order())
1618 {
1619 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1620 init(entity, intersection);
1621 }
1622
1624 void init(const EntityType& entity, const IntersectionType& intersection)
1625 {
1626 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1627
1628 BoundaryFunctionTraits<FunctionType>::init(localFunction_, entity, intersection);
1629 supported_ = BoundaryFunctionTraits<FunctionType>::supported(localFunction_);
1630
1631 intersection_ = &intersection;
1632 }
1633
1635 const IntersectionType &intersection() const
1636 {
1637 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1638 assert(intersection_);
1639 return *intersection_;
1640 }
1641
1642 bool supported() const
1643 {
1644 return globalSupport || (!emptySupport && supported_);
1645 }
1646
1647 private:
1648 const IntersectionType *intersection_;
1649 bool supported_;
1650
1652
1653 private:
1654 LocalFunctionType localFunction_;
1655 int order_;
1656 };
1657
1658 // UnaryGridFunctionExpression<IdentityOperation>
1659 // -------------------
1660
1662 template<class FunctionType>
1664 : public GridFunctionExpression<typename FunctionType::FunctionSpaceType,
1665 UnaryGridFunctionExpression<IdentityOperation, FunctionType> >,
1666 public UnaryGridFunctionExpressionTraits<IdentityOperation, FunctionType>::BoundarySupportTagType,
1667 public UnaryGridFunctionExpressionTraits<IdentityOperation, FunctionType>::PieceWiseConstantTagType
1668 {
1670 typedef Fem::Function<typename FunctionType::FunctionSpaceType, ThisType> BaseType;
1672 typedef UnaryFunctionExpressionOperation<OperationTagType> OperationType;
1673
1674 // Make sure the functions are discrete functions
1675 static_assert((std::is_base_of<Fem::HasLocalFunction, FunctionType>::value),
1676 "FunctionType must be a discrete function type.");
1677
1678 public:
1681
1682 enum {
1683 isPieceWiseConstant = TraitsType::isPieceWiseConstant
1684 };
1685
1686 typedef ThisType DiscreteFunctionType;
1687
1689 typedef typename TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
1690
1691 // type of discrete function space
1692 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
1693
1695 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1696
1698 typedef typename DiscreteFunctionSpaceType::GridType GridType;
1699
1701 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
1703 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
1705 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
1707 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
1709 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
1711 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
1712
1714 typedef typename TraitsType::EntityType EntityType;
1715 typedef typename GridPartType::IntersectionType IntersectionType;
1716
1717 private:
1718 class LocalFunction;
1720 using ExpressionBaseType::expressionName_;
1721
1722 public:
1724 typedef LocalFunction LocalFunctionType;
1725
1726 // reference to function this local belongs to
1727 UnaryGridFunctionExpression(const FunctionType &f)
1728 : space_(f.space().gridPart(), TraitsType::order),
1729 function_(f)
1730 {}
1731
1732 // reference to function this local belongs to
1733 UnaryGridFunctionExpression(const ThisType &other)
1734 : space_(other.space_),
1735 function_(other.function_)
1736 {}
1737
1739 void evaluate(const DomainType &global, RangeType &result) const
1740 {
1741 function().evaluate(global, result);
1742 }
1743
1745 void jacobian(const DomainType &global, JacobianRangeType &result) const
1746 {
1747 if (isPieceWiseConstant) {
1748 result = 0.;
1749 return;
1750 }
1751 function().jacobian(global, result);
1752 }
1753
1755 void hessian(const DomainType &global, HessianRangeType &result) const
1756 {
1757 if (isPieceWiseConstant) {
1758 result = 0.;
1759 return;
1760 }
1761 function().hessian(global, result);
1762 }
1763
1765 const LocalFunctionType localFunction(const EntityType &entity) const
1766 {
1767 return LocalFunctionType(entity, *this);
1768 }
1769
1772 {
1773 return LocalFunctionType(entity, *this);
1774 }
1775
1777 const std::string &expressionName() const
1778 {
1779 expressionName_ = operationName(OperationTagType(), function().name());
1780 return expressionName_;
1781 }
1782
1784 const DiscreteFunctionSpaceType &space() const { return space_; }
1785
1786 const GridPartType &gridPart() const { return space().gridPart(); }
1787
1794
1795 enum {
1799 globalSupport = std::is_same<IndicatorType, EntireBoundaryIndicatorType>::value
1801
1804 const IntersectionType &intersection) const
1805 {
1806 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1807 return LocalFunctionType(entity, intersection, *this);
1808 }
1809
1812 const IntersectionType &intersection)
1813 {
1814 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1815 return LocalFunctionType(entity, intersection, *this);
1816 }
1817
1819 {
1820 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1822 }
1823
1825
1826 private:
1827 const FunctionType& function() const { return function_(); }
1828 private:
1830 ExpressionStorage<FunctionType> function_;
1831 mutable std::string name_;
1832 };
1833
1834 // UnaryGridFunctionExpression::LocalFunction
1835 // ----------------------------------
1836 template<class FunctionType>
1837 class UnaryGridFunctionExpression<IdentityOperation, FunctionType>::LocalFunction
1838 {
1839 typedef LocalFunction ThisType;
1840 typedef IdentityOperation OperationTagType;
1841 typedef UnaryFunctionExpressionOperation<OperationTagType> OperationType;
1842 public:
1843 typedef UnaryGridFunctionExpression<IdentityOperation, FunctionType> DiscreteFunctionType;
1845 typedef typename DiscreteFunctionType::FunctionSpaceType FunctionSpaceType;
1846 private:
1847 typedef typename FunctionType::LocalFunctionType LocalFunctionType;
1848 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
1849 typedef typename EntityType::Geometry GeometryType;
1850 public:
1851
1852 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
1853 static const int dimDomain = GridPartType::GridType::dimensionworld;
1854 static const int dimLocal = GridPartType::GridType::dimension;
1855
1857 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
1859 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
1861 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
1863 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
1865 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
1867 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
1868
1870 LocalFunction(const EntityType &entity, const DiscreteFunctionType &df)
1871 : intersection_(0),
1872 supported_(true),
1873 localFunction_(df.function().localFunction(entity))
1874 {}
1875
1876 LocalFunction(const DiscreteFunctionType &df)
1877 : intersection_(0),
1878 supported_(true),
1879 localFunction_(df.function())
1880 {}
1881
1883 template<class PointType>
1884 void evaluate(const PointType &x, RangeType &ret) const
1885 {
1886 if (!supported()) {
1887 ret = 0.;
1888 return;
1889 }
1890 localFunction_.evaluate(x, ret);
1891 }
1892
1894 template<class PointType>
1895 void jacobian(const PointType &x, JacobianRangeType &ret) const
1896 {
1897 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1898 ret = 0.;
1899 return;
1900 }
1901 localFunction_.jacobian(x, ret);
1902 }
1903
1905 template<class PointType>
1906 void hessian(const PointType &x, HessianRangeType &ret) const
1907 {
1908 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
1909 ret = 0.;
1910 return;
1911 }
1912 localFunction_.hessian(x, ret);
1913 }
1914
1916 template<class QuadratureType, class VectorType>
1917 void evaluateQuadrature(const QuadratureType& quadrature, VectorType& values) const
1918 {
1919 assert(values.size() == quadrature.nop());
1920 DiscreteFunctionType::evaluateQuadratureImp(*this, quadrature, values, values[0]);
1921 }
1922
1924 int order() const { return localFunction_.order(); }
1925
1927 void init(const EntityType &entity)
1928 {
1929 localFunction_.init(entity);
1930 supported_ = true;
1931 }
1932
1933 const EntityType &entity() const
1934 {
1935 return localFunction_.entity();
1936 }
1937
1942 enum {
1947 };
1948
1949 LocalFunction(const EntityType& entity,
1950 const IntersectionType &intersection,
1951 const DiscreteFunctionType &df)
1952 : localFunction_(LocalFunctiontype(df.function()))
1953 {
1954 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1955 init(entity, intersection);
1956 }
1957
1959 void init(const EntityType& entity, const IntersectionType& intersection)
1960 {
1961 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1962
1963 BoundaryFunctionTraits<FunctionType>::init(localFunction_, entity, intersection);
1964 supported_ = BoundaryFunctionTraits<FunctionType>::supported(localFunction_);
1965
1966 intersection_ = &intersection;
1967 }
1968
1970 const IntersectionType &intersection() const
1971 {
1972 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
1973 assert(intersection_);
1974 return *intersection_;
1975 }
1976
1977 bool supported() const
1978 {
1979 return globalSupport || (!emptySupport && supported_);
1980 }
1981
1982 private:
1983 const IntersectionType *intersection_;
1984 bool supported_;
1985
1987
1988 private:
1989 LocalFunctionType localFunction_;
1990 };
1991
1992 // UnaryGridFunctionExpression<InvertOperation>
1993 // -------------------
1994
2003 template<class FunctionType>
2005 : public GridFunctionExpression<typename FunctionType::FunctionSpaceType,
2006 UnaryGridFunctionExpression<InvertOperation, FunctionType> >,
2007 public UnaryGridFunctionExpressionTraits<InvertOperation, FunctionType>::BoundarySupportTagType
2008 {
2010 typedef Fem::Function<typename FunctionType::FunctionSpaceType, ThisType> BaseType;
2012 typedef UnaryFunctionExpressionOperation<OperationTagType> OperationType;
2013
2014 // Make sure the functions are discrete functions
2015 static_assert((std::is_base_of<Fem::HasLocalFunction, FunctionType>::value),
2016 "FunctionType must be a discrete function type.");
2017
2018 public:
2021
2022 enum {
2023 isPieceWiseConstant = TraitsType::isPieceWiseConstant
2024 };
2025
2026 typedef ThisType DiscreteFunctionType;
2027
2029 typedef typename TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
2030
2031 // type of discrete function space
2032 typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
2033
2035 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
2036
2038 typedef typename DiscreteFunctionSpaceType::GridType GridType;
2039
2041 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
2043 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
2045 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
2047 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
2049 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
2051 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
2052
2054 typedef typename TraitsType::EntityType EntityType;
2055 typedef typename GridPartType::IntersectionType IntersectionType;
2056
2057 private:
2058 class LocalFunction;
2060 using ExpressionBaseType::expressionName_;
2061
2062 public:
2064 typedef LocalFunction LocalFunctionType;
2065
2066 // reference to function this local belongs to
2067 UnaryGridFunctionExpression(const FunctionType &f)
2068 : space_(f.space().gridPart(), TraitsType::order),
2069 function_(f)
2070 {}
2071
2072 // reference to function this local belongs to
2073 UnaryGridFunctionExpression(const ThisType &other)
2074 : space_(other.space_),
2075 function_(other.function_)
2076 {}
2077
2079 void evaluate(const DomainType &global, RangeType &result) const
2080 {
2081 OperationType::evaluate(function(), global, result);
2082 }
2083
2085 void jacobian(const DomainType &global, JacobianRangeType &result) const
2086 {
2087 if (isPieceWiseConstant) {
2088 result = 0.;
2089 return;
2090 }
2091 OperationType::jacobian(function(), global, result);
2092 }
2093
2095 void hessian(const DomainType &global, HessianRangeType &result) const
2096 {
2097 if (isPieceWiseConstant) {
2098 result = 0.;
2099 return;
2100 }
2101 OperationType::hessian(function(), global, result);
2102 }
2103
2105 const LocalFunctionType localFunction(const EntityType &entity) const
2106 {
2107 return LocalFunctionType(entity, *this);
2108 }
2109
2112 {
2113 return LocalFunctionType(entity, *this);
2114 }
2115
2117 const std::string &expressionName() const
2118 {
2119 expressionName_ = operationName(OperationTagType(), function().name());
2120 return expressionName_;
2121 }
2122
2124 const DiscreteFunctionSpaceType &space() const { return space_; }
2125
2126 const GridPartType &gridPart() const { return space().gridPart(); }
2127
2133 typedef EntireBoundaryIndicatorType IndicatorType;
2134
2135 enum {
2139 globalSupport = true
2141
2144 const IntersectionType &intersection) const
2145 {
2146 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2147 return LocalFunctionType(entity, intersection, *this);
2148 }
2149
2152 const IntersectionType &intersection)
2153 {
2154 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2155 return LocalFunctionType(entity, intersection, *this);
2156 }
2157
2159 {
2160 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2161 return IndicatorType();
2162 }
2163
2165
2166 private:
2167 const FunctionType& function() const { return function_(); }
2168 private:
2170 ExpressionStorage<FunctionType> function_;
2171 mutable std::string name_;
2172 };
2173
2174 // UnaryGridFunctionExpression::LocalFunction
2175 // ----------------------------------
2176 template<class FunctionType>
2177 class UnaryGridFunctionExpression<InvertOperation, FunctionType>::LocalFunction
2178 {
2179 typedef LocalFunction ThisType;
2180 public:
2181 typedef UnaryGridFunctionExpression<InvertOperation, FunctionType> DiscreteFunctionType;
2183 typedef typename DiscreteFunctionType::FunctionSpaceType FunctionSpaceType;
2184 private:
2185 typedef InvertOperation OperationTagType;
2186 typedef UnaryFunctionExpressionOperation<OperationTagType> OperationType;
2187 private:
2188 typedef typename FunctionType::LocalFunctionType LocalFunctionType;
2189 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
2190 typedef typename EntityType::Geometry GeometryType;
2191 public:
2192
2193 static const int dimRange = DiscreteFunctionSpaceType::dimRange;
2194 static const int dimDomain = GridPartType::GridType::dimensionworld;
2195 static const int dimLocal = GridPartType::GridType::dimension;
2196
2198 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
2200 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
2202 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
2204 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
2206 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
2208 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
2209
2211 LocalFunction(const EntityType &entity, const DiscreteFunctionType &df)
2212 : intersection_(0),
2213 supported_(true),
2214 localFunction_(df.function().localFunction(entity))
2215 {}
2216
2217 LocalFunction(const DiscreteFunctionType &df)
2218 : intersection_(0),
2219 supported_(true),
2220 localFunction_(df.function()) {}
2221
2223 template<class PointType>
2224 void evaluate(const PointType &x, RangeType &ret) const
2225 {
2226 if (!supported()) {
2227 ret = 0.;
2228 return;
2229 }
2230 OperationType::evaluate(localFunction_, x, ret);
2231 }
2232
2234 template<class PointType>
2235 void jacobian(const PointType &x, JacobianRangeType &ret) const
2236 {
2237 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
2238 ret = 0.;
2239 return;
2240 }
2241 OperationType::jacobian(localFunction_, x, ret);
2242 }
2243
2245 template<class PointType>
2246 void hessian(const PointType &x, HessianRangeType &ret) const
2247 {
2248 if (DiscreteFunctionType::isPieceWiseConstant || !supported()) {
2249 ret = 0.;
2250 return;
2251 }
2252 OperationType::hessian(localFunction_, x, ret);
2253 }
2254
2256 template<class QuadratureType, class VectorType>
2257 void evaluateQuadrature(const QuadratureType& quadrature, VectorType& values) const
2258 {
2259 assert(values.size() == quadrature.nop());
2260 DiscreteFunctionType::evaluateQuadratureImp(*this, quadrature, values, values[0]);
2261 }
2262
2264 int order() const
2265 {
2266 return DiscreteFunctionType::isPieceWiseConstant ? 0 : 111;
2267 }
2268
2270 void init(const EntityType &entity)
2271 {
2272 localFunction_.init(entity);
2273 supported_ = true;
2274 }
2275
2276 const EntityType &entity() const
2277 {
2278 return localFunction_.entity();
2279 }
2280
2285 enum {
2290 };
2291
2292 LocalFunction(const EntityType& entity,
2293 const IntersectionType &intersection,
2294 const DiscreteFunctionType &df)
2295 : localFunction_(LocalFunctiontype(df.function()))
2296 {
2297 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2298 init(entity, intersection);
2299 }
2300
2302 void init(const EntityType& entity, const IntersectionType& intersection)
2303 {
2304 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2305
2306 BoundaryFunctionTraits<FunctionType>::init(localFunction_, entity, intersection);
2307 supported_ = BoundaryFunctionTraits<FunctionType>::supported(localFunction_);
2308
2309 intersection_ = &intersection;
2310 }
2311
2313 const IntersectionType &intersection() const
2314 {
2315 static_assert(TraitsType::hasBoundarySupport, "This is not a BoundarySupportedFunction");
2316 assert(intersection_);
2317 return *intersection_;
2318 }
2319
2320 bool supported() const
2321 {
2322 return globalSupport || (!emptySupport && supported_);
2323 }
2324
2325 private:
2326 const IntersectionType *intersection_;
2327 bool supported_;
2328
2330
2331 private:
2332 LocalFunctionType localFunction_;
2333 };
2334
2336 //
2337 // Multiplication
2338
2339 template<class LeftFunction, class RightFunction>
2340 static inline
2341 BinaryGridFunctionExpression<MultiplyOperation, LeftFunction, RightFunction>
2342 operator*(const Fem::Function<typename LeftFunction::FunctionSpaceType, LeftFunction>& f_,
2343 const Fem::Function<typename RightFunction::FunctionSpaceType, RightFunction>& g_)
2344 {
2345 typedef LeftFunction LeftFunctionType;
2346 typedef RightFunction RightFunctionType;
2347
2348 static_assert(std::is_base_of<Fem::HasLocalFunction, LeftFunctionType>::value,
2349 "Left operand must be a GridFunction");
2350 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2351 "Right operand must be a GridFunction");
2352
2353 const LeftFunctionType& f(static_cast<const LeftFunctionType&>(f_));
2354 const RightFunctionType& g(static_cast<const RightFunctionType&>(g_));
2355
2356 typedef BinaryGridFunctionExpression<MultiplyOperation, LeftFunctionType, RightFunctionType> ExpressionType;
2357 return ExpressionType(f, g);
2358 }
2359
2360 // Multiplication by scalars
2361
2362 // s * f
2363 template<class Function>
2364 static inline
2365 BinaryGridFunctionExpression<SMultiplyOperation, typename Function::RangeFieldType, Function>
2366 operator*(const typename Function::RangeFieldType& s,
2367 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2368 {
2369 typedef Function FunctionType;
2370
2371 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2372 "Right operand must be a GridFunction");
2373
2374 const FunctionType& f(static_cast<const FunctionType&>(f_));
2375 typedef typename FunctionType::RangeFieldType RangeFieldType;
2376 typedef BinaryGridFunctionExpression<SMultiplyOperation, RangeFieldType, FunctionType> ExpressionType;
2377
2378 return ExpressionType(s, f);
2379 }
2380
2382 template<class Function>
2383 static inline
2384 auto operator*(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2385 const typename Function::RangeFieldType& s)
2386 -> decltype(s * asImp(f_))
2387 {
2388 return s * asImp(f_);
2389 }
2390
2392 template<class Parameter, class Function>
2393 static inline
2394 BinaryGridFunctionExpression<SMultiplyOperation, Parameter, Function>
2396 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2397 {
2398 typedef Function FunctionType;
2399 typedef Parameter ParameterType;
2400
2401 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2402 "Right operand must be a GridFunction");
2403
2404 const ParameterType& s(static_cast<const ParameterType&>(s_));
2405 const FunctionType& f(static_cast<const FunctionType&>(f_));
2407
2408 return ExpressionType(s, f);
2409 }
2410
2412 template<class Parameter, class Function>
2413 static inline
2414 auto operator*(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2416 -> decltype(asImp(s_) * asImp(f_))
2417
2418 {
2419 return asImp(s_) * asImp(f_);
2420 }
2421
2423 //
2424 // Division by multiplication with inverse.
2425
2427 template<class LeftFunction, class RightFunction>
2428 auto operator/(const Fem::Function<typename LeftFunction::FunctionSpaceType, LeftFunction>& l_,
2429 const Fem::Function<typename RightFunction::FunctionSpaceType::ScalarFunctionSpaceType, RightFunction>& r_)
2430 -> decltype(asImp(l_) * std::declval<UnaryGridFunctionExpression<InvertOperation, RightFunction> >())
2431 {
2432 typedef RightFunction RightFunctionType;
2433
2434 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2435 "Right operand must be a GridFunction");
2436
2438
2439 return asImp(l_) * InvertExpressionType(asImp(r_));
2440 }
2441
2443 template<class RightFunction>
2444 auto operator/(const typename RightFunction::RangeFieldType& l_,
2445 const Fem::Function<typename RightFunction::FunctionSpaceType::ScalarFunctionSpaceType, RightFunction>& r_)
2446 -> decltype(l_ * std::declval<UnaryGridFunctionExpression<InvertOperation, RightFunction> >())
2447 {
2448 typedef RightFunction RightFunctionType;
2449
2450 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2451 "Right operand must be a GridFunction");
2452
2454
2455 return l_ * InvertExpressionType(asImp(r_));
2456 }
2457
2459 template<class LeftParameter, class RightFunction>
2461 const Fem::Function<typename RightFunction::FunctionSpaceType::ScalarFunctionSpaceType, RightFunction>& r_)
2462 -> decltype(asImp(l_) * std::declval<UnaryGridFunctionExpression<InvertOperation, RightFunction> >())
2463 {
2464 typedef RightFunction RightFunctionType;
2465
2466 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2467 "Right operand must be a GridFunction");
2468
2470
2471 return asImp(l_) * InvertExpressionType(asImp(r_));
2472 }
2473
2475 template<class Function>
2476 static inline
2477 auto operator/(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2478 const typename Function::RangeFieldType& s)
2479 -> decltype((1.0/s) * asImp(f_))
2480 {
2481 return (1.0/s) * asImp(f_);
2482 }
2483
2484 // there is no operator/ yet for parameters, so what.
2485
2487 //
2488 // Addition
2489
2490 template<class LeftFunction, class RightFunction>
2491 static inline
2492 BinaryGridFunctionExpression<PlusOperation, LeftFunction, RightFunction>
2493 operator+(const Fem::Function<typename LeftFunction::FunctionSpaceType, LeftFunction>& f_,
2494 const Fem::Function<typename RightFunction::FunctionSpaceType, RightFunction>& g_)
2495 {
2496 typedef LeftFunction LeftFunctionType;
2497 typedef RightFunction RightFunctionType;
2498
2499 static_assert(std::is_base_of<Fem::HasLocalFunction, LeftFunctionType>::value,
2500 "Left operand must be a GridFunction");
2501 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2502 "Right operand must be a GridFunction");
2503
2504 const LeftFunctionType& f(static_cast<const LeftFunctionType&>(f_));
2505 const RightFunctionType& g(static_cast<const RightFunctionType&>(g_));
2506
2507 typedef BinaryGridFunctionExpression<PlusOperation, LeftFunctionType, RightFunctionType> ExpressionType;
2508
2509 return ExpressionType(f, g);
2510 }
2511
2517 template<class Function>
2518 static inline
2519 BinaryGridFunctionExpression<PlusOperation,
2520 ConstantGridFunction<typename Function::FunctionSpaceType, typename Function::GridPartType>,
2521 Function>
2522 operator+(const typename Function::RangeType& s,
2523 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2524 {
2525 typedef Function FunctionType;
2526
2527 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2528 "Right operand must be a GridFunction");
2529
2530 const FunctionType& f(static_cast<const FunctionType&>(f_));
2531
2532 typedef typename FunctionType::GridPartType GridPartType;
2533 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2534 typedef ConstantGridFunction<FunctionSpaceType, GridPartType> ConstantFunctionType;
2536 ExpressionType;
2537
2538 ConstantFunctionType constant(s, f.space().gridPart());
2539
2540 return ExpressionType(constant, f);
2541 }
2542
2548 template<class Function>
2549 static inline
2550 auto
2551 operator+(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2552 const typename Function::RangeType& s_)
2553 -> decltype(s_ + f_)
2554 {
2555 return s_ + f_;
2556 }
2557
2563 template<class Parameter, class Function>
2564 static inline
2565 BinaryGridFunctionExpression<PlusOperation,
2566 ParameterGridFunction<Parameter, typename Function::FunctionSpaceType,
2567 typename Function::GridPartType>,
2568 Function>
2570 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2571 {
2572 typedef Function FunctionType;
2573 typedef Parameter ParameterType;
2574
2575 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2576 "Right operand must be a GridFunction");
2577
2578 const ParameterType& s(static_cast<const ParameterType&>(s_));
2579 const FunctionType& f(static_cast<const FunctionType&>(f_));
2580
2581 typedef typename FunctionType::GridPartType GridPartType;
2582 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2585 ExpressionType;
2586
2587 ParameterFunctionType parameter(s, f.space().gridPart());
2588
2589 return ExpressionType(parameter, f);
2590 }
2591
2597 template<class Parameter, class Function>
2598 static inline
2599 auto
2600 operator+(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2602 -> decltype(s_ + f_)
2603 {
2604 return s_ + f_;
2605 }
2606
2608 //
2609 // Subtraction
2610
2611 template<class LeftFunction, class RightFunction>
2612 static inline
2613 BinaryGridFunctionExpression<MinusOperation, LeftFunction, RightFunction>
2614 operator-(const Fem::Function<typename LeftFunction::FunctionSpaceType, LeftFunction>& f_,
2615 const Fem::Function<typename RightFunction::FunctionSpaceType, RightFunction>& g_)
2616 {
2617 typedef LeftFunction LeftFunctionType;
2618 typedef RightFunction RightFunctionType;
2619
2620 static_assert(std::is_base_of<Fem::HasLocalFunction, LeftFunctionType>::value,
2621 "Left operand must be a GridFunction");
2622 static_assert(std::is_base_of<Fem::HasLocalFunction, RightFunctionType>::value,
2623 "Right operand must be a GridFunction");
2624
2625 const LeftFunctionType& f(static_cast<const LeftFunctionType&>(f_));
2626 const RightFunctionType& g(static_cast<const RightFunctionType&>(g_));
2627
2628 typedef BinaryGridFunctionExpression<MinusOperation, LeftFunctionType, RightFunctionType> ExpressionType;
2629
2630 return ExpressionType(f, g);
2631 }
2632
2633 template<class Function>
2634 static inline
2635 BinaryGridFunctionExpression<MinusOperation,
2636 ConstantGridFunction<typename Function::FunctionSpaceType::ScalarFunctionSpaceType, typename Function::GridPartType>,
2637 Function>
2638 operator-(const typename Function::RangeFieldType& s,
2639 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2640 {
2641 typedef Function FunctionType;
2642
2643 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2644 "Right operand must be a GridFunction");
2645
2646 const FunctionType& f(static_cast<const FunctionType&>(f_));
2647
2648 typedef typename FunctionType::GridPartType GridPartType;
2649 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2650 typedef typename FunctionSpaceType::ScalarFunctionSpaceType ScalarFunctionSpaceType;
2651 typedef ConstantGridFunction<ScalarFunctionSpaceType, GridPartType> ConstantFunctionType;
2652 typedef BinaryGridFunctionExpression<MinusOperation, ConstantFunctionType, FunctionType>
2653 ExpressionType;
2654
2655 ConstantFunctionType scalar(s, f.space().gridPart());
2656
2657 return ExpressionType(scalar, f);
2658 }
2659
2660 template<class Function>
2661 static inline
2662 BinaryGridFunctionExpression<MinusOperation,
2663 Function,
2664 ConstantGridFunction<typename Function::FunctionSpaceType::ScalarFunctionSpaceType, typename Function::GridPartType> >
2665 operator-(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2666 const typename Function::RangeFieldType& s)
2667 {
2668 typedef Function FunctionType;
2669
2670 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2671 "Left operand must be a GridFunction");
2672
2673 const FunctionType& f(static_cast<const FunctionType&>(f_));
2674
2675 typedef typename FunctionType::GridPartType GridPartType;
2676 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2677 typedef typename FunctionSpaceType::ScalarFunctionSpaceType ScalarFunctionSpaceType;
2678 typedef ConstantGridFunction<ScalarFunctionSpaceType, GridPartType> ConstantFunctionType;
2679 typedef BinaryGridFunctionExpression<MinusOperation, FunctionType, ConstantFunctionType>
2680 ExpressionType;
2681 ConstantFunctionType scalar(s, f.space().gridPart());
2682
2683 return ExpressionType(f, scalar);
2684 }
2685
2691 template<class Parameter, class Function>
2692 static inline
2693 BinaryGridFunctionExpression<MinusOperation,
2694 ParameterGridFunction<Parameter, typename Function::FunctionSpaceType,
2695 typename Function::GridPartType>,
2696 Function>
2698 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2699 {
2700 typedef Function FunctionType;
2701 typedef Parameter ParameterType;
2702
2703 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2704 "Right operand must be a GridFunction");
2705
2706 const ParameterType& s(static_cast<const ParameterType&>(s_));
2707 const FunctionType& f(static_cast<const FunctionType&>(f_));
2708
2709 typedef typename FunctionType::GridPartType GridPartType;
2710 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2713 ExpressionType;
2714
2715 ParameterFunctionType parameter(s, f.space().gridPart());
2716
2717 return ExpressionType(parameter, f);
2718 }
2719
2725 template<class Parameter, class Function>
2726 static inline
2727 BinaryGridFunctionExpression<MinusOperation,
2728 Function,
2729 ParameterGridFunction<Parameter, typename Function::FunctionSpaceType,
2730 typename Function::GridPartType> >
2731 operator-(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
2733
2734 {
2735 typedef Function FunctionType;
2736 typedef Parameter ParameterType;
2737
2738 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2739 "Right operand must be a GridFunction");
2740
2741 const ParameterType& s(static_cast<const ParameterType&>(s_));
2742 const FunctionType& f(static_cast<const FunctionType&>(f_));
2743
2744 typedef typename FunctionType::GridPartType GridPartType;
2745 typedef typename FunctionType::FunctionSpaceType FunctionSpaceType;
2748 ExpressionType;
2749
2750 ParameterFunctionType parameter(s, f.space().gridPart());
2751
2752 return ExpressionType(f, parameter);
2753 }
2754
2756 //
2757 // Identity
2758
2762 template<class Function>
2763 static inline
2764 UnaryGridFunctionExpression<IdentityOperation, Function>
2765 operator*(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2766 {
2767 typedef Function FunctionType;
2768
2769 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2770 "Operand must be a GridFunction");
2771
2772 const FunctionType& f(static_cast<const FunctionType&>(f_));
2773
2775
2776 return ExpressionType(f);
2777 }
2778
2780 //
2781 // unary plus is redirected to unary operator*()
2782
2784 template<class Function>
2785 static inline
2786 auto
2787 operator+(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2788 -> decltype(*asImp(f_))
2789 {
2790 return *asImp(f_);
2791 }
2792
2794 //
2795 // Unary minus
2796
2797 template<class Function>
2798 static inline
2799 UnaryGridFunctionExpression<MinusOperation, Function>
2800 operator-(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2801 {
2802 typedef Function FunctionType;
2803
2804 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2805 "Operand must be a GridFunction");
2806
2807 const FunctionType& f(static_cast<const FunctionType&>(f_));
2808
2809 typedef UnaryGridFunctionExpression<MinusOperation, FunctionType> ExpressionType;
2810
2811 return ExpressionType(f);
2812 }
2813
2815 //
2816 // Component-wise exponentiation
2817
2819 template<class Function>
2820 static inline
2821 UnaryGridFunctionExpression<ExpOperation, Function>
2822 exp(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2823 {
2824 typedef Function FunctionType;
2825
2826 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2827 "Operand must be a GridFunction");
2828
2829 const FunctionType& f(static_cast<const FunctionType&>(f_));
2830
2832
2833 return ExpressionType(f);
2834 }
2835
2837 //
2838 // Component-wise log
2839
2841 template<class Function>
2842 static inline
2843 UnaryGridFunctionExpression<LogOperation, Function>
2844 log(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2845 {
2846 typedef Function FunctionType;
2847
2848 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2849 "Operand must be a GridFunction");
2850
2851 const FunctionType& f(static_cast<const FunctionType&>(f_));
2852
2854
2855 return ExpressionType(f);
2856 }
2857
2859 //
2860 // Component-wise square-root
2861
2863 template<class Function>
2864 static inline
2865 UnaryGridFunctionExpression<SqrtOperation, Function>
2866 sqrt(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2867 {
2868 typedef Function FunctionType;
2869
2870 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2871 "Operand must be a GridFunction");
2872
2873 const FunctionType& f(static_cast<const FunctionType&>(f_));
2874
2876
2877 return ExpressionType(f);
2878 }
2879
2881 template<class Function>
2882 static inline
2883 UnaryGridFunctionExpression<SquareOperation, Function>
2884 sqr(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2885 {
2886 typedef Function FunctionType;
2887
2888 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2889 "Operand must be a GridFunction");
2890
2891 const FunctionType& f(static_cast<const FunctionType&>(f_));
2892
2894
2895 return ExpressionType(f);
2896 }
2897
2899 //
2900 // Component-wise sin
2901
2903 template<class Function>
2904 static inline
2905 UnaryGridFunctionExpression<SinOperation, Function>
2906 sin(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2907 {
2908 typedef Function FunctionType;
2909
2910 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2911 "Operand must be a GridFunction");
2912
2913 const FunctionType& f(static_cast<const FunctionType&>(f_));
2914
2916
2917 return ExpressionType(f);
2918 }
2919
2921 //
2922 // Component-wise cos
2923
2925 template<class Function>
2926 static inline
2927 UnaryGridFunctionExpression<CosOperation, Function>
2928 cos(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2929 {
2930 typedef Function FunctionType;
2931
2932 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2933 "Operand must be a GridFunction");
2934
2935 const FunctionType& f(static_cast<const FunctionType&>(f_));
2936
2938
2939 return ExpressionType(f);
2940 }
2941
2943 //
2944 // Component-wise tan
2945
2947 template<class Function>
2948 static inline
2949 UnaryGridFunctionExpression<TanOperation, Function>
2950 tan(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2951 {
2952 typedef Function FunctionType;
2953
2954 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2955 "Operand must be a GridFunction");
2956
2957 const FunctionType& f(static_cast<const FunctionType&>(f_));
2958
2960
2961 return ExpressionType(f);
2962 }
2963
2965 //
2966 // Component-wise arctan
2967
2969 template<class Function>
2970 static inline
2971 UnaryGridFunctionExpression<AtanOperation, Function>
2972 atan(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2973 {
2974 typedef Function FunctionType;
2975
2976 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2977 "Operand must be a GridFunction");
2978
2979 const FunctionType& f(static_cast<const FunctionType&>(f_));
2980
2982
2983 return ExpressionType(f);
2984 }
2985
2987 //
2988 // Component-wise arcsin
2989
2991 template<class Function>
2992 static inline
2993 UnaryGridFunctionExpression<AsinOperation, Function>
2994 asin(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
2995 {
2996 typedef Function FunctionType;
2997
2998 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
2999 "Operand must be a GridFunction");
3000
3001 const FunctionType& f(static_cast<const FunctionType&>(f_));
3002
3004
3005 return ExpressionType(f);
3006 }
3007
3009 //
3010 // Component-wise arccos
3011
3013 template<class Function>
3014 static inline
3015 UnaryGridFunctionExpression<AcosOperation, Function>
3016 acos(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
3017 {
3018 typedef Function FunctionType;
3019
3020 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
3021 "Operand must be a GridFunction");
3022
3023 const FunctionType& f(static_cast<const FunctionType&>(f_));
3024
3026
3027 return ExpressionType(f);
3028 }
3029
3040 template<class FunctionSpace, class Expression>
3041 static inline
3042 Expression
3044 {
3045 return f_.expression();
3046 }
3047
3049 template<class FunctionSpace, class GridPart>
3050 static inline
3051 ZeroGridFunction<FunctionSpace, GridPart>
3053 {
3054 return f_;
3055 }
3056
3063 template<class F>
3064 static inline
3066 -> decltype(+expr.function())
3067 {
3068 return +expr.function();
3069 };
3070
3072 template<class F1, class F2>
3073 static inline
3075 -> decltype(expr.rightFunction() - expr.leftFunction())
3076 {
3077 return expr.rightFunction() - expr.leftFunction();
3078 };
3079
3081 template<class F1, class F2>
3082 static inline
3085 -> decltype(expr.leftFunction().function() - expr.rightFunction())
3086 {
3087 return expr.leftFunction().function() - expr.rightFunction();
3088 };
3089
3091
3092
3099 template<class Field, class Function>
3100 static inline
3101 auto operator*(const Field& s_,
3102 const Fem::Function<typename Function::FunctionSpaceType,
3104 -> decltype((s_ * asImp(f_).scalar()) * asImp(f_).function())
3105 {
3106 return (s_ * asImp(f_).scalar()) * asImp(f_).function();
3107 }
3108
3110 template<class Field, class Function1, class Function2>
3111 static inline
3112 auto operator*(const Fem::Function<typename Function1::FunctionSpaceType,
3114 const Fem::Function<typename Function2::FunctionSpaceType, Function2>& f2_)
3115 -> decltype(asImp(f1_).scalar() * (asImp(f1_).function() * asImp(f2_)))
3116 {
3117 return asImp(f1_).scalar() * (asImp(f1_).function() * asImp(f2_));
3118 }
3119
3121 template<class Field, class Function1, class Function2>
3122 static inline
3123 auto operator*(const Fem::Function<typename Function1::FunctionSpaceType, Function1>& f1_,
3124 const Fem::Function<typename Function2::FunctionSpaceType,
3126 -> decltype(asImp(f2_).scalar() * (asImp(f1_) * asImp(f2_).function()))
3127 {
3128 return asImp(f2_).scalar() * (asImp(f1_) * asImp(f2_).function());
3129 }
3130
3132 template<class Field1, class Field2, class Function1, class Function2>
3133 static inline
3134 auto operator*(const Fem::Function<typename Function1::FunctionSpaceType,
3136 const Fem::Function<typename Function2::FunctionSpaceType,
3138 -> decltype((asImp(f1_).scalar() * asImp(f2_).scalar()) * (asImp(f1_).function() * asImp(f2_).function()))
3139 {
3140 return (asImp(f1_).scalar() * asImp(f2_).scalar()) * (asImp(f1_).function() * asImp(f2_).function());
3141 }
3142
3144 template<class Parameter, class Function>
3145 static inline
3147 const Fem::Function<typename Function::FunctionSpaceType,
3149 -> decltype(asImp(f_).scalar() * (asImp(p_) * asImp(f_).function()))
3150 {
3151 return asImp(f_).scalar() * (asImp(p_) * asImp(f_).function());
3152 }
3153
3155 template<class Field, class Parameter, class Function>
3156 static inline
3157 auto operator*(const BinaryParameterExpression<SMultiplyOperation, Field, Parameter>& p_,
3158 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
3159 -> decltype(p_.left() * (p_.right() * asImp(f_)))
3160 {
3161 return p_.left() * (p_.right() * asImp(f_));
3162 }
3163
3167 template<class Parameter, class ZeroExpression>
3168 static inline
3169 ZeroExpression
3170 operator*(const BinaryParameterExpression<SMultiplyOperation, typename ZeroExpression::RangeFieldType, Parameter>& p_,
3172 {
3173 return *z_;
3174 }
3175
3177
3183 template<class ZeroExpression>
3184 static inline
3185 ZeroExpression
3186 operator*(const typename ZeroExpression::RangeFieldType& s,
3188 {
3189 return *z;
3190 }
3191
3193 template<class Parameter, class ZeroExpression>
3194 static inline
3195 ZeroExpression
3198 {
3199 return *z;
3200 }
3201
3203 template<class Function, class ZeroExpression>
3204 static inline
3205 auto
3207 const Fem::Function<typename Function::FunctionSpaceType, Function>& g_)
3208 -> decltype(zeroProduct(*z_, g_))
3209 {
3210 return zeroProduct(*z_, g_);
3211 }
3212
3214 template<class Function, class ZeroExpression>
3215 static inline
3216 auto
3217 operator*(const Fem::Function<typename Function::FunctionSpaceType, Function>& g_,
3219 -> decltype(z_ * g_)
3220 {
3221 return z_ * g_;
3222 }
3223
3225 template<class ZeroExpression1, class ZeroExpression2>
3226 static inline
3227 auto
3230 -> decltype(zeroProduct(*z1, z2))
3231 {
3232 return zeroProduct(*z1, z2);
3233 }
3234
3237 template<class Function, class ZeroExpression>
3238 static inline
3239 auto
3241 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
3242 -> decltype(*asImp(f_))
3243 {
3244 typedef Function FunctionType;
3245
3246 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
3247 "Operand must be a GridFunction");
3248 // should also assert same space, maybe
3249
3250 return *asImp(f_);
3251 }
3252
3255 template<class Function, class ZeroExpression>
3256 static inline
3257 auto
3258 operator+(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
3260 -> decltype(*asImp(f_))
3261 {
3262 typedef Function FunctionType;
3263
3264 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
3265 "Operand must be a GridFunction");
3266
3267 return *asImp(f_);
3268 }
3269
3272 template<class ZeroExpression1, class ZeroExpression2>
3273 static inline
3274 ZeroExpression1
3277 {
3278 return *z1_; // just choose one
3279 }
3280
3283 template<class Function, class ZeroExpression>
3284 static inline
3285 auto
3286 operator-(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
3288 -> decltype(*asImp(f_))
3289 {
3290 typedef Function FunctionType;
3291
3292 static_assert(std::is_base_of<Fem::HasLocalFunction, FunctionType>::value,
3293 "Operand must be a GridFunction");
3294
3295 return *asImp(f_);
3296 }
3297
3300 template<class Function, class ZeroExpression>
3301 static inline
3302 auto
3304 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
3305 -> decltype(-asImp(f_))
3306 {
3307 return -asImp(f_);
3308 }
3309
3313 template<class ZeroExpression1, class ZeroExpression2>
3314 static inline
3315 ZeroExpression1
3318 {
3319 return *z1_; // just choose one
3320 }
3321
3323 template<class ZeroExpression>
3324 static inline
3325 ZeroExpression
3327 {
3328 return *z_;
3329 }
3330
3332
3340 template<class Function>
3341 static inline
3342 auto operator*(const decltype(oneFunction(std::declval<Function>()))& one,
3343 const Fem::Function<typename Function::FunctionSpaceType, Function>& f_)
3344 -> decltype(*asImp(f_))
3345 {
3346 return *asImp(f_);
3347 }
3348
3349 template<class Function>
3350 static inline
3351 auto operator*(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
3352 const decltype(oneFunction(std::declval<Function>()))& one)
3353 -> decltype(*asImp(f_))
3354 {
3355 return *asImp(f_);
3356 }
3357
3358 template<class FunctionSpace, class GridPart>
3359 static inline
3360 FractionGridFunction<typename FunctionSpace::ScalarFunctionSpaceType, GridPart, 1L, 1UL>
3361 operator*(const FractionGridFunction<FunctionSpace, GridPart, 1L, 1UL>& l,
3362 const FractionGridFunction<FunctionSpace, GridPart, 1L, 1UL>& r)
3363 {
3364 typedef
3365 FractionGridFunction<typename FunctionSpace::ScalarFunctionSpaceType, GridPart, 1L, 1UL>
3366 ExpressionType;
3367
3368 return ExpressionType(l.gridPart());
3369 }
3370
3371 template<class Function>
3372 static inline
3373 auto operator/(const Fem::Function<typename Function::FunctionSpaceType, Function>& f_,
3374 const decltype(oneFunction(std::declval<Function>()))& one)
3375 -> decltype(*asImp(f_))
3376 {
3377 return *asImp(f_);
3378 }
3379
3380 template<class FunctionSpace, class GridPart>
3381 static inline
3382 FractionGridFunction<typename FunctionSpace::ScalarFunctionSpaceType, GridPart, 1L, 1UL>
3383 operator/(const FractionGridFunction<FunctionSpace, GridPart, 1L, 1UL>& l,
3384 const FractionGridFunction<FunctionSpace, GridPart, 1L, 1UL>& r)
3385 {
3386 typedef
3387 FractionGridFunction<typename FunctionSpace::ScalarFunctionSpaceType, GridPart, 1L, 1UL>
3388 ExpressionType;
3389
3390 return ExpressionType(l.gridPart());
3391 }
3392
3394
3396
3398
3400
3401 } // namespace ACFem
3402
3403 namespace Fem {
3404 using ACFem::operator+;
3405 using ACFem::operator-;
3406 using ACFem::operator*;
3407 using ACFem::operator/;
3408 }
3409
3410} // namespace Dune
3411
3412#endif // __DUNE_ACFEM_GRIDFUNCTIONEXPRESSION_HH__
General local function object for binary grid-function expressions.
Definition: gridfunctionexpression.hh:613
const IntersectionType & intersection() const
Definition: gridfunctionexpression.hh:769
LocalFunction(const EntityType &entity, const DiscreteFunctionType &df)
constructor initializing local function
Definition: gridfunctionexpression.hh:649
LocalFunction(const EntityType &entity, const IntersectionType &intersection, const DiscreteFunctionType &df)
Construct from intersection.
Definition: gridfunctionexpression.hh:741
DiscreteFunctionSpaceType::FunctionSpaceType FunctionSpaceType
function space type
Definition: gridfunctionexpression.hh:633
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:646
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:640
int order() const
Return a bound on or suggestion for the piece-wise polynomial order.
Definition: gridfunctionexpression.hh:711
void init(const EntityType &entity, const IntersectionType &intersection)
Definition: gridfunctionexpression.hh:754
bool supported() const
Construct from intersection.
Definition: gridfunctionexpression.hh:776
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:644
void init(const EntityType &entity)
init local function
Definition: gridfunctionexpression.hh:714
void evaluateQuadrature(const QuadratureType &quadrature, VectorType &values) const
evaluate function or jacobian of function for given quadrature
Definition: gridfunctionexpression.hh:704
void evaluate(const PointType &x, RangeType &ret) const
evaluate local function
Definition: gridfunctionexpression.hh:671
void hessian(const PointType &x, HessianRangeType &ret) const
hessian of local function
Definition: gridfunctionexpression.hh:693
const EntityType & entity() const
Return a reference to the currently active entity.
Definition: gridfunctionexpression.hh:723
void jacobian(const PointType &x, JacobianRangeType &ret) const
jacobian of local function
Definition: gridfunctionexpression.hh:682
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:636
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:642
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:638
S-multiplication with RangeFieldType scalars.
Definition: gridfunctionexpression.hh:869
LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection)
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1029
void hessian(const DomainType &global, HessianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:968
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:918
TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
type of discrete function space
Definition: gridfunctionexpression.hh:898
BinaryGridFunctionExpressionTraits< SMultiplyOperation, Factor, ContainedFunctionType > TraitsType
type of traits
Definition: gridfunctionexpression.hh:887
Factor FactorType
type of factor
Definition: gridfunctionexpression.hh:881
DiscreteFunctionSpaceType::GridPartType GridPartType
type of gridPart
Definition: gridfunctionexpression.hh:904
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:912
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:910
const LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection) const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1022
const LocalFunctionType localFunction(const EntityType &entity) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity) const.
Definition: gridfunctionexpression.hh:981
void jacobian(const DomainType &global, JacobianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:957
ThisType DiscreteFunctionType
type of discrete function (self)
Definition: gridfunctionexpression.hh:895
DiscreteFunctionSpaceType::GridType GridType
type of grid
Definition: gridfunctionexpression.hh:907
const DiscreteFunctionSpaceType & space() const
See Fem::DiscreteFunctionInterface::space() const.
Definition: gridfunctionexpression.hh:1004
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:914
BoundaryFunctionTraits< ContainedFunctionType >::IndicatorType IndicatorType
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1013
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:916
void evaluate(const DomainType &global, RangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:950
const std::string & expressionName() const
See Fem::DiscreteFunctionInterface::name() const.
Definition: gridfunctionexpression.hh:993
Function ContainedFunctionType
type of function
Definition: gridfunctionexpression.hh:884
LocalFunctionType localFunction(const EntityType &entity)
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)
Definition: gridfunctionexpression.hh:987
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:920
TraitsType::EntityType EntityType
type of codim 0 entity
Definition: gridfunctionexpression.hh:923
LocalFunction LocalFunctionType
type of local function to export
Definition: gridfunctionexpression.hh:933
IndicatorType indicator() const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1036
BinaryGridFunctionExpression implements point-wise vector-space operations for GridFunction-types.
Definition: gridfunctionexpression.hh:423
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:460
const DiscreteFunctionSpaceType & space() const
See Fem::DiscreteFunctionInterface::space() const.
Definition: gridfunctionexpression.hh:550
void evaluate(const DomainType &global, RangeType &result) const
evaluate function at global coordinates
Definition: gridfunctionexpression.hh:504
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:470
BinaryGridFunctionExpressionTraits< BinOp, LeftFunctionType, RightFunctionType > TraitsType
type of traits
Definition: gridfunctionexpression.hh:428
GridPartType::IntersectionType IntersectionType
type of intersection
Definition: gridfunctionexpression.hh:476
const LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)> const.
Definition: gridfunctionexpression.hh:575
void hessian(const DomainType &global, HessianRangeType &result) const
hessian at global coordinates
Definition: gridfunctionexpression.hh:521
DiscreteFunctionSpaceType::GridType GridType
type of grid
Definition: gridfunctionexpression.hh:453
TraitsType::EntityType EntityType
type of codim 0 entity
Definition: gridfunctionexpression.hh:473
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:468
TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
type of discrete function space
Definition: gridfunctionexpression.hh:444
TraitsType::IndicatorType IndicatorType
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:565
ThisType DiscreteFunctionType
type of discrete function (self)
Definition: gridfunctionexpression.hh:441
const std::string & expressionName() const
See Fem::DiscreteFunctionInterface::name() const.
Definition: gridfunctionexpression.hh:543
BinaryGridFunctionExpression(const LeftFunctionType &f, const RightFunctionType &g)
Constructor.
Definition: gridfunctionexpression.hh:489
const LocalFunctionType localFunction(const EntityType &entity) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)> const.
Definition: gridfunctionexpression.hh:531
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:466
LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection)
See Fem::DiscreteFunctionInterface::localFunction(const IntersectionType &intersection)
Definition: gridfunctionexpression.hh:583
DiscreteFunctionSpaceType::GridPartType GridPartType
type of gridPart
Definition: gridfunctionexpression.hh:450
LocalFunction LocalFunctionType
type of local function to export
Definition: gridfunctionexpression.hh:486
LocalFunctionType localFunction(const EntityType &entity)
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)
Definition: gridfunctionexpression.hh:537
IndicatorType indicator() const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:589
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:464
@ emptySupport
Definition: gridfunctionexpression.hh:569
@ globalSupport
Definition: gridfunctionexpression.hh:571
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:462
void jacobian(const DomainType &global, JacobianRangeType &result) const
jacobian at global coordinates
Definition: gridfunctionexpression.hh:510
ConstantGridFunction implements a constant function.
Definition: constantfunction.hh:108
A grid-function which is constant with a scalar fractional value.
Definition: constantfunction.hh:598
A class providing some basic functionality common to all expressions.
Definition: gridfunctionexpressionbase.hh:35
ParameterGridFunction implements a constant function where the value is defined by something which fu...
Definition: parameterfunction.hh:112
Parameters are quasi-constant quantities, like the time-step size in one time-step when solving trans...
Definition: parameterinterface.hh:80
Specialize for the identity wrapper operation.
Definition: gridfunctionexpression.hh:1668
TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
type of discrete function space
Definition: gridfunctionexpression.hh:1689
const DiscreteFunctionSpaceType & space() const
See Fem::DiscreteFunctionInterface::space() const.
Definition: gridfunctionexpression.hh:1784
DiscreteFunctionSpaceType::GridType GridType
type of grid
Definition: gridfunctionexpression.hh:1698
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:1705
LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection)
See Fem::DiscreteFunctionInterface::localFunction(const IntersectionType &intersection)
Definition: gridfunctionexpression.hh:1811
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:1711
DiscreteFunctionSpaceType::GridPartType GridPartType
type of gridPart
Definition: gridfunctionexpression.hh:1695
const LocalFunctionType localFunction(const EntityType &entity) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity) const.
Definition: gridfunctionexpression.hh:1765
LocalFunctionType localFunction(const EntityType &entity)
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)
Definition: gridfunctionexpression.hh:1771
IndicatorType indicator() const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1818
void hessian(const DomainType &global, HessianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1755
LocalFunction LocalFunctionType
type of local function to export
Definition: gridfunctionexpression.hh:1724
const LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)> const.
Definition: gridfunctionexpression.hh:1803
BoundaryFunctionTraits< FunctionType >::IndicatorType IndicatorType
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1793
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:1707
const std::string & expressionName() const
See Fem::DiscreteFunctionInterface::name() const.
Definition: gridfunctionexpression.hh:1777
void evaluate(const DomainType &global, RangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1739
UnaryGridFunctionExpressionTraits< OperationTagType, FunctionType > TraitsType
type of traits
Definition: gridfunctionexpression.hh:1680
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:1703
TraitsType::EntityType EntityType
type of codim 0 entity
Definition: gridfunctionexpression.hh:1714
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:1709
void jacobian(const DomainType &global, JacobianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1745
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:1701
For the sake of BoundarySupportedFunction we need to specialize the InvertOperation.
Definition: gridfunctionexpression.hh:2008
void hessian(const DomainType &global, HessianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:2095
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:2043
const DiscreteFunctionSpaceType & space() const
See Fem::DiscreteFunctionInterface::space() const.
Definition: gridfunctionexpression.hh:2124
const LocalFunctionType localFunction(const EntityType &entity) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity) const.
Definition: gridfunctionexpression.hh:2105
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:2041
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:2049
DiscreteFunctionSpaceType::GridType GridType
type of grid
Definition: gridfunctionexpression.hh:2038
TraitsType::EntityType EntityType
type of codim 0 entity
Definition: gridfunctionexpression.hh:2054
UnaryGridFunctionExpressionTraits< OperationTagType, FunctionType > TraitsType
type of traits
Definition: gridfunctionexpression.hh:2020
LocalFunctionType localFunction(const EntityType &entity)
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)
Definition: gridfunctionexpression.hh:2111
const LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)> const.
Definition: gridfunctionexpression.hh:2143
LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection)
See Fem::DiscreteFunctionInterface::localFunction(const IntersectionType &intersection)
Definition: gridfunctionexpression.hh:2151
TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
type of discrete function space
Definition: gridfunctionexpression.hh:2029
IndicatorType indicator() const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:2158
LocalFunction LocalFunctionType
type of local function to export
Definition: gridfunctionexpression.hh:2064
EntireBoundaryIndicatorType IndicatorType
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:2133
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:2047
void evaluate(const DomainType &global, RangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:2079
const std::string & expressionName() const
See Fem::DiscreteFunctionInterface::name() const.
Definition: gridfunctionexpression.hh:2117
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:2045
void jacobian(const DomainType &global, JacobianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:2085
DiscreteFunctionSpaceType::GridPartType GridPartType
type of gridPart
Definition: gridfunctionexpression.hh:2035
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:2051
UnaryGridFunctionExpression implements point-wise unary operations for GridFunction-types.
Definition: gridfunctionexpression.hh:1328
const std::string & expressionName() const
See Fem::DiscreteFunctionInterface::name() const.
Definition: gridfunctionexpression.hh:1439
DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian type (from function space)
Definition: gridfunctionexpression.hh:1371
LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection)
See Fem::DiscreteFunctionInterface::localFunction(const IntersectionType &intersection)
Definition: gridfunctionexpression.hh:1473
DiscreteFunctionSpaceType::GridType GridType
type of grid
Definition: gridfunctionexpression.hh:1360
ThisType DiscreteFunctionType
export DiscreteFunctionType
Definition: gridfunctionexpression.hh:1348
UnaryGridFunctionExpressionTraits< UnOp, FunctionType > TraitsType
type of traits
Definition: gridfunctionexpression.hh:1331
void evaluate(const DomainType &global, RangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1401
LocalFunctionType localFunction(const EntityType &entity)
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)
Definition: gridfunctionexpression.hh:1433
DiscreteFunctionSpaceType::RangeType RangeType
range type (from function space)
Definition: gridfunctionexpression.hh:1369
@ globalSupport
Definition: gridfunctionexpression.hh:1461
@ emptySupport
Definition: gridfunctionexpression.hh:1459
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
range type (from function space)
Definition: gridfunctionexpression.hh:1365
void hessian(const DomainType &global, HessianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1417
DiscreteFunctionSpaceType::GridPartType GridPartType
type of gridPart
Definition: gridfunctionexpression.hh:1357
const DiscreteFunctionSpaceType & space() const
See Fem::DiscreteFunctionInterface::space() const.
Definition: gridfunctionexpression.hh:1446
TraitsType::EntityType EntityType
type of codim 0 entity
Definition: gridfunctionexpression.hh:1376
DiscreteFunctionSpaceType::DomainType DomainType
domain type (from function space)
Definition: gridfunctionexpression.hh:1367
DiscreteFunctionSpaceType::DomainFieldType DomainFieldType
domain type (from function space)
Definition: gridfunctionexpression.hh:1363
BoundaryFunctionTraits< FunctionType >::IndicatorType IndicatorType
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1455
void jacobian(const DomainType &global, JacobianRangeType &result) const
evaluate function on local coordinate local
Definition: gridfunctionexpression.hh:1407
IndicatorType indicator() const
Resulting Boundary indicator type.
Definition: gridfunctionexpression.hh:1480
TraitsType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
type of discrete function space
Definition: gridfunctionexpression.hh:1351
const LocalFunctionType localFunction(const EntityType &entity, const IntersectionType &intersection) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity)> const.
Definition: gridfunctionexpression.hh:1465
LocalFunction LocalFunctionType
type of local function to export
Definition: gridfunctionexpression.hh:1386
const LocalFunctionType localFunction(const EntityType &entity) const
See Fem::DiscreteFunctionInterface::localFunction(const EntityType &entity) const.
Definition: gridfunctionexpression.hh:1427
DiscreteFunctionSpaceType::HessianRangeType HessianRangeType
hessian type (from function space)
Definition: gridfunctionexpression.hh:1373
A base class for zero function expression.
Definition: gridfunctionexpressionbase.hh:110
A grid-function always returning 0.
Definition: constantfunction.hh:352
static auto operator/(const Fem::Function< typename Left::FunctionSpaceType, Left > &f_, const BoundarySupportedFunction< Right, RightInd > &g_) -> decltype(asBndryFct(asEssBndryFct(f_/asExprArg(g_))))
f / Wrapped(g) = Wrapped(f / g)
Definition: boundaryfunctionexpression.hh:657
const Implementation & asImp(const Fem::BartonNackmanInterface< Interface, Implementation > &arg)
Up-cast to the implementation for any Fem::BartonNackmanInterface.
Definition: expressionoperations.hh:71
auto constantFunction(const Fem::Function< typename F::FunctionSpaceType, F > &f_, const typename F::FunctionSpaceType::RangeType &c) -> ConstantGridFunction< typename F::FunctionSpaceType, typename F::GridPartType >
Generate a ConstantGridFunction on the same function space.
Definition: gridfunctionexpression.hh:191
static ZeroGridFunction< typename FunctionMultiplicationResultTraits< typename F1::FunctionSpaceType, typename F2::FunctionSpaceType >::FunctionSpaceType, typename F1::GridPartType > zeroProduct(const F1 &f1, const F2 &f2)
Generate the proper zero-functions for products evaluating to zero.
Definition: gridfunctionexpression.hh:116
static UnaryGridFunctionExpression< LogOperation, Function > log(const Fem::Function< typename Function::FunctionSpaceType, Function > &f_)
Component-wise logarithm, log(f)(x) = [log(f(x)_0),...,log(f(x)_N].
Definition: gridfunctionexpression.hh:2844
static const Fem::Function< typename GridFunction::FunctionSpaceType, GridFunction > & asExprFunction(const GridFunction &arg)
Cast to the proper underlying expression class; for the generic case this is just the underlying Fem:...
Definition: gridfunctionexpression.hh:61
auto zeroFunction(const Fem::Function< typename F::FunctionSpaceType, F > &f_) -> ZeroGridFunction< typename F::FunctionSpaceType, typename F::GridPartType >
Generate a ZeroGridFunction on the same function space.
Definition: gridfunctionexpression.hh:183
FractionGridFunction< typename FunctionSpace::FunctionSpaceType::ScalarFunctionSpaceType, GridPart, 1L, 1UL > oneFunction(const FunctionSpace &space, const GridPart &gridPart)
Generate a proper constant-one function from the given Fem::FunctionSpace and Fem::GridPart.
Definition: gridfunctionexpression.hh:157
ParameterValue< Value >::ResultType parameterValue(const Value &value)
Return the unaltered argument for non-parameters and otherwise the parameter value.
Definition: parameterinterface.hh:263
BinaryParameterExpression< MultiplyOperation, Left, Right > operator*(const ParameterInterface< Left > &left, const ParameterInterface< Right > &right)
Scalar product between parameters.
Definition: parameterexpression.hh:321
Multiplication from the left by a scalar constant.
Definition: functionoperations.hh:567
traits of BinaryFunctionExpression
Definition: functionexpression.hh:96
GridPartType::IndexSetType IndexSetType
type of IndexSet
Definition: gridfunctionexpression.hh:821
GridPartType::template Codim< 0 >::IteratorType IteratorType
type of iterator
Definition: gridfunctionexpression.hh:819
traits of BinaryGridFunctionExpression
Definition: gridfunctionexpression.hh:378
GridPartType::IndexSetType IndexSetType
type of IndexSet
Definition: gridfunctionexpression.hh:402
GridPartType::template Codim< 0 >::IteratorType IteratorType
type of iterator
Definition: gridfunctionexpression.hh:400
Helper traits in order to treat function w/o dedicated boundary support in the same way,...
Definition: boundarysupportedfunction.hh:52
const ExpressionType & expression() const
Return a const reference to the underlying expression.
Definition: expressionoperations.hh:42
A traits class in order to collect properties of expressions.
Definition: expressionoperations.hh:465
Identity, i.e. just wrap the object.
Definition: expressionoperations.hh:284
Inversion of an object.
Definition: expressionoperations.hh:275
Negation of IsPieceWiseConstant, in order to have something to branch to with std::conditional.
Definition: gridfunctionexpressionbase.hh:119
Tag type, consequences are zero Jacobian and Hessian.
Definition: gridfunctionexpressionbase.hh:114
Addition of two objects.
Definition: expressionoperations.hh:230
Multiplication by scalars from the left.
Definition: expressionoperations.hh:257
The default result type is void in order to trigger compilation errors.
Definition: expressionoperations.hh:589
traits of UnaryGridFunctionExpression
Definition: gridfunctionexpression.hh:1270
GridPartType::template Codim< 0 >::IteratorType IteratorType
type of iterator
Definition: gridfunctionexpression.hh:1287
GridPartType::IndexSetType IndexSetType
type of IndexSet
Definition: gridfunctionexpression.hh:1289
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)