Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
#include <config.h>
#include <fenv.h>
#include <iostream>
#include <iomanip>
#include <dune/common/fvector.hh>
#include <dune/geometry/quadraturerules.hh>
#include <dune/geometry/type.hh>
#include <dune/geometry/referenceelements.hh>
#include <dune/localfunctions/lagrange/pqkfactory.hh>
#include <dune/gfe/rotation.hh>
#include <dune/gfe/realtuple.hh>
#include <dune/gfe/unitvector.hh>
#include <dune/gfe/localprojectedfefunction.hh>
#include "multiindex.hh"
#include "valuefactory.hh"
const double eps = 1e-6;
using namespace Dune;
/** \brief Computes the diameter of a set */
template <class TargetSpace>
double diameter(const std::vector<TargetSpace>& v)
{
double d = 0;
for (size_t i=0; i<v.size(); i++)
for (size_t j=0; j<v.size(); j++)
d = std::max(d, TargetSpace::distance(v[i],v[j]));
return d;
}
template <int dim, class ctype, class LocalFunction>
auto
evaluateDerivativeFD(const LocalFunction& f, const Dune::FieldVector<ctype, dim>& local)
-> decltype(f.evaluateDerivative(local))
{
static const int embeddedDim = LocalFunction::TargetSpace::embeddedDim;
Dune::FieldMatrix<ctype, embeddedDim, dim> result;
for (int i=0; i<dim; i++) {
Dune::FieldVector<ctype, dim> forward = local;
Dune::FieldVector<ctype, dim> backward = local;
forward[i] += eps;
backward[i] -= eps;
auto fdDer = f.evaluate(forward).globalCoordinates() - f.evaluate(backward).globalCoordinates();
fdDer /= 2*eps;
for (int j=0; j<embeddedDim; j++)
result[j][i] = fdDer[j];
}
return result;
}

