version 3.10-dev
integrate.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_COMMON_INTEGRATE_HH
13#define DUMUX_COMMON_INTEGRATE_HH
14
15#include <cmath>
16#include <type_traits>
17
18#include <dune/common/typetraits.hh>
19#include <dune/geometry/quadraturerules.hh>
20#include <dune/common/concept.hh>
21
22#if HAVE_DUNE_FUNCTIONS
23#include <dune/functions/gridfunctions/gridfunction.hh>
24#endif
25
29
30namespace Dumux {
31
32// implementation details of the integrate functions
33#ifndef DOXYGEN
34namespace Detail {
35
36struct HasLocalFunction
37{
38 template<class F>
39 auto require(F&& f) -> decltype(
40 localFunction(f),
41 localFunction(f).unbind()
42 );
43};
44
45template<class F>
46static constexpr bool hasLocalFunction()
47{ return Dune::models<HasLocalFunction, F>(); }
48
49template<class Error,
50 typename std::enable_if_t<Dune::IsIndexable<Error>::value, int> = 0>
51Error sqrtNorm(const Error& error)
52{
53 using std::sqrt;
54 auto e = error;
55 for (int i = 0; i < error.size(); ++i)
56 e[i] = sqrt(error[i]);
57 return e;
58}
59
60template<class Error,
61 typename std::enable_if_t<!Dune::IsIndexable<Error>::value, int> = 0>
62Error sqrtNorm(const Error& error)
63{
64 using std::sqrt;
65 return sqrt(error);
66}
67
68template<class T, typename = int>
69struct FieldTypeImpl
70{
71 using type = T;
72};
73
74template<class T>
75struct FieldTypeImpl<T, typename std::enable_if<(sizeof(std::declval<T>()[0]) > 0), int>::type>
76{
77 using type = typename FieldTypeImpl<std::decay_t<decltype(std::declval<T>()[0])>>::type;
78};
79
80template<class T>
81using FieldType = typename FieldTypeImpl<T>::type;
82
83} // end namespace Detail
84#endif
85
92template<class GridGeometry, class SolutionVector,
93 typename std::enable_if_t<!Detail::hasLocalFunction<SolutionVector>(), int> = 0>
94auto integrateGridFunction(const GridGeometry& gg,
95 const SolutionVector& sol,
96 std::size_t order)
97{
98 using GridView = typename GridGeometry::GridView;
99 using Scalar = typename Detail::FieldType< std::decay_t<decltype(sol[0])> >;
100
101 Scalar integral(0.0);
102 for (const auto& element : elements(gg.gridView()))
103 {
104 const auto elemSol = elementSolution(element, sol, gg);
105 const auto geometry = element.geometry();
106 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
107 for (auto&& qp : quad)
108 {
109 auto value = evalSolution(element, geometry, gg, elemSol, geometry.global(qp.position()));
110 value *= qp.weight()*geometry.integrationElement(qp.position());
111 integral += value;
112 }
113 }
114 return integral;
115}
116
125template<class GridGeometry, class Sol1, class Sol2,
126 typename std::enable_if_t<!Detail::hasLocalFunction<Sol1>(), int> = 0>
127auto integrateL2Error(const GridGeometry& gg,
128 const Sol1& sol1,
129 const Sol2& sol2,
130 std::size_t order)
131{
132 using GridView = typename GridGeometry::GridView;
133 using Scalar = typename Detail::FieldType< std::decay_t<decltype(sol1[0])> >;
134
135 Scalar l2norm(0.0);
136 for (const auto& element : elements(gg.gridView()))
137 {
138 const auto elemSol1 = elementSolution(element, sol1, gg);
139 const auto elemSol2 = elementSolution(element, sol2, gg);
140
141 const auto geometry = element.geometry();
142 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
143 for (auto&& qp : quad)
144 {
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());
150 }
151 }
152
153 using std::sqrt;
154 return sqrt(l2norm);
155}
156
157#if HAVE_DUNE_FUNCTIONS
158
166template<class GridView, class F,
167 typename std::enable_if_t<Detail::hasLocalFunction<F>(), int> = 0>
168auto integrateGridFunction(const GridView& gv,
169 const F& f,
170 std::size_t order)
171{
172 auto fLocal = localFunction(f);
173
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>()))> >;
177
178 Scalar integral(0.0);
179 for (const auto& element : elements(gv))
180 {
181 fLocal.bind(element);
182
183 const auto geometry = element.geometry();
184 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
185 for (auto&& qp : quad)
186 {
187 auto value = fLocal(qp.position());
188 value *= qp.weight()*geometry.integrationElement(qp.position());
189 integral += value;
190 }
191
192 fLocal.unbind();
193 }
194 return integral;
195}
196
206template<class GridView, class F, class G,
207 typename std::enable_if_t<Detail::hasLocalFunction<F>(), int> = 0>
208auto integrateL2Error(const GridView& gv,
209 const F& f,
210 const G& g,
211 std::size_t order)
212{
213 auto fLocal = localFunction(f);
214 auto gLocal = localFunction(g);
215
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>()))> >;
219
220 Scalar l2norm(0.0);
221 for (const auto& element : elements(gv))
222 {
223 fLocal.bind(element);
224 gLocal.bind(element);
225
226 const auto geometry = element.geometry();
227 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
228 for (auto&& qp : quad)
229 {
230 const auto error = fLocal(qp.position()) - gLocal(qp.position());
231 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
232 }
233
234 gLocal.unbind();
235 fLocal.unbind();
236 }
237
238 using std::sqrt;
239 return sqrt(l2norm);
240}
241#endif
242
251template<class Scalar, class Function,
252 typename std::enable_if_t<std::is_invocable_r_v<Scalar, Function, Scalar>>...>
253Scalar integrateScalarFunction(const Function& f,
254 const Scalar lowerBound,
255 const Scalar upperBound,
256 const Scalar targetAbsoluteError = 1e-13)
257{
258 return DoubleExponentialIntegrator<Scalar>::integrate(f, lowerBound, upperBound, targetAbsoluteError);
259}
260
261} // end namespace Dumux
262
263#endif
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
FieldType
Identifier for vtk field types.
Definition: fieldtype.hh:22
Definition: adapt.hh:17
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