dune-pdelab 2.7-git
Loading...
Searching...
No Matches
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
13#include <dune/localfunctions/common/interfaceswitch.hh>
14
19
20namespace 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 // Define types
57 using FESwitch = FiniteElementInterfaceSwitch<typename LFSV::Traits::FiniteElementType>;
58 using BasisSwitch = BasisInterfaceSwitch<typename FESwitch::Basis>;
59 using Range = typename BasisSwitch::Range;
60
61 // Reference to cell
62 const auto& cell = eg.entity();
63
64 // get geometry
65 auto geo = eg.geometry();
66
67 // Initialize vectors outside for loop
68 std::vector<Range> phi(lfsv.size());
69 typename F::Traits::RangeType y(0.0);
70
71 // loop over quadrature points
72 for (const auto& ip : quadratureRule(geo,quadOrder_))
73 {
74 // evaluate shape functions
75 FESwitch::basis(lfsv.finiteElement()).
76 evaluateFunction(ip.position(),phi);
77
78 // evaluate right hand side parameter function
79 f_.evaluate(cell,ip.position(),y);
80
81 // integrate f
82 auto factor = r.weight() * ip.weight() * geo.integrationElement(ip.position());
83 for (size_t i=0; i<lfsv.size(); i++)
84 r.rawAccumulate(lfsv,i,y*phi[i]*factor);
85 }
86 }
87
88 protected:
89 const F& f_;
90
91 // Quadrature rule order
92 unsigned int quadOrder_;
93 };
94
96 } // namespace PDELab
97} // namespace Dune
98
99#endif // DUNE_PDELAB_LOCALOPERATOR_L2VOLUMEFUNCTIONAL_HH
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
QuadratureRuleWrapper< QuadratureRule< typename Geometry::ctype, Geometry::mydimension > > quadratureRule(const Geometry &geo, std::size_t order, QuadratureType::Enum quadrature_type=QuadratureType::GaussLegendre)
Returns a quadrature rule for the given geometry.
Definition: quadraturerules.hh:117
Default flags for all local operators.
Definition: flags.hh:19
A local operator that tests a function against a test function space defined on the entire grid.
Definition: l2volumefunctional.hh:38
L2VolumeFunctional(const F &f, unsigned int quadOrder)
Constructor.
Definition: l2volumefunctional.hh:48
@ doLambdaVolume
Definition: l2volumefunctional.hh:41
unsigned int quadOrder_
Definition: l2volumefunctional.hh:92
void lambda_volume(const EG &eg, const LFSV &lfsv, R &r) const
Definition: l2volumefunctional.hh:54
const F & f_
Definition: l2volumefunctional.hh:89