DUNE PDELab (git)

interpolatetest.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_FUNCTIONSPACEBASES_TEST_INTERPOLATETEST_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_TEST_INTERPOLATETEST_HH
5
6#include <tuple>
7#include <utility>
8
9#include <dune/common/test/testsuite.hh>
11#include <dune/common/hybridutilities.hh>
12
13#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
14#include <dune/functions/functionspacebases/interpolate.hh>
15
16
17double infinityDiff(const double& x, const double& y)
18{
19 return std::fabs(x-y);
20}
21
22double infinityDiff(const bool& x, const bool& y)
23{
24 return std::fabs(x-y);
25}
26
27template<class X, class Y>
28double infinityDiff(const X& x, const Y& y)
29{
30 if (x.size() != y.size())
31 return false;
32 double diff = 0;
33 Dune::Hybrid::forEach(Dune::range(Dune::Hybrid::size(x)), [&](auto i) {
34 auto&& xi = Dune::Hybrid::elementAt(x, i);
35 auto&& yi = Dune::Hybrid::elementAt(y, i);
36 diff = std::max(diff, infinityDiff(xi, yi));
37 });
38 return diff;
39}
40
41template<class Range, class Basis, class C>
42Dune::TestSuite checkInterpolateConsistency(Basis basis, C&& x)
43{
44 using Coefficients = std::decay_t<C>;
45
46 Dune::TestSuite suite("interpolate consistency check");
47 double coeffTol = 1e-10;
48
49 // generate a discrete function
50 auto fGridFunction = Dune::Functions::makeDiscreteGlobalBasisFunction<Range>(basis, x);
51
52 // Check both functions
53 {
54 const auto& f = fGridFunction;
55 Coefficients y;
57
58 suite.check(infinityDiff(x, y) < coeffTol)
59 << "Interpolation of DiscreteGlobalBasisFunction via local operator() differs from original coefficient vector" << std::endl;
60 }
61
62 return suite;
63}
64
65
66
67#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_TEST_INTERPOLATETEST_HH
A Simple helper class to organize your test suite.
Definition: testsuite.hh:31
Traits for type conversions and type information.
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:256
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:126
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)