DUNE-FEM (unstable)

gridfunctionview.hh
1#ifndef DUNE_FEMPY_FUNCTION_GRIDFUNCTIONVIEW_HH
2#define DUNE_FEMPY_FUNCTION_GRIDFUNCTIONVIEW_HH
3
4#include <dune/fem/function/common/discretefunction.hh>
5#include <dune/fem/function/localfunction/bindable.hh>
6#include <dune/fem/function/localfunction/const.hh>
7#include <dune/fem/common/intersectionside.hh>
8
9namespace Dune
10{
11
12 namespace Fem
13 {
14
15 // GridFunctionView
16 // ----------------
17 struct IsGridFunctionView {};
18
19 template< class GF, bool isDiscreteFunction = std::is_base_of< Fem::IsDiscreteFunction, GF >::value >
20 struct GridFunctionView;
21
22 template< class GF >
23 struct GridFunctionView< GF, false >
24 : public BindableGridFunction<typename GF::GridPartType, typename GF::RangeType>,
25 public IsGridFunctionView
26 {
27 using Base = BindableGridFunction<typename GF::GridPartType, typename GF::RangeType>;
28 typedef typename GF::EntityType Entity;
29 typedef typename GF::RangeType Value;
30
31 typedef typename Entity::Geometry::LocalCoordinate LocalCoordinate;
32
33 GridFunctionView ( const GF &gf )
34 : Base(gf.gridPart())
35 , localFunction_( gf ) {}
36
37 Value operator() ( const LocalCoordinate &x ) const
38 {
39 Value value;
40 localFunction_.evaluate( x, value );
41 return value;
42 }
43
44 void bind ( const Entity &entity ) { localFunction_.bind( entity ); }
45 void unbind () {}
46 template <class IntersectionType>
47 void bind(const IntersectionType &intersection, IntersectionSide side)
48 { localFunction_.bind(intersection, side); }
49
50 private:
51 Dune::Fem::ConstLocalFunction<GF> localFunction_;
52 };
53
54 template< class GF >
55 struct GridFunctionView< GF, true >
56 : public BindableGridFunctionWithSpace<typename GF::GridPartType, typename GF::RangeType>,
57 public IsGridFunctionView
58 {
59 using Base = BindableGridFunctionWithSpace<typename GF::GridPartType, typename GF::RangeType>;
60 typedef typename GF::EntityType Entity;
61 typedef typename GF::RangeType Value;
62
63 typedef typename Entity::Geometry::LocalCoordinate LocalCoordinate;
64
65 private:
66 typedef typename GF::DiscreteFunctionSpaceType DiscreteFunctionSpace;
67 typedef typename DiscreteFunctionSpace::BasisFunctionSetType BasisFunctionSet;
68
69 public:
70 GridFunctionView ( const GF &gf )
71 : Base(gf.gridPart(), gf.name(), gf.order())
72 , gf_( gf )
73 {
74 localDofVector_.reserve( DiscreteFunctionSpace::localBlockSize * space().blockMapper().maxNumDofs() );
75 }
76
77 Value operator() ( const LocalCoordinate &x ) const
78 {
79 Value value;
80 basisFunctionSet_.evaluateAll( x, localDofVector_, value );
81 return value;
82 }
83
84 void bind ( const Entity &entity )
85 {
86 basisFunctionSet_ = space().basisFunctionSet( entity );
87 localDofVector_.resize( basisFunctionSet_.size() );
88 gf_.getLocalDofs( entity, localDofVector_ );
89 }
90
91 void unbind ()
92 {
93 basisFunctionSet_ = BasisFunctionSet();
94 localDofVector_.resize( 0u );
95 }
96
97 template <class IntersectionType>
98 void bind(const IntersectionType &intersection, IntersectionSide side)
99 { defaultIntersectionBind(gf_, intersection, side); }
100
101 private:
102 const DiscreteFunctionSpace &space () const { return gf_.space(); }
103
104 const GF &gf_;
105 BasisFunctionSet basisFunctionSet_;
107 };
108
109
110
111 // localFunction
112 // -------------
113
114 template< class GF,
115 std::enable_if_t< !std::is_base_of< Fem::IsGridFunctionView, GF >::value &&
116 std::is_base_of< Fem::HasLocalFunction, GF >::value, int > = 0
117 >
118 inline static GridFunctionView< GF > localFunction ( const GF &gf )
119 {
120 return GridFunctionView< GF >( gf );
121 }
122
123 } // namespace Fem
124
125} // namespace Dune
126
127#endif // #ifndef DUNE_FEMPY_FUNCTION_GRIDFUNCTIONVIEW_HH
discrete function space
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)