24#ifndef DUMUX_COMMON_INTEGRATE_HH
25#define DUMUX_COMMON_INTEGRATE_HH
30#include <dune/common/typetraits.hh>
31#include <dune/geometry/quadraturerules.hh>
32#include <dune/common/concept.hh>
34#if HAVE_DUNE_FUNCTIONS
35#include <dune/functions/gridfunctions/gridfunction.hh>
48struct HasLocalFunction
51 auto require(F&& f) ->
decltype(
53 localFunction(f).unbind()
58static constexpr bool hasLocalFunction()
59{
return Dune::models<HasLocalFunction, F>(); }
62 typename std::enable_if_t<Dune::IsIndexable<Error>::value,
int> = 0>
63Error sqrtNorm(
const Error& error)
67 for (
int i = 0; i < error.size(); ++i)
68 e[i] = sqrt(error[i]);
73 typename std::enable_if_t<!Dune::IsIndexable<Error>::value,
int> = 0>
74Error sqrtNorm(
const Error& error)
80template<
class T,
typename =
int>
87struct FieldTypeImpl<T, typename std::enable_if<(sizeof(std::declval<T>()[0]) > 0), int>::type>
89 using type =
typename FieldTypeImpl<std::decay_t<decltype(std::declval<T>()[0])>>::type;
93using FieldType =
typename FieldTypeImpl<T>::type;
104template<
class GridGeometry,
class SolutionVector,
105 typename std::enable_if_t<!Detail::hasLocalFunction<SolutionVector>(),
int> = 0>
107 const SolutionVector& sol,
110 using GridView =
typename GridGeometry::GridView;
113 Scalar integral(0.0);
114 for (
const auto& element : elements(gg.gridView()))
117 const auto geometry = element.geometry();
118 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
119 for (
auto&& qp : quad)
121 auto value =
evalSolution(element, geometry, gg, elemSol, geometry.global(qp.position()));
122 value *= qp.weight()*geometry.integrationElement(qp.position());
137template<
class GridGeometry,
class Sol1,
class Sol2,
138 typename std::enable_if_t<!Detail::hasLocalFunction<Sol1>(),
int> = 0>
144 using GridView =
typename GridGeometry::GridView;
148 for (
const auto& element : elements(gg.gridView()))
153 const auto geometry = element.geometry();
154 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
155 for (
auto&& qp : quad)
157 const auto& globalPos = geometry.global(qp.position());
158 const auto value1 =
evalSolution(element, geometry, gg, elemSol1, globalPos);
159 const auto value2 =
evalSolution(element, geometry, gg, elemSol2, globalPos);
160 const auto error = (value1 - value2);
161 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
169#if HAVE_DUNE_FUNCTIONS
178template<
class GridView,
class F,
179 typename std::enable_if_t<Detail::hasLocalFunction<F>(),
int> = 0>
184 auto fLocal = localFunction(f);
186 using Element =
typename GridView::template Codim<0>::Entity;
187 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
188 using Scalar =
typename Detail::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
190 Scalar integral(0.0);
191 for (
const auto& element : elements(gv))
193 fLocal.bind(element);
195 const auto geometry = element.geometry();
196 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
197 for (
auto&& qp : quad)
199 auto value = fLocal(qp.position());
200 value *= qp.weight()*geometry.integrationElement(qp.position());
218template<
class GridView,
class F,
class G,
219 typename std::enable_if_t<Detail::hasLocalFunction<F>(),
int> = 0>
225 auto fLocal = localFunction(f);
226 auto gLocal = localFunction(g);
228 using Element =
typename GridView::template Codim<0>::Entity;
229 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
230 using Scalar =
typename Detail::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
233 for (
const auto& element : elements(gv))
235 fLocal.bind(element);
236 gLocal.bind(element);
238 const auto geometry =
element.geometry();
239 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
240 for (
auto&& qp : quad)
242 const auto error = fLocal(qp.position()) - gLocal(qp.position());
243 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
263template<
class Scalar,
class Function,
264 typename std::enable_if_t<std::is_invocable_r_v<Scalar, Function, Scalar>>...>
266 const Scalar lowerBound,
267 const Scalar upperBound,
268 const Scalar targetAbsoluteError = 1e-13)
A double exponential integrator.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
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==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
FieldType
Identifier for vtk field types.
Definition: fieldtype.hh:34
Scalar integrateScalarFunction(const Function &f, const Scalar lowerBound, const Scalar upperBound, const Scalar targetAbsoluteError=1e-13)
Integrate a scalar function.
Definition: integrate.hh:265
auto integrateGridFunction(const GridGeometry &gg, const SolutionVector &sol, std::size_t order)
Integrate a grid function over a grid view.
Definition: integrate.hh:106
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:139
static Scalar integrate(const Function &f, const Scalar a, const Scalar b, const Scalar targetAbsoluteError, int &numFunctionEvaluations, Scalar &errorEstimate)
Integrate an analytic function over a finite interval.
Definition: doubleexpintegrator.hh:86