12#ifndef DUMUX_COMMON_INTEGRATE_HH
13#define DUMUX_COMMON_INTEGRATE_HH
18#include <dune/common/typetraits.hh>
19#include <dune/geometry/quadraturerules.hh>
20#include <dune/common/concept.hh>
22#if HAVE_DUNE_FUNCTIONS
23#include <dune/functions/gridfunctions/gridfunction.hh>
36struct HasLocalFunction
39 auto require(F&& f) ->
decltype(
41 localFunction(f).unbind()
46static constexpr bool hasLocalFunction()
47{
return Dune::models<HasLocalFunction, F>(); }
50 typename std::enable_if_t<Dune::IsIndexable<Error>::value,
int> = 0>
51Error sqrtNorm(
const Error& error)
55 for (
int i = 0; i < error.size(); ++i)
56 e[i] = sqrt(error[i]);
61 typename std::enable_if_t<!Dune::IsIndexable<Error>::value,
int> = 0>
62Error sqrtNorm(
const Error& error)
68template<
class T,
typename =
int>
75struct FieldTypeImpl<T, typename std::enable_if<(sizeof(std::declval<T>()[0]) > 0), int>::type>
77 using type =
typename FieldTypeImpl<std::decay_t<decltype(std::declval<T>()[0])>>::type;
81using FieldType =
typename FieldTypeImpl<T>::type;
92template<
class GridGeometry,
class SolutionVector,
93 typename std::enable_if_t<!Detail::hasLocalFunction<SolutionVector>(),
int> = 0>
95 const SolutionVector& sol,
98 using GridView =
typename GridGeometry::GridView;
101 Scalar integral(0.0);
102 for (
const auto& element : elements(gg.gridView()))
105 const auto geometry = element.geometry();
106 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
107 for (
auto&& qp : quad)
109 auto value =
evalSolution(element, geometry, gg, elemSol, geometry.global(qp.position()));
110 value *= qp.weight()*geometry.integrationElement(qp.position());
125template<
class GridGeometry,
class Sol1,
class Sol2,
126 typename std::enable_if_t<!Detail::hasLocalFunction<Sol1>(),
int> = 0>
132 using GridView =
typename GridGeometry::GridView;
136 for (
const auto& element : elements(gg.gridView()))
141 const auto geometry = element.geometry();
142 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
143 for (
auto&& qp : quad)
145 const auto& globalPos = geometry.global(qp.position());
146 const auto value1 =
evalSolution(element, geometry, gg, elemSol1, globalPos);
147 const auto value2 =
evalSolution(element, geometry, gg, elemSol2, globalPos);
148 const auto error = (value1 - value2);
149 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
157#if HAVE_DUNE_FUNCTIONS
166template<
class GridView,
class F,
167 typename std::enable_if_t<Detail::hasLocalFunction<F>(),
int> = 0>
172 auto fLocal = localFunction(f);
174 using Element =
typename GridView::template Codim<0>::Entity;
175 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
176 using Scalar =
typename Detail::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
178 Scalar integral(0.0);
179 for (
const auto& element : elements(gv))
181 fLocal.bind(element);
183 const auto geometry = element.geometry();
184 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
185 for (
auto&& qp : quad)
187 auto value = fLocal(qp.position());
188 value *= qp.weight()*geometry.integrationElement(qp.position());
206template<
class GridView,
class F,
class G,
207 typename std::enable_if_t<Detail::hasLocalFunction<F>(),
int> = 0>
213 auto fLocal = localFunction(f);
214 auto gLocal = localFunction(g);
216 using Element =
typename GridView::template Codim<0>::Entity;
217 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
218 using Scalar =
typename Detail::FieldType< std::decay_t<decltype(fLocal(std::declval<LocalPosition>()))> >;
221 for (
const auto& element : elements(gv))
223 fLocal.bind(element);
224 gLocal.bind(element);
226 const auto geometry =
element.geometry();
227 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
228 for (
auto&& qp : quad)
230 const auto error = fLocal(qp.position()) - gLocal(qp.position());
231 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
251template<
class Scalar,
class Function,
252 typename std::enable_if_t<std::is_invocable_r_v<Scalar, Function, Scalar>>...>
254 const Scalar lowerBound,
255 const Scalar upperBound,
256 const Scalar targetAbsoluteError = 1e-13)
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:80
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 CVFEElementSolution< 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:152
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
Scalar integrateScalarFunction(const Function &f, const Scalar lowerBound, const Scalar upperBound, const Scalar targetAbsoluteError=1e-13)
Integrate a scalar function.
Definition: integrate.hh:253
auto integrateGridFunction(const GridGeometry &gg, const SolutionVector &sol, std::size_t order)
Integrate a grid function over a grid view.
Definition: integrate.hh:94
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:127