11#ifndef DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
12#define DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
39 using RangeType =
typename FE::Traits::Basis::Traits::Range;
40 using DomainType =
typename FE::Traits::Basis::Traits::DomainLocal;
43 using RangeType =
typename FE::Traits::Basis::Traits::Range;
44 using DomainType =
typename FE::Traits::Basis::Traits::DomainLocal;
47 typedef typename FE::Traits::Basis::Traits::RangeField CT;
49 std::vector<CT> coeff;
51 FEFunction(
const FE& fe_) : fe(fe_) { resetCoefficients(); }
53 void resetCoefficients() {
54 coeff.resize(fe.basis().size());
55 for(std::size_t i=0; i<coeff.size(); ++i)
59 void setRandom(
double max) {
60 coeff.resize(fe.basis().size());
61 for(std::size_t i=0; i<coeff.size(); ++i)
62 coeff[i] = ((1.0*std::rand()) / RAND_MAX - 0.5)*2.0*
max;
65 void evaluate (
const DomainType& x, RangeType& y)
const {
66 std::vector<RangeType> yy;
67 fe.basis().evaluateFunction(x, yy);
69 for (std::size_t i=0; i<yy.size(); ++i)
70 y.axpy(coeff[i], yy[i]);
88bool testInterpolation(
const FE& fe,
double eps,
int n=5)
93 std::vector<typename FEFunction<FE>::CT> coeff;
94 std::vector<typename FEFunction<FE>::CT> coeff2;
95 for(
int i=0; i<n && success; ++i) {
106 fe.interpolation().interpolate(f, coeff);
112 auto callableF = [&](
const auto& x) {
113 using Range =
typename FE::Traits::Basis::Traits::Range;
118 fe.interpolation().interpolate(callableF, coeff2);
119 if (coeff != coeff2) {
120 std::cout <<
"Passing callable and function with evaluate yields different "
121 <<
"interpolation weights for finite element type "
122 << Dune::className<FE>() <<
":" << std::endl;
127 if (coeff.size() != fe.basis().size()) {
128 std::cout <<
"Bug in LocalInterpolation for finite element type "
129 << Dune::className<FE>() <<
":" << std::endl;
130 std::cout <<
" Interpolation vector has size " << coeff.size()
132 std::cout <<
" Basis has size " << fe.basis().size() << std::endl;
133 std::cout << std::endl;
141 for(std::size_t j=0; j<coeff.size() && success; ++j) {
142 if ( std::abs(coeff[j]-f.coeff[j]) >
143 (2*coeff.size()*eps)*(std::max(std::abs(f.coeff[j]), 1.0)) )
145 std::cout << std::setprecision(16);
146 std::cout <<
"Bug in LocalInterpolation for finite element type "
147 << Dune::className<FE>() <<
":" << std::endl;
148 std::cout <<
" Interpolation weight " << j <<
" differs by "
149 << std::abs(coeff[j]-f.coeff[j]) <<
" from coefficient of "
150 <<
"linear combination." << std::endl;
151 std::cout << std::endl;
171template<
class Geo,
class FE>
172bool testJacobian(
const Geo &geo,
const FE& fe,
double eps,
double delta,
173 std::size_t order = 2)
175 typedef typename FE::Traits::Basis Basis;
177 typedef typename Basis::Traits::DomainField DF;
178 static const std::size_t dimDLocal = Basis::Traits::dimDomainLocal;
179 typedef typename Basis::Traits::DomainLocal DomainLocal;
180 static const std::size_t dimDGlobal = Basis::Traits::dimDomainGlobal;
182 static const std::size_t dimR = Basis::Traits::dimRange;
183 typedef typename Basis::Traits::Range Range;
185 typedef typename Basis::Traits::Jacobian Jacobian;
199 for (std::size_t i=0; i < quad.size(); i++) {
202 const DomainLocal& testPoint = quad[i].position();
205 std::vector<Jacobian> jacobians;
206 fe.basis().evaluateJacobian(testPoint, jacobians);
207 if(jacobians.size() != fe.basis().size()) {
208 std::cout <<
"Bug in evaluateJacobianGlobal() for finite element type "
209 << Dune::className<FE>() <<
":" << std::endl;
210 std::cout <<
" Jacobian vector has size " << jacobians.size()
212 std::cout <<
" Basis has size " << fe.basis().size() << std::endl;
213 std::cout << std::endl;
218 geo.jacobianTransposed(testPoint);
221 for (std::size_t j=0; j<fe.basis().size(); ++j) {
227 for(std::size_t k = 0; k < dimR; ++k)
228 for(std::size_t l = 0; l < dimDGlobal; ++l)
229 for(std::size_t m = 0; m < dimDLocal; ++m)
230 localJacobian[k][m] += jacobians[j][k][l] * geoJT[m][l];
233 for (std::size_t m = 0; m < dimDLocal; ++m) {
236 DomainLocal upPos = testPoint;
237 DomainLocal downPos = testPoint;
242 std::vector<Range> upValues, downValues;
244 fe.basis().evaluateFunction(upPos, upValues);
245 fe.basis().evaluateFunction(downPos, downValues);
248 for(std::size_t k = 0; k < dimR; ++k) {
253 double finiteDiff = (upValues[j][k] - downValues[j][k]) / (2*delta);
257 eps/delta*(std::max(std::abs(finiteDiff), 1.0)) )
259 std::cout << std::setprecision(16);
260 std::cout <<
"Bug in evaluateJacobian() for finite element type "
261 << Dune::className<FE>() <<
":" << std::endl;
262 std::cout <<
" Shape function derivative does not agree with "
263 <<
"FD approximation" << std::endl;
264 std::cout <<
" Shape function " << j <<
" component " << k
265 <<
" at position " << testPoint <<
": derivative in "
266 <<
"local direction " << m <<
" is "
267 <<
derivative <<
", but " << finiteDiff <<
" is "
268 <<
"expected." << std::endl;
269 std::cout << std::endl;
281template<
class Geo,
class FE>
282bool testFE(
const Geo &geo,
const FE& fe,
double eps,
double delta,
287 success = testInterpolation(fe, eps) and success;
288 success = testJacobian(geo, fe, eps, delta, order) and success;
A dense n x m matrix.
Definition: fmatrix.hh:117
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.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
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:81