DUNE PDELab (2.8)

test-fe.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4// This header is not part of the official Dune API and might be subject
5// to change. You can use this header to test external finite element
6// implementations, but be warned that your tests might break with future
7// Dune versions.
8
9#ifndef DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
10#define DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
11
12#include <algorithm>
13#include <cmath>
14#include <cstddef>
15#include <cstdlib>
16#include <iomanip>
17#include <iostream>
18#include <ostream>
19#include <vector>
20
23
25
26// This class defines a local finite element function.
27// It is determined by a local finite element and
28// representing the local basis and a coefficient vector.
29// This provides the evaluate method needed by the interpolate()
30// method.
31template<class FE>
32class FEFunction
33{
34 const FE& fe;
35
36public:
37 using RangeType = typename FE::Traits::Basis::Traits::Range;
38 using DomainType = typename FE::Traits::Basis::Traits::DomainLocal;
39
40 struct Traits {
41 using RangeType = typename FE::Traits::Basis::Traits::Range;
42 using DomainType = typename FE::Traits::Basis::Traits::DomainLocal;
43 };
44
45 typedef typename FE::Traits::Basis::Traits::RangeField CT;
46
47 std::vector<CT> coeff;
48
49 FEFunction(const FE& fe_) : fe(fe_) { resetCoefficients(); }
50
51 void resetCoefficients() {
52 coeff.resize(fe.basis().size());
53 for(std::size_t i=0; i<coeff.size(); ++i)
54 coeff[i] = 0;
55 }
56
57 void setRandom(double max) {
58 coeff.resize(fe.basis().size());
59 for(std::size_t i=0; i<coeff.size(); ++i)
60 coeff[i] = ((1.0*std::rand()) / RAND_MAX - 0.5)*2.0*max;
61 }
62
63 void evaluate (const DomainType& x, RangeType& y) const {
64 std::vector<RangeType> yy;
65 fe.basis().evaluateFunction(x, yy);
66 y = 0.0;
67 for (std::size_t i=0; i<yy.size(); ++i)
68 y.axpy(coeff[i], yy[i]);
69 }
70};
71
72
73// Check if interpolation is consistent with basis evaluation.
85template<class FE>
86bool testInterpolation(const FE& fe, double eps, int n=5)
87{
88 bool success = true;
89 FEFunction<FE> f(fe);
90
91 std::vector<typename FEFunction<FE>::CT> coeff;
92 std::vector<typename FEFunction<FE>::CT> coeff2;
93 for(int i=0; i<n && success; ++i) {
94 // Set random coefficient vector
95 f.setRandom(100);
96
97 // Compute interpolation weights
98
100 // Part A: Feed the function to the 'interpolate' method in form of
101 // a class providing an evaluate() method.
102 // This way is deprecated since dune-localfunctions 2.7.
104 fe.interpolation().interpolate(f, coeff);
105
107 // Part B: Redo the same test, but feed the function to the
108 // 'interpolate' method in form of a callable.
110 auto callableF = [&](const auto& x) {
111 using Range = typename FE::Traits::Basis::Traits::Range;
112 Range y(0);
113 f.evaluate(x,y);
114 return y;
115 };
116 fe.interpolation().interpolate(callableF, coeff2);
117 if (coeff != coeff2) {
118 std::cout << "Passing callable and function with evaluate yields different "
119 << "interpolation weights for finite element type "
120 << Dune::className<FE>() << ":" << std::endl;
121 success = false;
122 }
123
124 // Check size of weight vector
125 if (coeff.size() != fe.basis().size()) {
126 std::cout << "Bug in LocalInterpolation for finite element type "
127 << Dune::className<FE>() << ":" << std::endl;
128 std::cout << " Interpolation vector has size " << coeff.size()
129 << std::endl;
130 std::cout << " Basis has size " << fe.basis().size() << std::endl;
131 std::cout << std::endl;
132 success = false;
133
134 // skip rest of loop since that depends on matching sizes
135 continue;
136 }
137
138 // Check if interpolation weights are equal to coefficients
139 for(std::size_t j=0; j<coeff.size() && success; ++j) {
140 if ( std::abs(coeff[j]-f.coeff[j]) >
141 eps*(std::max(std::abs(f.coeff[j]), 1.0)) )
142 {
143 std::cout << std::setprecision(16);
144 std::cout << "Bug in LocalInterpolation for finite element type "
145 << Dune::className<FE>() << ":" << std::endl;
146 std::cout << " Interpolation weight " << j << " differs by "
147 << std::abs(coeff[j]-f.coeff[j]) << " from coefficient of "
148 << "linear combination." << std::endl;
149 std::cout << std::endl;
150 success = false;
151 }
152 }
153 }
154 return success;
155}
156
157// check whether Jacobian agrees with FD approximation
169template<class Geo, class FE>
170bool testJacobian(const Geo &geo, const FE& fe, double eps, double delta,
171 std::size_t order = 2)
172{
173 typedef typename FE::Traits::Basis Basis;
174
175 typedef typename Basis::Traits::DomainField DF;
176 static const std::size_t dimDLocal = Basis::Traits::dimDomainLocal;
177 typedef typename Basis::Traits::DomainLocal DomainLocal;
178 static const std::size_t dimDGlobal = Basis::Traits::dimDomainGlobal;
179
180 static const std::size_t dimR = Basis::Traits::dimRange;
181 typedef typename Basis::Traits::Range Range;
182
183 typedef typename Basis::Traits::Jacobian Jacobian;
184
185 bool success = true;
186
187 // ////////////////////////////////////////////////////////////
188 // Check the partial derivatives by comparing them
189 // to finite difference approximations
190 // ////////////////////////////////////////////////////////////
191
192 // A set of test points
195
196 // Loop over all quadrature points
197 for (std::size_t i=0; i < quad.size(); i++) {
198
199 // Get a test point
200 const DomainLocal& testPoint = quad[i].position();
201
202 // Get the shape function derivatives there
203 std::vector<Jacobian> jacobians;
204 fe.basis().evaluateJacobian(testPoint, jacobians);
205 if(jacobians.size() != fe.basis().size()) {
206 std::cout << "Bug in evaluateJacobianGlobal() for finite element type "
207 << Dune::className<FE>() << ":" << std::endl;
208 std::cout << " Jacobian vector has size " << jacobians.size()
209 << std::endl;
210 std::cout << " Basis has size " << fe.basis().size() << std::endl;
211 std::cout << std::endl;
212 return false;
213 }
214
216 geo.jacobianTransposed(testPoint);
217
218 // Loop over all shape functions in this set
219 for (std::size_t j=0; j<fe.basis().size(); ++j) {
220
221 // basis.evaluateJacobian returns global derivatives, however we can
222 // only do local derivatives, so transform the derivatives back into
223 // local coordinates
225 for(std::size_t k = 0; k < dimR; ++k)
226 for(std::size_t l = 0; l < dimDGlobal; ++l)
227 for(std::size_t m = 0; m < dimDLocal; ++m)
228 localJacobian[k][m] += jacobians[j][k][l] * geoJT[m][l];
229
230 // Loop over all local directions
231 for (std::size_t m = 0; m < dimDLocal; ++m) {
232
233 // Compute an approximation to the derivative by finite differences
234 DomainLocal upPos = testPoint;
235 DomainLocal downPos = testPoint;
236
237 upPos[m] += delta;
238 downPos[m] -= delta;
239
240 std::vector<Range> upValues, downValues;
241
242 fe.basis().evaluateFunction(upPos, upValues);
243 fe.basis().evaluateFunction(downPos, downValues);
244
245 //Loop over all components
246 for(std::size_t k = 0; k < dimR; ++k) {
247
248 // The current partial derivative, just for ease of notation
249 double derivative = localJacobian[k][m];
250
251 double finiteDiff = (upValues[j][k] - downValues[j][k]) / (2*delta);
252
253 // Check
254 if ( std::abs(derivative-finiteDiff) >
255 eps/delta*(std::max(std::abs(finiteDiff), 1.0)) )
256 {
257 std::cout << std::setprecision(16);
258 std::cout << "Bug in evaluateJacobian() for finite element type "
259 << Dune::className<FE>() << ":" << std::endl;
260 std::cout << " Shape function derivative does not agree with "
261 << "FD approximation" << std::endl;
262 std::cout << " Shape function " << j << " component " << k
263 << " at position " << testPoint << ": derivative in "
264 << "local direction " << m << " is "
265 << derivative << ", but " << finiteDiff << " is "
266 << "expected." << std::endl;
267 std::cout << std::endl;
268 success = false;
269 }
270 } //Loop over all components
271 } // Loop over all local directions
272 } // Loop over all shape functions in this set
273 } // Loop over all quadrature points
274
275 return success;
276}
277
278// call tests for given finite element
279template<class Geo, class FE>
280bool testFE(const Geo &geo, const FE& fe, double eps, double delta,
281 unsigned order = 2)
282{
283 bool success = true;
284
285 success = testInterpolation(fe, eps) and success;
286 success = testJacobian(geo, fe, eps, delta, order) and success;
287
288 return success;
289}
290
291#endif // DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
A dense n x m matrix.
Definition: fmatrix.hh:69
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:152
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:280
A free function to provide the demangled class name of a given object or type as a string.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:39
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)