9#ifndef DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
10#define DUNE_LOCALFUNCTIONS_TEST_TEST_FE_HH
33 typedef typename FE::Traits::Basis::Traits::DomainLocal DomainLocal;
34 typedef typename FE::Traits::Basis::Traits::Range Range;
39 typedef typename FE::Traits::Basis::Traits::RangeField CT;
41 std::vector<CT> coeff;
43 FEFunction(
const FE& fe_) : fe(fe_) { resetCoefficients(); }
45 void resetCoefficients() {
46 coeff.resize(fe.basis().size());
47 for(std::size_t i=0; i<coeff.size(); ++i)
51 void setRandom(
double max) {
52 coeff.resize(fe.basis().size());
53 for(std::size_t i=0; i<coeff.size(); ++i)
54 coeff[i] = ((1.0*std::rand()) / RAND_MAX - 0.5)*2.0*max;
57 void evaluate (
const DomainLocal& x, Range& y)
const {
58 std::vector<Range> yy;
59 fe.basis().evaluateFunction(x, yy);
62 for (std::size_t i=0; i<yy.size(); ++i)
63 y.axpy(coeff[i], yy[i]);
81bool testInterpolation(
const FE& fe,
double eps,
int n=5)
86 std::vector<typename FEFunction<FE>::CT> coeff;
87 for(
int i=0; i<n && success; ++i) {
92 fe.interpolation().interpolate(f, coeff);
95 if (coeff.size() != fe.basis().size()) {
96 std::cout <<
"Bug in LocalInterpolation for finite element type "
97 << Dune::className<FE>() <<
":" << std::endl;
98 std::cout <<
" Interpolation vector has size " << coeff.size()
100 std::cout <<
" Basis has size " << fe.basis().size() << std::endl;
101 std::cout << std::endl;
109 for(std::size_t j=0; j<coeff.size() && success; ++j) {
110 if ( std::abs(coeff[j]-f.coeff[j]) >
111 eps*(std::max(std::abs(f.coeff[j]), 1.0)) )
113 std::cout << std::setprecision(16);
114 std::cout <<
"Bug in LocalInterpolation for finite element type "
115 << Dune::className<FE>() <<
":" << std::endl;
116 std::cout <<
" Interpolation weight " << j <<
" differs by "
117 << std::abs(coeff[j]-f.coeff[j]) <<
" from coefficient of "
118 <<
"linear combination." << std::endl;
119 std::cout << std::endl;
139template<
class Geo,
class FE>
140bool testJacobian(
const Geo &geo,
const FE& fe,
double eps,
double delta,
141 std::size_t order = 2)
143 typedef typename FE::Traits::Basis Basis;
145 typedef typename Basis::Traits::DomainField DF;
146 static const std::size_t dimDLocal = Basis::Traits::dimDomainLocal;
147 typedef typename Basis::Traits::DomainLocal DomainLocal;
148 static const std::size_t dimDGlobal = Basis::Traits::dimDomainGlobal;
150 static const std::size_t dimR = Basis::Traits::dimRange;
151 typedef typename Basis::Traits::Range Range;
153 typedef typename Basis::Traits::Jacobian Jacobian;
167 for (std::size_t i=0; i < quad.size(); i++) {
170 const DomainLocal& testPoint = quad[i].position();
173 std::vector<Jacobian> jacobians;
174 fe.basis().evaluateJacobian(testPoint, jacobians);
175 if(jacobians.size() != fe.basis().size()) {
176 std::cout <<
"Bug in evaluateJacobianGlobal() for finite element type "
177 << Dune::className<FE>() <<
":" << std::endl;
178 std::cout <<
" Jacobian vector has size " << jacobians.size()
180 std::cout <<
" Basis has size " << fe.basis().size() << std::endl;
181 std::cout << std::endl;
186 geo.jacobianTransposed(testPoint);
189 for (std::size_t j=0; j<fe.basis().size(); ++j) {
195 for(std::size_t k = 0; k < dimR; ++k)
196 for(std::size_t l = 0; l < dimDGlobal; ++l)
197 for(std::size_t m = 0; m < dimDLocal; ++m)
198 localJacobian[k][m] += jacobians[j][k][l] * geoJT[m][l];
201 for (std::size_t m = 0; m < dimDLocal; ++m) {
204 DomainLocal upPos = testPoint;
205 DomainLocal downPos = testPoint;
210 std::vector<Range> upValues, downValues;
212 fe.basis().evaluateFunction(upPos, upValues);
213 fe.basis().evaluateFunction(downPos, downValues);
216 for(std::size_t k = 0; k < dimR; ++k) {
219 double derivative = localJacobian[k][m];
221 double finiteDiff = (upValues[j][k] - downValues[j][k]) / (2*delta);
224 if ( std::abs(derivative-finiteDiff) >
225 eps/delta*(std::max(std::abs(finiteDiff), 1.0)) )
227 std::cout << std::setprecision(16);
228 std::cout <<
"Bug in evaluateJacobian() for finite element type "
229 << Dune::className<FE>() <<
":" << std::endl;
230 std::cout <<
" Shape function derivative does not agree with "
231 <<
"FD approximation" << std::endl;
232 std::cout <<
" Shape function " << j <<
" component " << k
233 <<
" at position " << testPoint <<
": derivative in "
234 <<
"local direction " << m <<
" is "
235 << derivative <<
", but " << finiteDiff <<
" is "
236 <<
"expected." << std::endl;
237 std::cout << std::endl;
249template<
class Geo,
class FE>
250bool testFE(
const Geo &geo,
const FE& fe,
double eps,
double delta,
255 success = testInterpolation(fe, eps) and success;
256 success = testJacobian(geo, fe, eps, delta, order) and success;
A dense n x m matrix.
Definition: fmatrix.hh:68
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:97
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:225
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 ...