Lisa Julia Nebel
committed
template <int domainDim, int dim>
void testDerivativeTangentiality(const RealTuple<double,dim>& x,
const FieldMatrix<double,dim,domainDim>& derivative)
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
{
// By construction, derivatives of RealTuples are always tangent
}
// the columns of the derivative must be tangential to the manifold
template <int domainDim, int vectorDim>
void testDerivativeTangentiality(const UnitVector<double,vectorDim>& x,
const FieldMatrix<double,vectorDim,domainDim>& derivative)
{
for (int i=0; i<domainDim; i++) {
// The i-th column is a tangent vector if its scalar product with the global coordinates
// of x vanishes.
double sp = 0;
for (int j=0; j<vectorDim; j++)
sp += x.globalCoordinates()[j] * derivative[j][i];
if (std::fabs(sp) > 1e-8)
DUNE_THROW(Dune::Exception, "Derivative is not tangential: Column: " << i << ", product: " << sp);
}
}
// the columns of the derivative must be tangential to the manifold
template <int domainDim, int vectorDim>
void testDerivativeTangentiality(const Rotation<double,vectorDim-1>& x,
const FieldMatrix<double,vectorDim,domainDim>& derivative)
{
}
// the columns of the derivative must be tangential to the manifold
template <int domainDim, int vectorDim>
void testDerivativeTangentiality(const RigidBodyMotion<double,3>& x,
const FieldMatrix<double,vectorDim,domainDim>& derivative)
{
}
/** \brief Test whether interpolation is invariant under permutation of the simplex vertices
* \todo Implement this for all dimensions
*/
template <int domainDim, class TargetSpace>
void testPermutationInvariance(const std::vector<TargetSpace>& corners)
{
// works only for 2d domains
if (domainDim!=2)
return;
PQkLocalFiniteElementCache<double,double,domainDim,1> feCache;
typedef typename PQkLocalFiniteElementCache<double,double,domainDim,1>::FiniteElementType LocalFiniteElement;
GeometryType simplex = GeometryTypes::simplex(domainDim);
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
//
std::vector<TargetSpace> cornersRotated1(domainDim+1);
std::vector<TargetSpace> cornersRotated2(domainDim+1);
cornersRotated1[0] = cornersRotated2[2] = corners[1];
cornersRotated1[1] = cornersRotated2[0] = corners[2];
cornersRotated1[2] = cornersRotated2[1] = corners[0];
GFE::LocalProjectedFEFunction<2,double,LocalFiniteElement,TargetSpace> f0(feCache.get(simplex), corners);
GFE::LocalProjectedFEFunction<2,double,LocalFiniteElement,TargetSpace> f1(feCache.get(simplex), cornersRotated1);
GFE::LocalProjectedFEFunction<2,double,LocalFiniteElement,TargetSpace> f2(feCache.get(simplex), cornersRotated2);
// A quadrature rule as a set of test points
int quadOrder = 3;
const Dune::QuadratureRule<double, domainDim>& quad
= Dune::QuadratureRules<double, domainDim>::rule(simplex, quadOrder);
for (size_t pt=0; pt<quad.size(); pt++) {
const Dune::FieldVector<double,domainDim>& quadPos = quad[pt].position();
Dune::FieldVector<double,domainDim> l0 = quadPos;
Dune::FieldVector<double,domainDim> l1, l2;
l1[0] = quadPos[1];
l1[1] = 1-quadPos[0]-quadPos[1];
l2[0] = 1-quadPos[0]-quadPos[1];
l2[1] = quadPos[0];
// evaluate the three functions
TargetSpace v0 = f0.evaluate(l0);
TargetSpace v1 = f1.evaluate(l1);
TargetSpace v2 = f2.evaluate(l2);
// Check that they are all equal
assert(TargetSpace::distance(v0,v1) < eps);
assert(TargetSpace::distance(v0,v2) < eps);
}
}
template <int domainDim, class TargetSpace>
void testDerivative(const GFE::LocalProjectedFEFunction<domainDim,double,typename PQkLocalFiniteElementCache<double,double,domainDim,1>::FiniteElementType, TargetSpace>& f)
{
static const int embeddedDim = TargetSpace::EmbeddedTangentVector::dimension;
// A quadrature rule as a set of test points
int quadOrder = 3;
const auto& quad = Dune::QuadratureRules<double, domainDim>::rule(f.type(), quadOrder);
for (size_t pt=0; pt<quad.size(); pt++) {
const Dune::FieldVector<double,domainDim>& quadPos = quad[pt].position();
// evaluate actual derivative
Dune::FieldMatrix<double, embeddedDim, domainDim> derivative = f.evaluateDerivative(quadPos);
// evaluate fd approximation of derivative
Dune::FieldMatrix<double, embeddedDim, domainDim> fdDerivative = evaluateDerivativeFD(f,quadPos);
Dune::FieldMatrix<double, embeddedDim, domainDim> diff = derivative;
diff -= fdDerivative;
if ( diff.infinity_norm() > 100*eps ) {
std::cout << className<TargetSpace>() << ": Analytical gradient does not match fd approximation." << std::endl;
std::cout << "Analytical: " << derivative << std::endl;
std::cout << "FD : " << fdDerivative << std::endl;

Lisa Julia Nebel
committed
assert(false);
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
}
testDerivativeTangentiality(f.evaluate(quadPos), derivative);
}
}
template <class TargetSpace, int domainDim>
void test(const GeometryType& element)
{
std::cout << " --- Testing " << className<TargetSpace>() << ", domain dimension: " << element.dim() << " ---" << std::endl;
std::vector<TargetSpace> testPoints;
ValueFactory<TargetSpace>::get(testPoints);
int nTestPoints = testPoints.size();
size_t nVertices = Dune::ReferenceElements<double,domainDim>::general(element).size(domainDim);
// Set up elements of the target space
std::vector<TargetSpace> corners(nVertices);
MultiIndex index(nVertices, nTestPoints);
int numIndices = index.cycle();
for (int i=0; i<numIndices; i++, ++index) {
for (size_t j=0; j<nVertices; j++)
corners[j] = testPoints[index[j]];
if (diameter(corners) > 0.5*M_PI)
continue;
// Make local gfe function to be tested
PQkLocalFiniteElementCache<double,double,domainDim,1> feCache;
typedef typename PQkLocalFiniteElementCache<double,double,domainDim,1>::FiniteElementType LocalFiniteElement;
GFE::LocalProjectedFEFunction<domainDim,double,LocalFiniteElement,TargetSpace> f(feCache.get(element),corners);
//testPermutationInvariance(corners);
testDerivative<domainDim>(f);
}
}
int main()
{
// choke on NaN -- don't enable this by default, as there are
// a few harmless NaN in the loopsolver
//feenableexcept(FE_INVALID);
std::cout << std::setw(15) << std::setprecision(12);
////////////////////////////////////////////////////////////////
// Test functions on 1d elements
////////////////////////////////////////////////////////////////
test<RealTuple<double,1>,1>(GeometryTypes::line);
test<UnitVector<double,2>,1>(GeometryTypes::line);
test<UnitVector<double,3>,1>(GeometryTypes::line);
test<Rotation<double,3>,1>(GeometryTypes::line);
test<RigidBodyMotion<double,3>,1>(GeometryTypes::line);
////////////////////////////////////////////////////////////////
// Test functions on 2d simplex elements
////////////////////////////////////////////////////////////////
test<RealTuple<double,1>,2>(GeometryTypes::triangle);
test<UnitVector<double,2>,2>(GeometryTypes::triangle);

Lisa Julia Nebel
committed
test<RealTuple<double,3>,2>(GeometryTypes::triangle);
test<UnitVector<double,3>,2>(GeometryTypes::triangle);
test<Rotation<double,3>,2>(GeometryTypes::triangle);
test<RigidBodyMotion<double,3>,2>(GeometryTypes::triangle);
////////////////////////////////////////////////////////////////
// Test functions on 2d quadrilateral elements
////////////////////////////////////////////////////////////////
test<RealTuple<double,1>,2>(GeometryTypes::quadrilateral);
test<UnitVector<double,2>,2>(GeometryTypes::quadrilateral);
test<UnitVector<double,3>,2>(GeometryTypes::quadrilateral);
test<Rotation<double,3>,2>(GeometryTypes::quadrilateral);
test<RigidBodyMotion<double,3>,2>(GeometryTypes::quadrilateral);