24#ifndef DUMUX_COMMON_INTEGRATE_HH
25#define DUMUX_COMMON_INTEGRATE_HH
30#include <dune/geometry/quadraturerules.hh>
31#include <dune/common/concept.hh>
33#if HAVE_DUNE_FUNCTIONS
34#include <dune/functions/gridfunctions/gridfunction.hh>
47struct HasLocalFunction
50 auto require(F&& f) ->
decltype(
52 localFunction(f).unbind()
57static constexpr bool hasLocalFunction()
58{
return Dune::models<HasLocalFunction, F>(); }
61 typename std::enable_if_t<IsIndexable<Error>::value,
int> = 0>
62Error sqrtNorm(
const Error& error)
66 for (
int i = 0; i < error.size(); ++i)
67 e[i] = sqrt(error[i]);
72 typename std::enable_if_t<!IsIndexable<Error>::value,
int> = 0>
73Error sqrtNorm(
const Error& error)
79template<
class T,
typename =
int>
86struct FieldTypeImpl<T, typename std::enable_if<(sizeof(std::declval<T>()[0]) > 0), int>::type>
88 using type =
typename FieldTypeImpl<std::decay_t<decltype(std::declval<T>()[0])>>::type;
92using FieldType =
typename FieldTypeImpl<T>::type;
103template<
class GridGeometry,
class SolutionVector,
104 typename std::enable_if_t<!Impl::hasLocalFunction<SolutionVector>(),
int> = 0>
106 const SolutionVector& sol,
109 using GridView =
typename GridGeometry::GridView;
110 using Scalar =
typename Impl::FieldType< std::decay_t<
decltype(sol[0])> >;
112 Scalar integral(0.0);
113 for (
const auto& element : elements(gg.gridView()))
116 const auto geometry = element.geometry();
117 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
118 for (
auto&& qp : quad)
120 auto value =
evalSolution(element, geometry, gg, elemSol, geometry.global(qp.position()));
121 value *= qp.weight()*geometry.integrationElement(qp.position());
136template<
class GridGeometry,
class Sol1,
class Sol2,
137 typename std::enable_if_t<!Impl::hasLocalFunction<Sol1>(),
int> = 0>
143 using GridView =
typename GridGeometry::GridView;
144 using Scalar =
typename Impl::FieldType< std::decay_t<
decltype(sol1[0])> >;
147 for (
const auto& element : elements(gg.gridView()))
152 const auto geometry = element.geometry();
153 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
154 for (
auto&& qp : quad)
156 const auto& globalPos = geometry.global(qp.position());
157 const auto value1 =
evalSolution(element, geometry, gg, elemSol1, globalPos);
158 const auto value2 =
evalSolution(element, geometry, gg, elemSol2, globalPos);
159 const auto error = (value1 - value2);
160 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
168#if HAVE_DUNE_FUNCTIONS
177template<
class GridView,
class F,
178 typename std::enable_if_t<Impl::hasLocalFunction<F>(),
int> = 0>
183 auto fLocal = localFunction(f);
185 using Element =
typename GridView::template Codim<0>::Entity;
186 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
187 using Scalar =
typename Impl::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
189 Scalar integral(0.0);
190 for (
const auto& element : elements(gv))
192 fLocal.bind(element);
194 const auto geometry = element.geometry();
195 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
196 for (
auto&& qp : quad)
198 auto value = fLocal(qp.position());
199 value *= qp.weight()*geometry.integrationElement(qp.position());
217template<
class GridView,
class F,
class G,
218 typename std::enable_if_t<Impl::hasLocalFunction<F>(),
int> = 0>
224 auto fLocal = localFunction(f);
225 auto gLocal = localFunction(g);
227 using Element =
typename GridView::template Codim<0>::Entity;
228 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
229 using Scalar =
typename Impl::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
232 for (
const auto& element : elements(gv))
234 fLocal.bind(element);
235 gLocal.bind(element);
237 const auto geometry = element.geometry();
238 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
239 for (
auto&& qp : quad)
241 const auto error = fLocal(qp.position()) - gLocal(qp.position());
242 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
free functions for the evaluation of primary variables inside elements.
Element solution classes and factory functions.
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:95
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
auto integrateGridFunction(const GridGeometry &gg, const SolutionVector &sol, std::size_t order)
Integrate a grid function over a grid view.
Definition: integrate.hh:105
auto integrateL2Error(const GridGeometry &gg, const Sol1 &sol1, const Sol2 &sol2, std::size_t order)
Integrate a function over a grid view.
Definition: integrate.hh:138