Dune Core Modules (2.9.1)

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