Dune Core Modules (2.6.0)

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 typedef typename FE::Traits::Basis::Traits::DomainLocal DomainLocal;
34 typedef typename FE::Traits::Basis::Traits::Range Range;
35
36 const FE& fe;
37
38public:
39 typedef typename FE::Traits::Basis::Traits::RangeField CT;
40
41 std::vector<CT> coeff;
42
43 FEFunction(const FE& fe_) : fe(fe_) { resetCoefficients(); }
44
45 void resetCoefficients() {
46 coeff.resize(fe.basis().size());
47 for(std::size_t i=0; i<coeff.size(); ++i)
48 coeff[i] = 0;
49 }
50
51 void setRandom(double max) {
52 coeff.resize(fe.basis().size());
53 for(std::size_t i=0; i<coeff.size(); ++i)
54 coeff[i] = ((1.0*std::rand()) / RAND_MAX - 0.5)*2.0*max;
55 }
56
57 void evaluate (const DomainLocal& x, Range& y) const {
58 std::vector<Range> yy;
59 fe.basis().evaluateFunction(x, yy);
60
61 y = 0.0;
62 for (std::size_t i=0; i<yy.size(); ++i)
63 y.axpy(coeff[i], yy[i]);
64 }
65};
66
67
68// Check if interpolation is consistent with basis evaluation.
80template<class FE>
81bool testInterpolation(const FE& fe, double eps, int n=5)
82{
83 bool success = true;
84 FEFunction<FE> f(fe);
85
86 std::vector<typename FEFunction<FE>::CT> coeff;
87 for(int i=0; i<n && success; ++i) {
88 // Set random coefficient vector
89 f.setRandom(100);
90
91 // Compute interpolation weights
92 fe.interpolation().interpolate(f, coeff);
93
94 // Check size of weight vector
95 if (coeff.size() != fe.basis().size()) {
96 std::cout << "Bug in LocalInterpolation for finite element type "
97 << Dune::className<FE>() << ":" << std::endl;
98 std::cout << " Interpolation vector has size " << coeff.size()
99 << std::endl;
100 std::cout << " Basis has size " << fe.basis().size() << std::endl;
101 std::cout << std::endl;
102 success = false;
103
104 // skip rest of loop since that depends on matching sizes
105 continue;
106 }
107
108 // Check if interpolation weights are equal to coefficients
109 for(std::size_t j=0; j<coeff.size() && success; ++j) {
110 if ( std::abs(coeff[j]-f.coeff[j]) >
111 eps*(std::max(std::abs(f.coeff[j]), 1.0)) )
112 {
113 std::cout << std::setprecision(16);
114 std::cout << "Bug in LocalInterpolation for finite element type "
115 << Dune::className<FE>() << ":" << std::endl;
116 std::cout << " Interpolation weight " << j << " differs by "
117 << std::abs(coeff[j]-f.coeff[j]) << " from coefficient of "
118 << "linear combination." << std::endl;
119 std::cout << std::endl;
120 success = false;
121 }
122 }
123 }
124 return success;
125}
126
127// check whether Jacobian agrees with FD approximation
139template<class Geo, class FE>
140bool testJacobian(const Geo &geo, const FE& fe, double eps, double delta,
141 std::size_t order = 2)
142{
143 typedef typename FE::Traits::Basis Basis;
144
145 typedef typename Basis::Traits::DomainField DF;
146 static const std::size_t dimDLocal = Basis::Traits::dimDomainLocal;
147 typedef typename Basis::Traits::DomainLocal DomainLocal;
148 static const std::size_t dimDGlobal = Basis::Traits::dimDomainGlobal;
149
150 static const std::size_t dimR = Basis::Traits::dimRange;
151 typedef typename Basis::Traits::Range Range;
152
153 typedef typename Basis::Traits::Jacobian Jacobian;
154
155 bool success = true;
156
157 // ////////////////////////////////////////////////////////////
158 // Check the partial derivatives by comparing them
159 // to finite difference approximations
160 // ////////////////////////////////////////////////////////////
161
162 // A set of test points
165
166 // Loop over all quadrature points
167 for (std::size_t i=0; i < quad.size(); i++) {
168
169 // Get a test point
170 const DomainLocal& testPoint = quad[i].position();
171
172 // Get the shape function derivatives there
173 std::vector<Jacobian> jacobians;
174 fe.basis().evaluateJacobian(testPoint, jacobians);
175 if(jacobians.size() != fe.basis().size()) {
176 std::cout << "Bug in evaluateJacobianGlobal() for finite element type "
177 << Dune::className<FE>() << ":" << std::endl;
178 std::cout << " Jacobian vector has size " << jacobians.size()
179 << std::endl;
180 std::cout << " Basis has size " << fe.basis().size() << std::endl;
181 std::cout << std::endl;
182 return false;
183 }
184
186 geo.jacobianTransposed(testPoint);
187
188 // Loop over all shape functions in this set
189 for (std::size_t j=0; j<fe.basis().size(); ++j) {
190
191 // basis.evaluateJacobian returns global derivatives, however we can
192 // only do local derivatives, so transform the derivatives back into
193 // local coordinates
195 for(std::size_t k = 0; k < dimR; ++k)
196 for(std::size_t l = 0; l < dimDGlobal; ++l)
197 for(std::size_t m = 0; m < dimDLocal; ++m)
198 localJacobian[k][m] += jacobians[j][k][l] * geoJT[m][l];
199
200 // Loop over all local directions
201 for (std::size_t m = 0; m < dimDLocal; ++m) {
202
203 // Compute an approximation to the derivative by finite differences
204 DomainLocal upPos = testPoint;
205 DomainLocal downPos = testPoint;
206
207 upPos[m] += delta;
208 downPos[m] -= delta;
209
210 std::vector<Range> upValues, downValues;
211
212 fe.basis().evaluateFunction(upPos, upValues);
213 fe.basis().evaluateFunction(downPos, downValues);
214
215 //Loop over all components
216 for(std::size_t k = 0; k < dimR; ++k) {
217
218 // The current partial derivative, just for ease of notation
219 double derivative = localJacobian[k][m];
220
221 double finiteDiff = (upValues[j][k] - downValues[j][k]) / (2*delta);
222
223 // Check
224 if ( std::abs(derivative-finiteDiff) >
225 eps/delta*(std::max(std::abs(finiteDiff), 1.0)) )
226 {
227 std::cout << std::setprecision(16);
228 std::cout << "Bug in evaluateJacobian() for finite element type "
229 << Dune::className<FE>() << ":" << std::endl;
230 std::cout << " Shape function derivative does not agree with "
231 << "FD approximation" << std::endl;
232 std::cout << " Shape function " << j << " component " << k
233 << " at position " << testPoint << ": derivative in "
234 << "local direction " << m << " is "
235 << derivative << ", but " << finiteDiff << " is "
236 << "expected." << std::endl;
237 std::cout << std::endl;
238 success = false;
239 }
240 } //Loop over all components
241 } // Loop over all local directions
242 } // Loop over all shape functions in this set
243 } // Loop over all quadrature points
244
245 return success;
246}
247
248// call tests for given finite element
249template<class Geo, class FE>
250bool testFE(const Geo &geo, const FE& fe, double eps, double delta,
251 unsigned order = 2)
252{
253 bool success = true;
254
255 success = testInterpolation(fe, eps) and success;
256 success = testJacobian(geo, fe, eps, delta, order) and success;
257
258 return success;
259}
260
261#endif // DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
A dense n x m matrix.
Definition: fmatrix.hh:68
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:97
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:225
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 ...
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)