3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_COMMON_INTEGRATE_HH
25#define DUMUX_COMMON_INTEGRATE_HH
26
27#include <cmath>
28#include <type_traits>
29
30#include <dune/geometry/quadraturerules.hh>
31#include <dune/common/concept.hh>
32
33#if HAVE_DUNE_FUNCTIONS
34#include <dune/functions/gridfunctions/gridfunction.hh>
35#endif
36
40
41namespace Dumux {
42
43// implementation details of the integrate functions
44#ifndef DOXYGEN
45namespace Impl {
46
47struct HasLocalFunction
48{
49 template<class F>
50 auto require(F&& f) -> decltype(
51 localFunction(f),
52 localFunction(f).unbind()
53 );
54};
55
56template<class F>
57static constexpr bool hasLocalFunction()
58{ return Dune::models<HasLocalFunction, F>(); }
59
60template<class Error,
61 typename std::enable_if_t<IsIndexable<Error>::value, int> = 0>
62Error sqrtNorm(const Error& error)
63{
64 using std::sqrt;
65 auto e = error;
66 for (int i = 0; i < error.size(); ++i)
67 e[i] = sqrt(error[i]);
68 return e;
69}
70
71template<class Error,
72 typename std::enable_if_t<!IsIndexable<Error>::value, int> = 0>
73Error sqrtNorm(const Error& error)
74{
75 using std::sqrt;
76 return sqrt(error);
77}
78
79template<class T, typename = int>
80struct FieldTypeImpl
81{
82 using type = T;
83};
84
85template<class T>
86struct FieldTypeImpl<T, typename std::enable_if<(sizeof(std::declval<T>()[0]) > 0), int>::type>
87{
88 using type = typename FieldTypeImpl<std::decay_t<decltype(std::declval<T>()[0])>>::type;
89};
90
91template<class T>
92using FieldType = typename FieldTypeImpl<T>::type;
93
94} // end namespace Impl
95#endif
96
103template<class GridGeometry, class SolutionVector,
104 typename std::enable_if_t<!Impl::hasLocalFunction<SolutionVector>(), int> = 0>
105auto integrateGridFunction(const GridGeometry& gg,
106 const SolutionVector& sol,
107 std::size_t order)
108{
109 using GridView = typename GridGeometry::GridView;
110 using Scalar = typename Impl::FieldType< std::decay_t<decltype(sol[0])> >;
111
112 Scalar integral(0.0);
113 for (const auto& element : elements(gg.gridView()))
114 {
115 const auto elemSol = elementSolution(element, sol, gg);
116 const auto geometry = element.geometry();
117 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
118 for (auto&& qp : quad)
119 {
120 auto value = evalSolution(element, geometry, gg, elemSol, geometry.global(qp.position()));
121 value *= qp.weight()*geometry.integrationElement(qp.position());
122 integral += value;
123 }
124 }
125 return integral;
126}
127
136template<class GridGeometry, class Sol1, class Sol2,
137 typename std::enable_if_t<!Impl::hasLocalFunction<Sol1>(), int> = 0>
138auto integrateL2Error(const GridGeometry& gg,
139 const Sol1& sol1,
140 const Sol2& sol2,
141 std::size_t order)
142{
143 using GridView = typename GridGeometry::GridView;
144 using Scalar = typename Impl::FieldType< std::decay_t<decltype(sol1[0])> >;
145
146 Scalar l2norm(0.0);
147 for (const auto& element : elements(gg.gridView()))
148 {
149 const auto elemSol1 = elementSolution(element, sol1, gg);
150 const auto elemSol2 = elementSolution(element, sol2, gg);
151
152 const auto geometry = element.geometry();
153 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
154 for (auto&& qp : quad)
155 {
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());
161 }
162 }
163
164 using std::sqrt;
165 return sqrt(l2norm);
166}
167
168#if HAVE_DUNE_FUNCTIONS
169
177template<class GridView, class F,
178 typename std::enable_if_t<Impl::hasLocalFunction<F>(), int> = 0>
179auto integrateGridFunction(const GridView& gv,
180 const F& f,
181 std::size_t order)
182{
183 auto fLocal = localFunction(f);
184
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>()))> >;
188
189 Scalar integral(0.0);
190 for (const auto& element : elements(gv))
191 {
192 fLocal.bind(element);
193
194 const auto geometry = element.geometry();
195 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
196 for (auto&& qp : quad)
197 {
198 auto value = fLocal(qp.position());
199 value *= qp.weight()*geometry.integrationElement(qp.position());
200 integral += value;
201 }
202
203 fLocal.unbind();
204 }
205 return integral;
206}
207
217template<class GridView, class F, class G,
218 typename std::enable_if_t<Impl::hasLocalFunction<F>(), int> = 0>
219auto integrateL2Error(const GridView& gv,
220 const F& f,
221 const G& g,
222 std::size_t order)
223{
224 auto fLocal = localFunction(f);
225 auto gLocal = localFunction(g);
226
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>()))> >;
230
231 Scalar l2norm(0.0);
232 for (const auto& element : elements(gv))
233 {
234 fLocal.bind(element);
235 gLocal.bind(element);
236
237 const auto geometry = element.geometry();
238 const auto& quad = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(geometry.type(), order);
239 for (auto&& qp : quad)
240 {
241 const auto error = fLocal(qp.position()) - gLocal(qp.position());
242 l2norm += (error*error)*qp.weight()*geometry.integrationElement(qp.position());
243 }
244
245 gLocal.unbind();
246 fLocal.unbind();
247 }
248
249 using std::sqrt;
250 return sqrt(l2norm);
251}
252#endif
253
254} // end namespace Dumux
255
256#endif
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==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