dune-fem  2.4.1-rc
dgl2projection.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_DGL2PROJECTION_HH
2 #define DUNE_FEM_DGL2PROJECTION_HH
3 
9 
10 namespace Dune
11 {
12 
13  namespace Fem
14  {
15 
16  // DGL2ProjectionImpl
17  // ------------------
18 
19  // implementation of L2 projection for discontinuous spaces
21  {
22  template <int dummy, bool hasLocalFunction>
23  struct ProjectChooser
24  {
25  template <class FunctionImp, class FunctionSpace>
27  {
28  const FunctionImp& function_;
29  public:
33  FunctionAdapter(const FunctionImp& f) : function_(f) {}
34 
35  void evaluate(const DomainType& local,
36  RangeType& ret) const
37  {
38  function_.evaluate( local , ret );
39  }
40  };
41 
42  template <class FunctionImp, class DiscreteFunctionImp>
43  static void project(const FunctionImp& f,
44  DiscreteFunctionImp& discFunc,
45  const int polOrd)
46  {
47  // some typedefs
48  typedef typename DiscreteFunctionImp :: DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
49  typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
51  // create function adapter in case of incorrect implementation
52  FunctionAdapterType af( f );
53  // create discrete function adapter
55  "DGL2projection::adapter" , f , discFunc.space().gridPart());
56  DGL2ProjectionImpl::projectFunction(adapter, discFunc, polOrd);
57  }
58  };
59 
60  template <int dummy>
61  struct ProjectChooser<dummy,true>
62  {
63  template <class FunctionImp, class DiscreteFunctionImp>
64  static void project(const FunctionImp& f,
65  DiscreteFunctionImp& discFunc,
66  const int polOrd )
67  {
68  DGL2ProjectionImpl::projectFunction(f, discFunc, polOrd);
69  }
70  };
71 
72  public:
80  template <class FunctionImp, class DiscreteFunctionImp>
81  static void project(const FunctionImp& f, DiscreteFunctionImp& discFunc,
82  const int quadOrd = -1, const bool communicate = true )
83  {
84  ProjectChooser<0, Conversion<FunctionImp, HasLocalFunction> ::exists > :: project(f,discFunc,quadOrd);
85 
86  // do communication in parallel cases
87  if( communicate )
88  discFunc.communicate();
89  }
90 
92  //- make interface equal to LagrangeInterpolation
93  template <class FunctionImp, class DiscreteFunctionImp>
94  static void apply(const FunctionImp& f, DiscreteFunctionImp& discFunc )
95  {
96  project( f, discFunc );
97  }
98 
99  protected:
100  template <class FunctionImp, class DiscreteFunctionImp>
101  static void projectFunction(const FunctionImp& func,
102  DiscreteFunctionImp& discFunc,
103  int polOrd = -1)
104  {
105  typedef typename DiscreteFunctionImp::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
106  typedef typename DiscreteFunctionImp::LocalFunctionType LocalFuncType;
107  typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
108  typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
109 
110  typedef typename FunctionImp::LocalFunctionType LocalFType;
111 
112  const DiscreteFunctionSpaceType& space = discFunc.space();
113 
114  // type of quadrature
115  typedef CachingQuadrature<GridPartType,0> QuadratureType;
116  // type of local mass matrix
118 
119  const int quadOrd = (polOrd == -1) ? (2 * space.order()) : polOrd;
120 
121  // create local mass matrix object
122  LocalMassMatrixType massMatrix( space, quadOrd );
123 
124  // clear destination
125  discFunc.clear();
126 
127  // extract type from grid part
128  typedef typename GridPartType::template Codim<0>::GeometryType Geometry;
129 
130  // create storage for values
131  std::vector< RangeType > values;
132 
133  for(const auto & en : space)
134  {
135  // get geometry
136  const Geometry& geo = en.geometry();
137 
138  // get quadrature
139  QuadratureType quad(en, quadOrd);
140 
141  // get local function of destination
142  LocalFuncType lf = discFunc.localFunction(en);
143 
144  // get local function of argument
145  const LocalFType f = func.localFunction(en);
146 
147  const int quadNop = quad.nop();
148 
149  // adjust size of values
150  values.resize( quadNop );
151 
152  for(int qP = 0; qP < quadNop ; ++qP)
153  {
154  const double intel =
155  quad.weight(qP) * geo.integrationElement( quad.point(qP) );
156 
157  // evaluate function at quadrature point
158  f.evaluate(quad[ qP ], values[ qP ] );
159 
160  // apply weight
161  values[ qP ] *= intel;
162  }
163 
164  // add values to local function
165  lf.axpyQuadrature( quad, values );
166 
167  // apply inverse of mass matrix to local function
168  massMatrix.applyInverse( en, lf );
169  }
170  }
171  };
172 
173  } // namespace Fem
174 
175 } // namespace Dune
176 
177 #endif // #ifndef DUNE_FEM_DGL2PROJECTION_HH
FunctionAdapter(const FunctionImp &f)
Definition: dgl2projection.hh:33
A vector valued function space.
Definition: functionspace.hh:16
void evaluate(const DomainType &local, RangeType &ret) const
Definition: dgl2projection.hh:35
FunctionSpaceType::RangeType RangeType
Definition: dgl2projection.hh:31
static void apply(const FunctionImp &f, DiscreteFunctionImp &discFunc)
project function to discrete space
Definition: dgl2projection.hh:94
FunctionSpace FunctionSpaceType
Definition: dgl2projection.hh:30
static void project(const FunctionImp &f, DiscreteFunctionImp &discFunc, const int quadOrd=-1, const bool communicate=true)
Definition: dgl2projection.hh:81
Definition: coordinate.hh:4
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:66
FunctionSpaceType::DomainType DomainType
Definition: dgl2projection.hh:32
Definition: dgl2projection.hh:20
GridFunctionAdapter provides local functions for a Function.
Definition: gridfunctionadapter.hh:36
static void projectFunction(const FunctionImp &func, DiscreteFunctionImp &discFunc, int polOrd=-1)
Definition: dgl2projection.hh:101
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:70