dune-pdelab  2.4-dev
l2volumefunctional.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_PDELAB_LOCALOPERATOR_L2VOLUMEFUNCTIONAL_HH
4 #define DUNE_PDELAB_LOCALOPERATOR_L2VOLUMEFUNCTIONAL_HH
5 
6 #include <vector>
7 
8 #include <dune/common/exceptions.hh>
9 #include <dune/common/fvector.hh>
10 
11 #include <dune/geometry/type.hh>
12 #include <dune/geometry/quadraturerules.hh>
13 
14 #include <dune/localfunctions/common/interfaceswitch.hh>
15 
19 
20 namespace Dune {
21  namespace PDELab {
25 
35  template<typename F>
38  {
39  public:
40  // residual assembly flags
41  enum { doLambdaVolume = true };
42 
48  L2VolumeFunctional (const F& f, unsigned int quadOrder)
49  : f_(f), quadOrder_(quadOrder)
50  {}
51 
52  // volume integral depending only on test functions
53  template<typename EG, typename LFSV, typename R>
54  void lambda_volume (const EG& eg, const LFSV& lfsv, R& r) const
55  {
56  // domain and range field type
57  typedef FiniteElementInterfaceSwitch<typename LFSV::Traits::FiniteElementType> FESwitch;
58  typedef BasisInterfaceSwitch<typename FESwitch::Basis> BasisSwitch;
59 
60  typedef typename BasisSwitch::DomainField DF;
61  typedef typename BasisSwitch::RangeField RF;
62  typedef typename BasisSwitch::Range Range;
63 
64  // dimensions
65  static const int dimLocal = EG::Geometry::mydimension;
66 
67  // select quadrature rule
68  Dune::GeometryType gt = eg.geometry().type();
69  const Dune::QuadratureRule<DF,dimLocal>& rule =
70  Dune::QuadratureRules<DF,dimLocal>::rule(gt,quadOrder_);
71 
72  // loop over quadrature points
73  for (typename Dune::QuadratureRule<DF,dimLocal>::const_iterator it =
74  rule.begin(); it!=rule.end(); ++it)
75  {
76  // evaluate shape functions
77  std::vector<Range> phi(lfsv.size());
78  FESwitch::basis(lfsv.finiteElement()).
79  evaluateFunction(it->position(),phi);
80 
81  // evaluate right hand side parameter function
82  typename F::Traits::RangeType y(0.0);
83  f_.evaluate(eg.entity(),it->position(),y);
84 
85  // integrate f
86  RF factor = r.weight() * it->weight() * eg.geometry().integrationElement(it->position());
87  for (size_t i=0; i<lfsv.size(); i++)
88  r.rawAccumulate(lfsv,i,y*phi[i]*factor);
89  }
90  }
91 
92  protected:
93  const F& f_;
94 
95  // Quadrature rule order
96  unsigned int quadOrder_;
97  };
98 
100  } // namespace PDELab
101 } // namespace Dune
102 
103 #endif
A local operator that tests a function against a test function space defined on the entire grid...
Definition: l2volumefunctional.hh:36
void lambda_volume(const EG &eg, const LFSV &lfsv, R &r) const
Definition: l2volumefunctional.hh:54
L2VolumeFunctional(const F &f, unsigned int quadOrder)
Constructor.
Definition: l2volumefunctional.hh:48
Default flags for all local operators.
Definition: flags.hh:18
const F & f_
Definition: l2volumefunctional.hh:93
Definition: adaptivity.hh:27
unsigned int quadOrder_
Definition: l2volumefunctional.hh:96
Definition: l2volumefunctional.hh:41
const EG & eg
Definition: constraints.hh:280