Dune Core Modules (2.9.0)

gridfunctiontest.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_TEST_GRIDFUNCTIONTEST_HH
4#define DUNE_FUNCTIONS_GRIDFUNCTIONS_TEST_GRIDFUNCTIONTEST_HH
5
6#include <type_traits>
7
9
10#include <dune/functions/common/functionconcepts.hh>
11#include <dune/functions/gridfunctions/gridfunction.hh>
12#include <dune/functions/gridfunctions/gridviewfunction.hh>
13#include <dune/functions/gridfunctions/gridviewentityset.hh>
14
15
16namespace Dune {
17namespace Functions {
18namespace Test {
19
20
21
22bool checkTrue(bool value, std::string error)
23{
24 if (not(value))
25 std::cout << "TEST FAILURE:" << error << std::endl;
26 return value;
27}
28
29
30
31template<class GridView, class F>
32double integrateGridViewFunction(const GridView& gridView, const F& f, unsigned int quadOrder)
33{
34 static const int dim = GridView::dimension;
35
36 double integral = 0;
37
38 auto fLocal = localFunction(f);
39
40 // Loop over elements and integrate over the function
41 for (const auto& e : elements(gridView))
42 {
43 auto geometry = e.geometry();
44
45 fLocal.bind(e);
46
47 if (not fLocal.bound())
48 DUNE_THROW(Exception, "Function is not bound after bind()");
49
50 // A quadrature rule
51 const auto& quad = QuadratureRules<double, dim>::rule(e.type(), quadOrder);
52
53 // Loop over all quadrature points
54 for ( size_t pt=0; pt < quad.size(); pt++ ) {
55
56 // Position of the current quadrature point in the reference element
57 auto quadPos = quad[pt].position();
58
59 // The multiplicative factor in the integral transformation formula
60 auto integrationElement = geometry.integrationElement(quadPos);
61
62 integral += fLocal(quadPos) * quad[pt].weight() * integrationElement;
63 }
64 fLocal.unbind();
65 }
66 return integral;
67}
68
69
70template<class GridView, class F>
71bool checkGridViewFunction(const GridView& gridView, const F& f, double exactIntegral, unsigned int quadOrder=1)
72{
73 bool passed = true;
74 double integral;
75
77
78 using EntitySet = typename F::EntitySet;
79 using Domain = typename EntitySet::GlobalCoordinate;
80 using Range = typename std::invoke_result<F, Domain>::type;
81
82
83
84 std::cout << "Checking raw function f on grid view" << std::endl;
85 passed = checkTrue(isGridFunction<decltype(f), Range(Domain), EntitySet>(), "Function does not model GridFunction concept");
86 integral = integrateGridViewFunction(gridView, f, quadOrder);
87 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
88
89 std::cout << "Checking GridFunction<Range(Domain), EntitySet>(f) on grid view" << std::endl;
90 GridFunction<Range(Domain), EntitySet> f2 = f;
91 passed = checkTrue(isGridFunction<decltype(f2), Range(Domain), EntitySet>(), "Function does not model GridFunction concept");
92 integral = integrateGridViewFunction(gridView, f2, quadOrder);
93 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
94
95 std::cout << "Checking GridViewFunction<Range(Domain), GridView>(f) on grid view" << std::endl;
96 GridViewFunction<Range(Domain), GridView> f3 = f;
97 passed = checkTrue(isGridFunction<decltype(f3), Range(Domain), EntitySet>(), "Function does not model GridFunction concept");
98 integral = integrateGridViewFunction(gridView, f3, quadOrder);
99 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
100
101 std::cout << "Checking makeGridFunction(f) on grid view" << std::endl;
102 auto f4 = makeGridViewFunction(f, gridView);
103 passed = checkTrue(isGridFunction<decltype(f4), Range(Domain), EntitySet>(), "Function does not model GridFunction concept");
104 integral = integrateGridViewFunction(gridView, f4, quadOrder);
105 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
106
107 std::cout << "Checking GridViewFunction<Range(Domain), GridView>(makeGridFunction(f)) on grid view" << std::endl;
108 GridViewFunction<Range(Domain), GridView> f5 = f4;
109 integral = integrateGridViewFunction(gridView, f5, quadOrder);
110 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
111
112 std::cout << "Checking default constructed and assigned GridViewFunction<Range(Domain), GridView> on grid view" << std::endl;
113 GridViewFunction<Range(Domain), GridView> f6;
114 f6 = f5;
115 integral = integrateGridViewFunction(gridView, f6, quadOrder);
116 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
117
118 std::cout << "Checking reassigned GridViewFunction<Range(Domain), GridView> on grid view" << std::endl;
119 f6 = f3;
120 integral = integrateGridViewFunction(gridView, f6, quadOrder);
121 passed = checkTrue(std::abs(integral-exactIntegral) < 1e-10, "Integral is "+std::to_string(integral)+" but should be "+std::to_string(exactIntegral));
122
123 return passed;
124}
125
126
127
128
129} // namespace Test
130} // namespace Functions
131} // namespace Dune
132
133
134
135
136#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_TEST_GRIDFUNCTIONTEST_HH
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:266
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
static constexpr bool isGridFunction()
Check if F models the GridFunction concept with given signature and entity set.
Definition: functionconcepts.hh:269
static constexpr int dimension
The dimension of the grid.
Definition: gridview.hh:148
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)