Dune Core Modules (2.9.0)

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