3 #ifndef DUNE_PDELAB_LOCALOPERATOR_L2VOLUMEFUNCTIONAL_HH
4 #define DUNE_PDELAB_LOCALOPERATOR_L2VOLUMEFUNCTIONAL_HH
8 #include <dune/common/exceptions.hh>
9 #include <dune/common/fvector.hh>
11 #include <dune/geometry/type.hh>
13 #include <dune/localfunctions/common/interfaceswitch.hh>
53 template<
typename EG,
typename LFSV,
typename R>
57 using FESwitch = FiniteElementInterfaceSwitch<typename LFSV::Traits::FiniteElementType>;
58 using BasisSwitch = BasisInterfaceSwitch<typename FESwitch::Basis>;
59 using Range =
typename BasisSwitch::Range;
62 const auto& cell = eg.entity();
65 auto geo = eg.geometry();
68 std::vector<Range> phi(lfsv.size());
69 typename F::Traits::RangeType y(0.0);
75 FESwitch::basis(lfsv.finiteElement()).
76 evaluateFunction(ip.position(),phi);
79 f_.evaluate(cell,ip.position(),y);
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);
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