DUNE PDELab (2.8)

test-localfe.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_LOCALFUNCTIONS_TEST_TEST_LOCALFE_HH
4#define DUNE_LOCALFUNCTIONS_TEST_TEST_LOCALFE_HH
5
14#include <iomanip>
15#include <iostream>
16#include <typeinfo>
17
20#include <dune/geometry/referenceelements.hh>
21
22#include <dune/localfunctions/common/virtualinterface.hh>
23#include <dune/localfunctions/common/virtualwrappers.hh>
24
25double TOL = 1e-9;
26// The FD approximation used for checking the Jacobian uses half of the
27// precision -- so we have to be a little bit more tolerant here.
28double jacobianTOL = 1e-5; // sqrt(TOL)
29
30// This class wraps one shape function of a local finite element as a function
31// that can be feed to the LocalInterpolation::interpolate method.
32template<class FE>
33class ShapeFunctionAsFunction
34{
35public:
36 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
37 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
38
39 struct Traits {
40 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
41 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
42 };
43
44 typedef typename FE::Traits::LocalBasisType::Traits::RangeFieldType CT;
45
46 ShapeFunctionAsFunction(const FE& fe, int shapeFunction) :
47 fe_(fe),
48 shapeFunction_(shapeFunction)
49 {}
50
51 void evaluate (const DomainType& x, RangeType& y) const
52 {
53 std::vector<RangeType> yy;
54 fe_.localBasis().evaluateFunction(x, yy);
55 y = yy[shapeFunction_];
56 }
57
58private:
59 const FE& fe_;
60 int shapeFunction_;
61};
62
63// This class wraps one shape function of a local finite element as a callable
64// that can be fed to the LocalInterpolation::interpolate method.
65template<class FE>
66class ShapeFunctionAsCallable
67{
68 // These types are deliberately private: They are not part of a Callable API
69 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
70 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
71public:
72
73 ShapeFunctionAsCallable(const FE& fe, int shapeFunction) :
74 fe_(fe),
75 shapeFunction_(shapeFunction)
76 {}
77
78 RangeType operator() (DomainType x) const
79 {
80 std::vector<RangeType> yy;
81 fe_.localBasis().evaluateFunction(x, yy);
82 return yy[shapeFunction_];
83 }
84
85private:
86 const FE& fe_;
87 int shapeFunction_;
88};
89
90
91// Check whether the degrees of freedom computed by LocalInterpolation
92// are dual to the shape functions. See Ciarlet, "The Finite Element Method
93// for Elliptic Problems", 1978, for details.
94template<class FE>
95bool testLocalInterpolation(const FE& fe)
96{
97 std::vector<typename ShapeFunctionAsFunction<FE>::CT> coeff;
98 for(size_t i=0; i<fe.size(); ++i)
99 {
101 // Part A: Feed the shape functions to the 'interpolate' method in form of
102 // a class providing an evaluate() method.
103 // This way is deprecated since dune-localfunctions 2.7.
105
106 // The i-th shape function as a function that 'interpolate' can deal with
107 ShapeFunctionAsFunction<FE> f(fe, i);
108
109 // Compute degrees of freedom for that shape function
110 // We expect the result to be the i-th unit vector
111 fe.localInterpolation().interpolate(f, coeff);
112
113 // Check size of weight vector
114 if (coeff.size() != fe.localBasis().size())
115 {
116 std::cout << "Bug in LocalInterpolation for finite element type "
117 << Dune::className(fe) << std::endl;
118 std::cout << " Interpolation produces " << coeff.size() << " degrees of freedom" << std::endl;
119 std::cout << " Basis has size " << fe.localBasis().size() << std::endl;
120 std::cout << std::endl;
121 return false;
122 }
123
124 // Check if interpolation weights are equal to coefficients
125 for(std::size_t j=0; j<coeff.size(); ++j)
126 {
127 if ( std::abs(coeff[j] - (i==j)) > TOL)
128 {
129 std::cout << std::setprecision(16);
130 std::cout << "Bug in LocalInterpolation for finite element type "
131 << Dune::className(fe) << std::endl;
132 std::cout << " Degree of freedom " << j << " applied to shape function " << i
133 << " yields value " << coeff[j] << ", not the expected value " << (i==j) << std::endl;
134 std::cout << std::endl;
135 return false;
136 }
137 }
138
139
141 // Part B: Redo the same test, but feed the shape functions to the
142 // 'interpolate' method in form of a callable.
144
145 // The i-th shape function as a function that 'interpolate' can deal with
146 ShapeFunctionAsCallable<FE> sfAsCallable(fe, i);
147
148 // Compute degrees of freedom for that shape function
149 // We expect the result to be the i-th unit vector
150 fe.localInterpolation().interpolate(sfAsCallable, coeff);
151
152 // Check size of weight vector
153 if (coeff.size() != fe.localBasis().size())
154 {
155 std::cout << "Bug in LocalInterpolation for finite element type "
156 << Dune::className(fe) << std::endl;
157 std::cout << " Interpolation produces " << coeff.size() << " degrees of freedom" << std::endl;
158 std::cout << " Basis has size " << fe.localBasis().size() << std::endl;
159 std::cout << std::endl;
160 return false;
161 }
162
163 // Check if interpolation weights are equal to coefficients
164 for(std::size_t j=0; j<coeff.size(); ++j)
165 {
166 if ( std::abs(coeff[j] - (i==j)) > TOL)
167 {
168 std::cout << std::setprecision(16);
169 std::cout << "Bug in LocalInterpolation for finite element type "
170 << Dune::className(fe) << std::endl;
171 std::cout << " Degree of freedom " << j << " applied to shape function " << i
172 << " yields value " << coeff[j] << ", not the expected value " << (i==j) << std::endl;
173 std::cout << std::endl;
174 return false;
175 }
176 }
177 }
178 return true;
179}
180
181
182// Check whether the space spanned by the shape functions
183// contains the constant functions
184template<class FE>
185bool testCanRepresentConstants(const FE& fe,
186 unsigned order = 5)
187{
188 typedef typename FE::Traits::LocalBasisType LB;
189 using RangeType = typename LB::Traits::RangeType;
190
191 bool success = true;
192
193 // Construct the constant '1' function
194 auto constantOne = [](const typename LB::Traits::DomainType& xi) { return RangeType(1.0); };
195
196 // Project the constant function onto the FE space
197 std::vector<double> coefficients;
198 fe.localInterpolation().interpolate(constantOne, coefficients);
199
200 // A set of test points
201 const auto& quad = Dune::QuadratureRules<double,LB::Traits::dimDomain>::rule(fe.type(),order);
202
203 // Loop over all quadrature points
204 for (size_t i=0; i<quad.size(); i++) {
205
206 // Get a test point
207 const auto& testPoint = quad[i].position();
208
209 // Compute value of the representation of constantOne at the test point
210 std::vector<RangeType> values;
211 fe.localBasis().evaluateFunction(testPoint, values);
212
213 RangeType sum(0);
214 for (size_t j=0; j<values.size(); j++)
215 sum += coefficients[j] * values[j];
216
217 if ((RangeType(1.0)-sum).two_norm() > TOL)
218 {
219 std::cout << "Finite element type " << Dune::className(fe)
220 << " cannot represent constant functions!" << std::endl;
221 std::cout << " At position: " << testPoint << ","
222 << std::endl;
223 std::cout << " discrete approximation of the '1' function has value " << sum
224 << std::endl;
225 std::cout << std::endl;
226 success = false;
227 }
228
229 } // Loop over all quadrature points
230
231 return success;
232}
233
234// check whether Jacobian agrees with FD approximation
235template<class FE>
236bool testJacobian(const FE& fe,
237 unsigned order = 2,
238 const std::function<bool(const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip = nullptr)
239{
240 typedef typename FE::Traits::LocalBasisType LB;
241
242 bool success = true;
243
244 // ////////////////////////////////////////////////////////////
245 // Check the partial derivatives by comparing them
246 // to finite difference approximations
247 // ////////////////////////////////////////////////////////////
248
249 // A set of test points
252
253 // Loop over all quadrature points
254 for (size_t i=0; i<quad.size(); i++) {
255
256 // Get a test point
258 quad[i].position();
259
260 // Get the shape function derivatives there
261 std::vector<typename LB::Traits::JacobianType> jacobians;
262 fe.localBasis().evaluateJacobian(testPoint, jacobians);
263 if(jacobians.size() != fe.localBasis().size()) {
264 std::cout << "Bug in evaluateJacobian() for finite element type "
265 << Dune::className(fe) << std::endl;
266 std::cout << " Jacobian vector has size " << jacobians.size()
267 << std::endl;
268 std::cout << " Basis has size " << fe.localBasis().size()
269 << std::endl;
270 std::cout << std::endl;
271 return false;
272 }
273
274 // Skip this test point if we are supposed to
275 if (derivativePointSkip && derivativePointSkip(quad[i].position()))
276 continue;
277
278 // Loop over all directions
279 for (int k=0; k<LB::Traits::dimDomain; k++) {
280
281 // Compute an approximation to the derivative by finite differences
284
285 upPos[k] += jacobianTOL;
286 downPos[k] -= jacobianTOL;
287
288 std::vector<typename LB::Traits::RangeType> upValues, downValues;
289
290 fe.localBasis().evaluateFunction(upPos, upValues);
291 fe.localBasis().evaluateFunction(downPos, downValues);
292
293 // Loop over all shape functions in this set
294 for (unsigned int j=0; j<fe.localBasis().size(); ++j) {
295 //Loop over all components
296 for(int l=0; l < LB::Traits::dimRange; ++l) {
297
298 // The current partial derivative, just for ease of notation
299 double derivative = jacobians[j][l][k];
300
301 double finiteDiff = (upValues[j][l] - downValues[j][l])
302 / (2*jacobianTOL);
303
304 // Check
305 if ( std::abs(derivative-finiteDiff) >
306 TOL/jacobianTOL*((std::abs(finiteDiff)>1) ? std::abs(finiteDiff) : 1.) )
307 {
308 std::cout << std::setprecision(16);
309 std::cout << "Bug in evaluateJacobian() for finite element type "
310 << Dune::className(fe) << std::endl;
311 std::cout << " Shape function derivative does not agree with "
312 << "FD approximation" << std::endl;
313 std::cout << " Shape function " << j << " component " << l
314 << " at position " << testPoint << ": derivative in "
315 << "direction " << k << " is " << derivative << ", but "
316 << finiteDiff << " is expected." << std::endl;
317 std::cout << std::endl;
318 success = false;
319 }
320 } //Loop over all components
321 } // Loop over all shape functions in this set
322 } // Loop over all directions
323 } // Loop over all quadrature points
324
325 return success;
326}
327
333{
334 template <class FE>
335 static bool test(const FE& fe,
336 double eps, double delta, unsigned int diffOrder, std::size_t order = 2,
337 const std::function<bool(const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip = nullptr)
338 {
339 bool success = true;
340
341 if (diffOrder > 2)
342 std::cout << "No test for differentiability orders larger than 2!" << std::endl;
343
344 if (diffOrder >= 2)
345 success = success and testOrder2(fe, eps, delta, order, derivativePointSkip);
346
347 if (diffOrder >= 1)
348 success = success and testOrder1(fe, eps, delta, order, derivativePointSkip);
349
350 success = success and testOrder0(fe, eps, delta, order);
351
352 return success;
353 }
354
356 template <class FE>
357 static bool testOrder0(const FE& fe,
358 double eps,
359 double delta,
360 std::size_t order = 2)
361 {
362 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
363 constexpr auto dimDomain = FE::Traits::LocalBasisType::Traits::dimDomain;
364
365 bool success = true;
366
368 // Check the partial derivatives by comparing them
369 // to finite difference approximations
371
372 // A set of test points
373 const auto& quad = Dune::QuadratureRules<double, dimDomain>::rule(fe.type(),
374 order);
375
376 // Loop over all quadrature points
377 for (size_t i = 0; i < quad.size(); i++)
378 {
379 // Get a test point
380 const Dune::FieldVector<double, dimDomain>& testPoint = quad[i].position();
381
382 // Get the shape function values there using the 'partial' method
383 std::vector<RangeType> partialValues;
384 std::array<unsigned int, dimDomain> multiIndex;
385 std::fill(multiIndex.begin(), multiIndex.end(), 0);
386
387 fe.localBasis().partial(multiIndex, testPoint, partialValues);
388
389 if (partialValues.size() != fe.localBasis().size())
390 {
391 std::cout << "Bug in partial() for finite element type "
392 << Dune::className(fe) << std::endl;
393 std::cout << " values vector has size "
394 << partialValues.size() << std::endl;
395 std::cout << " Basis has size " << fe.localBasis().size()
396 << std::endl;
397 std::cout << std::endl;
398 return false;
399 }
400
401 // Get reference values
402 std::vector<RangeType> referenceValues;
403 fe.localBasis().evaluateFunction(testPoint, referenceValues);
404
405 // Loop over all shape functions in this set
406 for (unsigned int j = 0; j < fe.localBasis().size(); ++j)
407 {
408 // Loop over all components
409 for (int l = 0; l < FE::Traits::LocalBasisType::Traits::dimRange; ++l)
410 {
411 // Check the 'partial' method
412 if (std::abs(partialValues[j][l] - referenceValues[j][l])
413 > TOL / jacobianTOL
414 * ((std::abs(referenceValues[j][l]) > 1) ? std::abs(referenceValues[j][l]) : 1.))
415 {
416 std::cout << std::setprecision(16);
417 std::cout << "Bug in partial() for finite element type "
418 << Dune::className(fe) << std::endl;
419 std::cout << " Shape function value does not agree with "
420 << "output of method evaluateFunction." << std::endl;
421 std::cout << " Shape function " << j << " component " << l
422 << " at position " << testPoint << ": value is " << partialValues[j][l]
423 << ", but " << referenceValues[j][l] << " is expected." << std::endl;
424 std::cout << std::endl;
425 success = false;
426 }
427
428 } //Loop over all components
429 } // Loop over all shape functions in this set
430 } // Loop over all quadrature points
431
432 return success;
433 }
434
436 template <class FE>
437 static bool testOrder1(const FE& fe,
438 double eps,
439 double delta,
440 std::size_t order = 2,
441 const std::function<bool(const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip = nullptr)
442 {
443 typedef typename FE::Traits::LocalBasisType LB;
444 typedef typename LB::Traits::RangeFieldType RangeField;
445
446 bool success = true;
447
449 // Check the partial derivatives by comparing them
450 // to finite difference approximations
452
453 // A set of test points
456 order);
457
458 // Loop over all quadrature points
459 for (size_t i = 0; i < quad.size(); i++)
460 {
461 // Get a test point
462 const Dune::FieldVector<double, LB::Traits::dimDomain>& testPoint = quad[i].position();
463
464 // Skip the test points we are supposed to skip
465 if (derivativePointSkip && derivativePointSkip(testPoint))
466 continue;
467
468 // Loop over all directions
469 for (int k = 0; k < LB::Traits::dimDomain; k++)
470 {
471 // Get the shape function derivatives there using the 'partial' method
472 std::vector<typename LB::Traits::RangeType> firstPartialDerivatives;
473 std::array<unsigned int, LB::Traits::dimDomain> multiIndex;
474 std::fill(multiIndex.begin(), multiIndex.end(), 0);
475 multiIndex[k]++;
476 fe.localBasis().partial(multiIndex, testPoint, firstPartialDerivatives);
477 if (firstPartialDerivatives.size() != fe.localBasis().size())
478 {
479 std::cout << "Bug in partial() for finite element type "
480 << Dune::className(fe) << std::endl;
481 std::cout << " firstPartialDerivatives vector has size "
482 << firstPartialDerivatives.size() << std::endl;
483 std::cout << " Basis has size " << fe.localBasis().size()
484 << std::endl;
485 std::cout << std::endl;
486 return false;
487 }
488
489 // Compute an approximation to the derivative by finite differences
492
493 upPos[k] += jacobianTOL;
494 downPos[k] -= jacobianTOL;
495
496 std::vector<typename LB::Traits::RangeType> upValues, downValues;
497
498 fe.localBasis().evaluateFunction(upPos, upValues);
499 fe.localBasis().evaluateFunction(downPos, downValues);
500
501 // Loop over all shape functions in this set
502 for (unsigned int j = 0; j < fe.localBasis().size(); ++j)
503 {
504 // Loop over all components
505 for (int l = 0; l < LB::Traits::dimRange; ++l)
506 {
507 RangeField finiteDiff = (upValues[j][l] - downValues[j][l])
508 / (2 * jacobianTOL);
509
510 // Check the 'partial' method
511 RangeField partialDerivative = firstPartialDerivatives[j][l];
512 if (std::abs(partialDerivative - finiteDiff)
513 > TOL / jacobianTOL
514 * ((std::abs(finiteDiff) > 1) ? std::abs(finiteDiff) : 1.))
515 {
516 std::cout << std::setprecision(16);
517 std::cout << "Bug in partial() for finite element type "
518 << Dune::className(fe) << std::endl;
519 std::cout << " Shape function derivative does not agree with "
520 << "FD approximation" << std::endl;
521 std::cout << " Shape function " << j << " component " << l
522 << " at position " << testPoint << ": derivative in "
523 << "direction " << k << " is " << partialDerivative << ", but "
524 << finiteDiff << " is expected." << std::endl;
525 std::cout << std::endl;
526 success = false;
527 }
528
529 } // Loop over all directions
530 } // Loop over all shape functions in this set
531 } //Loop over all components
532 } // Loop over all quadrature points
533
534 return success;
535 }
536
538 template <class FE>
539 static bool testOrder2(const FE& fe,
540 double eps,
541 double delta,
542 std::size_t order = 2,
543 const std::function<bool(const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip = nullptr)
544 {
545 typedef typename FE::Traits::LocalBasisType LocalBasis;
546 typedef typename LocalBasis::Traits::DomainFieldType DF;
547 typedef typename LocalBasis::Traits::DomainType Domain;
548 static const int dimDomain = LocalBasis::Traits::dimDomain;
549
550 static const std::size_t dimR = LocalBasis::Traits::dimRange;
551 typedef typename LocalBasis::Traits::RangeType Range;
552 typedef typename LocalBasis::Traits::RangeFieldType RangeField;
553
554 bool success = true;
555
557 // Check the partial derivatives by comparing them
558 // to finite difference approximations
560
561 // A set of test points
564
565 // Loop over all quadrature points
566 for (std::size_t i = 0; i < quad.size(); i++)
567 {
568 // Get a test point
569 const Domain& testPoint = quad[i].position();
570
571 // Skip the test points we are supposed to skip
572 if (derivativePointSkip && derivativePointSkip(testPoint))
573 continue;
574
575 // For testing the 'partial' method
576 std::array<std::vector<Dune::FieldMatrix<RangeField, dimDomain, dimDomain> >, dimR> partialHessians;
577 for (size_t k = 0; k < dimR; k++)
578 partialHessians[k].resize(fe.size());
579
580 //loop over all local directions
581 for (int dir0 = 0; dir0 < dimDomain; dir0++)
582 {
583 for (int dir1 = 0; dir1 < dimDomain; dir1++)
584 {
585 // Get the shape function derivatives there using the 'partial' method
586 std::vector<Range> secondPartialDerivative;
587 std::array<unsigned int,dimDomain> multiIndex;
588 std::fill(multiIndex.begin(), multiIndex.end(), 0);
589 multiIndex[dir0]++;
590 multiIndex[dir1]++;
591 fe.localBasis().partial(multiIndex, testPoint, secondPartialDerivative);
592
593 if (secondPartialDerivative.size() != fe.localBasis().size())
594 {
595 std::cout << "Bug in partial() for finite element type "
596 << Dune::className<FE>() << ":" << std::endl;
597 std::cout << " return vector has size " << secondPartialDerivative.size()
598 << std::endl;
599 std::cout << " Basis has size " << fe.localBasis().size()
600 << std::endl;
601 std::cout << std::endl;
602 return false;
603 }
604
605 //combine to Hesse matrices
606 for (size_t k = 0; k < dimR; k++)
607 for (std::size_t j = 0; j < fe.localBasis().size(); ++j)
608 partialHessians[k][j][dir0][dir1] = secondPartialDerivative[j][k];
609
610 }
611 } //loop over all directions
612
613 // Loop over all local directions
614 for (std::size_t dir0 = 0; dir0 < dimDomain; ++dir0)
615 {
616 for (unsigned int dir1 = 0; dir1 < dimDomain; dir1++)
617 {
618 // Compute an approximation to the derivative by finite differences
619 std::array<Domain,4> neighbourPos;
620 std::fill(neighbourPos.begin(), neighbourPos.end(), testPoint);
621
622 neighbourPos[0][dir0] += delta;
623 neighbourPos[0][dir1] += delta;
624 neighbourPos[1][dir0] -= delta;
625 neighbourPos[1][dir1] += delta;
626 neighbourPos[2][dir0] += delta;
627 neighbourPos[2][dir1] -= delta;
628 neighbourPos[3][dir0] -= delta;
629 neighbourPos[3][dir1] -= delta;
630
631 std::array<std::vector<Range>, 4> neighbourValues;
632 for (int k = 0; k < 4; k++)
633 fe.localBasis().evaluateFunction(neighbourPos[k],
634 neighbourValues[k]);
635
636 // Loop over all shape functions in this set
637 for (std::size_t j = 0; j < fe.localBasis().size(); ++j)
638 {
639 //Loop over all components
640 for (std::size_t k = 0; k < dimR; ++k)
641 {
642 RangeField finiteDiff = (neighbourValues[0][j][k]
643 - neighbourValues[1][j][k] - neighbourValues[2][j][k]
644 + neighbourValues[3][j][k]) / (4 * delta * delta);
645
646 // The current partial derivative, just for ease of notation, evaluated by the 'partial' method
647 RangeField partialDerivative = partialHessians[k][j][dir0][dir1];
648
649 // Check
650 if (std::abs(partialDerivative - finiteDiff)
651 > eps / delta * (std::max(std::abs(finiteDiff), 1.0)))
652 {
653 std::cout << std::setprecision(16);
654 std::cout << "Bug in partial() for finite element type "
655 << Dune::className<FE>() << ":" << std::endl;
656 std::cout << " Second shape function derivative does not agree with "
657 << "FD approximation" << std::endl;
658 std::cout << " Shape function " << j << " component " << k
659 << " at position " << testPoint << ": derivative in "
660 << "local direction (" << dir0 << ", " << dir1 << ") is "
661 << partialDerivative << ", but " << finiteDiff
662 << " is expected." << std::endl;
663 std::cout << std::endl;
664 success = false;
665 }
666
667 } //Loop over all components
668 }
669 } // Loop over all local directions
670 } // Loop over all shape functions in this set
671 } // Loop over all quadrature points
672
673 return success;
674 }
675
676};
677
678// Flags for disabling parts of testFE
679enum {
680 DisableNone = 0,
681 DisableLocalInterpolation = 1,
682 DisableVirtualInterface = 2,
683 DisableJacobian = 4,
684 DisableEvaluate = 8,
685 DisableRepresentConstants = 16
686};
687
698template<class FE>
699bool testFE(const FE& fe,
700 char disabledTests = DisableNone,
701 unsigned int diffOrder = 0,
702 const std::function<bool(const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip = nullptr)
703{
704 // Order of the quadrature rule used to generate test points
705 unsigned int quadOrder = 2;
706
707 bool success = true;
708
709 if (FE::Traits::LocalBasisType::Traits::dimDomain != fe.type().dim())
710 {
711 std::cout << "Bug in type() for finite element type "
712 << Dune::className(fe) << std::endl;
713 std::cout << " Coordinate dimension is " << FE::Traits::LocalBasisType::Traits::dimDomain << std::endl;
714 std::cout << " but GeometryType is " << fe.type() << " with dimension " << fe.type().dim() << std::endl;
715 success = false;
716 }
717
718 if (fe.size() != fe.localBasis().size())
719 {
720 std::cout << "Bug in finite element type "
721 << Dune::className(fe) << std::endl;
722 std::cout << " Size reported by LocalFiniteElement is " << fe.size() << std::endl;
723 std::cout << " but size reported by LocalBasis is " << fe.localBasis().size() << std::endl;
724 success = false;
725 }
726
727 // Make sure evaluateFunction returns the correct number of values
728 std::vector<typename FE::Traits::LocalBasisType::Traits::RangeType> values;
729 fe.localBasis().evaluateFunction(Dune::ReferenceElements<double,FE::Traits::LocalBasisType::Traits::dimDomain>::general(fe.type()).position(0,0), values);
730
731 if (values.size() != fe.size())
732 {
733 std::cout << "Bug in finite element type "
734 << Dune::className(fe) << std::endl;
735 std::cout << " LocalFiniteElement.size() returns " << fe.size() << "," << std::endl;
736 std::cout << " but LocalBasis::evaluateFunction returns " << values.size() << " values!" << std::endl;
737 success = false;
738 }
739
740 if (fe.size() != fe.localCoefficients().size())
741 {
742 std::cout << "Bug in finite element type "
743 << Dune::className(fe) << std::endl;
744 std::cout << " Size reported by LocalFiniteElement is " << fe.size() << std::endl;
745 std::cout << " but size reported by LocalCoefficients is " << fe.localCoefficients().size() << std::endl;
746 success = false;
747 }
748
749 const auto& lc = fe.localCoefficients();
750 for(size_t i=0; i<lc.size(); i++)
751 {
752 const auto& lk = lc.localKey(i);
753 if (lk.codim() > fe.type().dim())
754 {
755 std::cout << "Bug in finite element type "
756 << Dune::className(fe) << std::endl;
757 std::cout << " Codimension reported by localKey(" << i << ") is " << lk.codim() << std::endl;
758 std::cout << " but geometry is " << fe.type() << " with dimension " << fe.type().dim() << std::endl;
759 success = false;
760 }
761 }
762
763 if (not (disabledTests & DisableLocalInterpolation))
764 {
765 success = testLocalInterpolation<FE>(fe) and success;
766 }
767
768 if (not (disabledTests & DisableRepresentConstants))
769 {
770 success = testCanRepresentConstants<FE>(fe) and success;
771 }
772
773 if (not (disabledTests & DisableJacobian))
774 {
775 success = testJacobian<FE>(fe, quadOrder, derivativePointSkip) and success;
776 }
777 else
778 {
779 // make sure diffOrder is 0
780 success = (diffOrder == 0) and success;
781 }
782
783 if (not (disabledTests & DisableEvaluate))
784 {
785 success = TestPartial::test(fe, TOL, jacobianTOL, diffOrder, quadOrder, derivativePointSkip) and success;
786 }
787
788 if (not (disabledTests & DisableVirtualInterface))
789 {
790 typedef typename FE::Traits::LocalBasisType::Traits ImplementationLBTraits;
792 typedef typename Dune::LocalFiniteElementVirtualImp<FE> VirtualFEImp;
793
794 const VirtualFEImp virtualFE(fe);
795 if (not (disabledTests & DisableLocalInterpolation))
796 success = testLocalInterpolation<VirtualFEInterface>(virtualFE) and success;
797 if (not (disabledTests & DisableJacobian))
798 {
799 success = testJacobian<VirtualFEInterface>(virtualFE, quadOrder, derivativePointSkip) and success;
800 }
801 else
802 {
803 // make sure diffOrder is 0
804 success = (diffOrder == 0) and success;
805 }
806 }
807
808 return success;
809}
810
811#define TEST_FE(A) { bool b = testFE(A); std::cout << "testFE(" #A ") " << (b?"succeeded\n":"failed\n"); success &= b; }
812#define TEST_FE2(A,B) { bool b = testFE(A, B); std::cout << "testFE(" #A ", " #B ") " << (b?"succeeded\n":"failed\n"); success &= b; }
813#define TEST_FE3(A,B,C) { bool b = testFE(A, B, C); std::cout << "testFE(" #A ", " #B ", " #C ") " << (b?"succeeded\n":"failed\n"); success &= b; }
814#define TEST_FE4(A,B,C,D) { bool b = testFE(A, B, C, D); std::cout << "testFE(" #A ", " #B ", " #C ", " #D ") " << (b?"succeeded\n":"failed\n"); success &= b; }
815
816#endif // DUNE_LOCALFUNCTIONS_TEST_TEST_LOCALFE_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:95
class for wrapping a finite element using the virtual interface
Definition: virtualwrappers.hh:238
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:284
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.
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
std::string className()
Provide the demangled class name of a type T as a string.
Definition: classname.hh:45
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:168
Helper class to test the 'partial' method.
Definition: test-localfe.hh:333
static bool testOrder0(const FE &fe, double eps, double delta, std::size_t order=2)
Test the 'partial' method for zero-order partial derivatives, i.e., values.
Definition: test-localfe.hh:357
static bool testOrder1(const FE &fe, double eps, double delta, std::size_t order=2, const std::function< bool(const typename FE::Traits::LocalBasisType::Traits::DomainType &)> derivativePointSkip=nullptr)
Test the 'partial' method for first-order partial derivatives.
Definition: test-localfe.hh:437
static bool testOrder2(const FE &fe, double eps, double delta, std::size_t order=2, const std::function< bool(const typename FE::Traits::LocalBasisType::Traits::DomainType &)> derivativePointSkip=nullptr)
Test second-order partial derivatives.
Definition: test-localfe.hh:539
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)