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