5#ifndef DUNE_LOCALFUNCTIONS_TEST_TEST_LOCALFE_HH
6#define DUNE_LOCALFUNCTIONS_TEST_TEST_LOCALFE_HH
22#include <dune/geometry/referenceelements.hh>
24#include <dune/localfunctions/common/virtualinterface.hh>
25#include <dune/localfunctions/common/virtualwrappers.hh>
30double jacobianTOL = 1e-5;
35class ShapeFunctionAsFunction
38 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
39 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
42 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
43 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
46 typedef typename FE::Traits::LocalBasisType::Traits::RangeFieldType CT;
48 ShapeFunctionAsFunction(
const FE& fe,
int shapeFunction) :
50 shapeFunction_(shapeFunction)
53 void evaluate (
const DomainType& x, RangeType& y)
const
55 std::vector<RangeType> yy;
56 fe_.localBasis().evaluateFunction(x, yy);
57 y = yy[shapeFunction_];
68class ShapeFunctionAsCallable
71 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
72 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
75 ShapeFunctionAsCallable(
const FE& fe,
int shapeFunction) :
77 shapeFunction_(shapeFunction)
80 RangeType operator() (DomainType x)
const
82 std::vector<RangeType> yy;
83 fe_.localBasis().evaluateFunction(x, yy);
84 return yy[shapeFunction_];
97bool testLocalInterpolation(
const FE& fe)
99 std::vector<typename ShapeFunctionAsFunction<FE>::CT> coeff;
100 for(
size_t i=0; i<fe.size(); ++i)
109 ShapeFunctionAsFunction<FE> f(fe, i);
113 fe.localInterpolation().interpolate(f, coeff);
116 if (coeff.size() != fe.localBasis().size())
118 std::cout <<
"Bug in LocalInterpolation for finite element type "
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;
127 for(std::size_t j=0; j<coeff.size(); ++j)
129 if ( std::abs(coeff[j] - (i==j)) > TOL)
131 std::cout << std::setprecision(16);
132 std::cout <<
"Bug in LocalInterpolation for finite element type "
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;
148 ShapeFunctionAsCallable<FE> sfAsCallable(fe, i);
152 fe.localInterpolation().interpolate(sfAsCallable, coeff);
155 if (coeff.size() != fe.localBasis().size())
157 std::cout <<
"Bug in LocalInterpolation for finite element type "
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;
166 for(std::size_t j=0; j<coeff.size(); ++j)
168 if ( std::abs(coeff[j] - (i==j)) > TOL)
170 std::cout << std::setprecision(16);
171 std::cout <<
"Bug in LocalInterpolation for finite element type "
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;
187bool testCanRepresentConstants(
const FE& fe,
190 typedef typename FE::Traits::LocalBasisType LB;
191 using RangeType =
typename LB::Traits::RangeType;
196 auto constantOne = [](
const typename LB::Traits::DomainType& xi) {
return RangeType(1.0); };
199 std::vector<double> coefficients;
200 fe.localInterpolation().interpolate(constantOne, coefficients);
206 for (
size_t i=0; i<quad.size(); i++) {
209 const auto& testPoint = quad[i].position();
212 std::vector<RangeType> values;
213 fe.localBasis().evaluateFunction(testPoint, values);
216 for (
size_t j=0; j<values.size(); j++)
217 sum += coefficients[j] * values[j];
219 if ((RangeType(1.0)-sum).two_norm() > TOL)
222 <<
" cannot represent constant functions!" << std::endl;
223 std::cout <<
" At position: " << testPoint <<
","
225 std::cout <<
" discrete approximation of the '1' function has value " << sum
227 std::cout << std::endl;
238bool testJacobian(
const FE& fe,
240 const std::function<
bool(
const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip =
nullptr)
242 typedef typename FE::Traits::LocalBasisType LB;
256 for (
size_t i=0; i<quad.size(); i++) {
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 "
268 std::cout <<
" Jacobian vector has size " << jacobians.size()
270 std::cout <<
" Basis has size " << fe.localBasis().size()
272 std::cout << std::endl;
277 if (derivativePointSkip && derivativePointSkip(quad[i].position()))
281 for (
int k=0; k<LB::Traits::dimDomain; k++) {
287 upPos[k] += jacobianTOL;
288 downPos[k] -= jacobianTOL;
290 std::vector<typename LB::Traits::RangeType> upValues, downValues;
292 fe.localBasis().evaluateFunction(upPos, upValues);
293 fe.localBasis().evaluateFunction(downPos, downValues);
296 for (
unsigned int j=0; j<fe.localBasis().size(); ++j) {
298 for(
int l=0; l < LB::Traits::dimRange; ++l) {
303 double finiteDiff = (upValues[j][l] - downValues[j][l])
308 TOL/jacobianTOL*((std::abs(finiteDiff)>1) ? std::abs(finiteDiff) : 1.) )
310 std::cout << std::setprecision(16);
311 std::cout <<
"Bug in evaluateJacobian() for finite element type "
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;
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)
344 std::cout <<
"No test for differentiability orders larger than 2!" << std::endl;
347 success = success and
testOrder2(fe, eps, delta, order, derivativePointSkip);
350 success = success and
testOrder1(fe, eps, delta, order, derivativePointSkip);
352 success = success and
testOrder0(fe, eps, delta, order);
362 std::size_t order = 2)
364 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
365 constexpr auto dimDomain = FE::Traits::LocalBasisType::Traits::dimDomain;
379 for (
size_t i = 0; i < quad.size(); i++)
385 std::vector<RangeType> partialValues;
386 std::array<unsigned int, dimDomain> multiIndex;
387 std::fill(multiIndex.begin(), multiIndex.end(), 0);
389 fe.localBasis().partial(multiIndex, testPoint, partialValues);
391 if (partialValues.size() != fe.localBasis().size())
393 std::cout <<
"Bug in partial() for finite element type "
395 std::cout <<
" values vector has size "
396 << partialValues.size() << std::endl;
397 std::cout <<
" Basis has size " << fe.localBasis().size()
399 std::cout << std::endl;
404 std::vector<RangeType> referenceValues;
405 fe.localBasis().evaluateFunction(testPoint, referenceValues);
408 for (
unsigned int j = 0; j < fe.localBasis().size(); ++j)
411 for (
int l = 0; l < FE::Traits::LocalBasisType::Traits::dimRange; ++l)
414 if (std::abs(partialValues[j][l] - referenceValues[j][l])
416 * ((std::abs(referenceValues[j][l]) > 1) ? std::abs(referenceValues[j][l]) : 1.))
418 std::cout << std::setprecision(16);
419 std::cout <<
"Bug in partial() for finite element type "
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;
442 std::size_t order = 2,
443 const std::function<
bool(
const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip =
nullptr)
445 typedef typename FE::Traits::LocalBasisType LB;
446 typedef typename LB::Traits::RangeFieldType RangeField;
461 for (
size_t i = 0; i < quad.size(); i++)
467 if (derivativePointSkip && derivativePointSkip(testPoint))
471 for (
int k = 0; k < LB::Traits::dimDomain; k++)
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);
478 fe.localBasis().partial(multiIndex, testPoint, firstPartialDerivatives);
479 if (firstPartialDerivatives.size() != fe.localBasis().size())
481 std::cout <<
"Bug in partial() for finite element type "
483 std::cout <<
" firstPartialDerivatives vector has size "
484 << firstPartialDerivatives.size() << std::endl;
485 std::cout <<
" Basis has size " << fe.localBasis().size()
487 std::cout << std::endl;
495 upPos[k] += jacobianTOL;
496 downPos[k] -= jacobianTOL;
498 std::vector<typename LB::Traits::RangeType> upValues, downValues;
500 fe.localBasis().evaluateFunction(upPos, upValues);
501 fe.localBasis().evaluateFunction(downPos, downValues);
504 for (
unsigned int j = 0; j < fe.localBasis().size(); ++j)
507 for (
int l = 0; l < LB::Traits::dimRange; ++l)
509 RangeField finiteDiff = (upValues[j][l] - downValues[j][l])
513 RangeField partialDerivative = firstPartialDerivatives[j][l];
514 if (std::abs(partialDerivative - finiteDiff)
516 * ((std::abs(finiteDiff) > 1) ? std::abs(finiteDiff) : 1.))
518 std::cout << std::setprecision(16);
519 std::cout <<
"Bug in partial() for finite element type "
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;
544 std::size_t order = 2,
545 const std::function<
bool(
const typename FE::Traits::LocalBasisType::Traits::DomainType&)> derivativePointSkip =
nullptr)
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;
552 static const std::size_t dimR = LocalBasis::Traits::dimRange;
553 typedef typename LocalBasis::Traits::RangeType Range;
554 typedef typename LocalBasis::Traits::RangeFieldType RangeField;
568 for (std::size_t i = 0; i < quad.size(); i++)
571 const Domain& testPoint = quad[i].position();
574 if (derivativePointSkip && derivativePointSkip(testPoint))
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());
583 for (
int dir0 = 0; dir0 < dimDomain; dir0++)
585 for (
int dir1 = 0; dir1 < dimDomain; dir1++)
588 std::vector<Range> secondPartialDerivative;
589 std::array<unsigned int,dimDomain> multiIndex;
590 std::fill(multiIndex.begin(), multiIndex.end(), 0);
593 fe.localBasis().partial(multiIndex, testPoint, secondPartialDerivative);
595 if (secondPartialDerivative.size() != fe.localBasis().size())
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()
601 std::cout <<
" Basis has size " << fe.localBasis().size()
603 std::cout << std::endl;
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];
616 for (std::size_t dir0 = 0; dir0 < dimDomain; ++dir0)
618 for (
unsigned int dir1 = 0; dir1 < dimDomain; dir1++)
621 std::array<Domain,4> neighbourPos;
622 std::fill(neighbourPos.begin(), neighbourPos.end(), testPoint);
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;
633 std::array<std::vector<Range>, 4> neighbourValues;
634 for (
int k = 0; k < 4; k++)
635 fe.localBasis().evaluateFunction(neighbourPos[k],
639 for (std::size_t j = 0; j < fe.localBasis().size(); ++j)
642 for (std::size_t k = 0; k < dimR; ++k)
644 RangeField finiteDiff = (neighbourValues[0][j][k]
645 - neighbourValues[1][j][k] - neighbourValues[2][j][k]
646 + neighbourValues[3][j][k]) / (4 * delta * delta);
649 RangeField partialDerivative = partialHessians[k][j][dir0][dir1];
652 if (std::abs(partialDerivative - finiteDiff)
653 > eps / delta * (std::max(std::abs(finiteDiff), 1.0)))
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;
683 DisableLocalInterpolation = 1,
684 DisableVirtualInterface = 2,
687 DisableRepresentConstants = 16
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)
707 unsigned int quadOrder = 2;
711 if (FE::Traits::LocalBasisType::Traits::dimDomain != fe.type().dim())
713 std::cout <<
"Bug in type() for finite element type "
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;
720 if (fe.size() != fe.localBasis().size())
722 std::cout <<
"Bug in finite element type "
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;
730 std::vector<typename FE::Traits::LocalBasisType::Traits::RangeType> values;
733 if (values.size() != fe.size())
735 std::cout <<
"Bug in finite element type "
737 std::cout <<
" LocalFiniteElement.size() returns " << fe.size() <<
"," << std::endl;
738 std::cout <<
" but LocalBasis::evaluateFunction returns " << values.size() <<
" values!" << std::endl;
742 if (fe.size() != fe.localCoefficients().size())
744 std::cout <<
"Bug in finite element type "
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;
751 const auto& lc = fe.localCoefficients();
752 for(
size_t i=0; i<lc.size(); i++)
754 const auto& lk = lc.localKey(i);
755 if (lk.codim() > fe.type().dim())
757 std::cout <<
"Bug in finite element type "
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;
765 if (not (disabledTests & DisableLocalInterpolation))
767 success = testLocalInterpolation<FE>(fe) and success;
770 if (not (disabledTests & DisableRepresentConstants))
772 success = testCanRepresentConstants<FE>(fe) and success;
775 if (not (disabledTests & DisableJacobian))
777 success = testJacobian<FE>(fe, quadOrder, derivativePointSkip) and success;
782 success = (diffOrder == 0) and success;
785 if (not (disabledTests & DisableEvaluate))
787 success = TestPartial::test(fe, TOL, jacobianTOL, diffOrder, quadOrder, derivativePointSkip) and success;
790 if (not (disabledTests & DisableVirtualInterface))
792 typedef typename FE::Traits::LocalBasisType::Traits ImplementationLBTraits;
796 const VirtualFEImp virtualFE(fe);
797 if (not (disabledTests & DisableLocalInterpolation))
798 success = testLocalInterpolation<VirtualFEInterface>(virtualFE) and success;
799 if (not (disabledTests & DisableJacobian))
801 success = testJacobian<VirtualFEInterface>(virtualFE, quadOrder, derivativePointSkip) and success;
806 success = (diffOrder == 0) and success;
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; }
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