DUNE-FUNCTIONS (unstable)

basistest.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
5// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
6
7#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_TEST_BASISTEST_HH
8#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_TEST_BASISTEST_HH
9
10#include <set>
11#include <algorithm>
12#include <string>
13#include <sstream>
14#include <map>
15#include <optional>
16
17#include <dune/common/test/testsuite.hh>
18#include <dune/common/concept.hh>
19#include <dune/common/typetraits.hh>
20#include <dune/common/hybridutilities.hh>
21
22#include <dune/geometry/quadraturerules.hh>
23
24#include <dune/functions/functionspacebases/concepts.hh>
25
26struct CheckBasisFlag {};
27struct AllowZeroBasisFunctions {};
28
29template<class T, class... S>
30struct IsContained : public std::disjunction<std::is_same<T,S>...>
31{};
32
33
34
35/*
36 * Get string identifier of element
37 */
38template<class Element, class GridView>
39std::string elementStr(const Element& element, const GridView& gridView)
40{
41 std::stringstream s;
42 s << element.type() << "#" << gridView.indexSet().index(element);
43 return s.str();
44}
45
46/*
47 * Check if two multi-indices are consecutive.
48 * This is a used by checkBasisIndexTreeConsistency()
49 */
50template<class MultiIndex>
51bool multiIndicesConsecutive(const MultiIndex& a, const MultiIndex& b)
52{
53 std::size_t i = 0;
54
55 // find largest common prefix
56 for (; (i<a.size()) and (i<b.size()) and (a[i] == b[i]); ++i)
57 {};
58
59 // if b is exhausted but a is not, then b is a strict prefix of a and does not succeed a
60 if ((i<a.size()) and (i==b.size()))
61 return false;
62
63 // if a and b are not exhausted, then the first non-common index must be an increment
64 if ((i<a.size()) and (i<b.size()))
65 {
66 if (b[i] != a[i]+1)
67 return false;
68 ++i;
69 }
70
71 // if b is not exhausted, then the following indices should be zero
72 if (i<b.size())
73 {
74 for (; i<b.size(); ++i)
75 {
76 if (b[i] != 0)
77 return false;
78 }
79 }
80 return true;
81}
82
83
84
85/*
86 * Check if given set of multi-indices is consistent, i.e.,
87 * if it induces a consistent ordered tree. This is used
88 * by checkBasisIndices()
89 */
90template<class MultiIndexSet>
91Dune::TestSuite checkBasisIndexTreeConsistency(const MultiIndexSet& multiIndexSet)
92{
93 Dune::TestSuite test("index tree consistency check");
94
95 using namespace Dune;
96
97 auto it = multiIndexSet.begin();
98 auto end = multiIndexSet.end();
99
100 // get first multi-index
101 auto lastMultiIndex = *it;
102
103 // assert that index is non-empty
104 test.require(lastMultiIndex.size()>0, "multi-index size check")
105 << "empty multi-index found";
106
107 // check if first multi-index is [0,...,0]
108 for (decltype(lastMultiIndex.size()) i = 0; i<lastMultiIndex.size(); ++i)
109 {
110 test.require(lastMultiIndex[i] == 0, "smallest index check")
111 << "smallest index contains non-zero entry " << lastMultiIndex[i] << " in position " << i;
112 }
113
114 ++it;
115 for(; it != end; ++it)
116 {
117 auto multiIndex = *it;
118
119 // assert that index is non-empty
120 test.require(multiIndex.size()>0, "multi-index size check")
121 << "empty multi-index found";
122
123 // assert that indices are consecutive
124 test.check(multiIndicesConsecutive(lastMultiIndex, multiIndex), "consecutive index check")
125 << "multi-indices " << lastMultiIndex << " and " << multiIndex << " are subsequent but not consecutive";
126
127 lastMultiIndex = multiIndex;
128 }
129
130 return test;
131}
132
133
134
135/*
136 * Check consistency of basis.size(prefix)
137 */
138template<class Basis, class MultiIndexSet>
139Dune::TestSuite checkBasisSizeConsistency(const Basis& basis, const MultiIndexSet& multiIndexSet)
140{
141 Dune::TestSuite test("index size consistency check");
142
143 // Based on the index tree, build a map that contains all possible
144 // prefixes and maps each prefix to the size (of the subsequent digit).
145 using Prefix = typename Basis::SizePrefix;
146 auto prefixSet = std::map<Prefix, std::size_t>();
147 for(const auto& index : multiIndexSet)
148 {
149 auto prefix = Prefix();
150 for (const auto& i: index)
151 {
152 prefixSet[prefix] = std::max(prefixSet[prefix], i+1);
153 prefix.push_back(i);
154 }
155 prefixSet[prefix] = 0;
156 }
157
158 // Now check for all prefixes, if the size computed from the
159 // index tree is consistent with basis.size(prefix).
160 for(const auto& [prefix, size] : prefixSet)
161 {
162 auto prefixSize = basis.size(prefix);
163 test.check(prefixSize == size, "basis.size(prefix) check")
164 << "basis.size(" << prefix << ")=" << prefixSize << ", but should be " << size;
165 }
166
167 return test;
168}
169
170
171
172/*
173 * Check indices of basis:
174 * - First store the whole index tree in a set
175 * - Check if this corresponds to a consistent index tree
176 * - Check if index tree is consistent with basis.size(prefix) and basis.dimension()
177 */
178template<class Basis>
179Dune::TestSuite checkBasisIndices(const Basis& basis)
180{
181 Dune::TestSuite test("basis index check");
182
183 using MultiIndex = typename Basis::MultiIndex;
184
185 static_assert(Dune::IsIndexable<MultiIndex>(), "MultiIndex must support operator[]");
186
187 auto compare = [](const auto& a, const auto& b) {
188 return std::lexicographical_compare(a.begin(), a.end(), b.begin(), b.end());
189 };
190
191 auto multiIndexSet = std::set<MultiIndex, decltype(compare)>{compare};
192
193 auto localView = basis.localView();
194 for (const auto& e : elements(basis.gridView()))
195 {
196 localView.bind(e);
197
198 for (decltype(localView.size()) i=0; i< localView.size(); ++i)
199 {
200 auto multiIndex = localView.index(i);
201 for(auto mi: multiIndex)
202 test.check(mi>=0)
203 << "Global multi-index contains negative entry for shape function " << i
204 << " in element " << elementStr(localView.element(), basis.gridView());
205 multiIndexSet.insert(multiIndex);
206 }
207 }
208
209 test.subTest(checkBasisIndexTreeConsistency(multiIndexSet));
210 test.subTest(checkBasisSizeConsistency(basis, multiIndexSet));
211 test.check(basis.dimension() == multiIndexSet.size())
212 << "basis.dimension() does not equal the total number of basis functions.";
213
214 return test;
215}
216
217
218
219/*
220 * Check if shape functions are not constant zero.
221 * This is called by checkLocalView().
222 */
223template<class LocalFiniteElement>
224Dune::TestSuite checkNonZeroShapeFunctions(const LocalFiniteElement& fe, std::size_t order = 5, double tol = 1e-10)
225{
226 Dune::TestSuite test;
227 static const int dimension = LocalFiniteElement::Traits::LocalBasisType::Traits::dimDomain;
228
229 auto quadRule = Dune::QuadratureRules<double, dimension>::rule(fe.type(), order);
230
231 std::vector<typename LocalFiniteElement::Traits::LocalBasisType::Traits::RangeType> values;
232 std::vector<bool> isNonZero;
233 isNonZero.resize(fe.size(), false);
234 for (const auto& qp : quadRule)
235 {
236 fe.localBasis().evaluateFunction(qp.position(), values);
237 for(std::size_t i=0; i<fe.size(); ++i)
238 isNonZero[i] = (isNonZero[i] or (values[i].infinity_norm() > tol));
239 }
240 for(std::size_t i=0; i<fe.size(); ++i)
241 test.check(isNonZero[i])
242 << "Found a constant zero basis function";
243 return test;
244}
245
246
247
248/*
249 * Check localView. This especially checks for
250 * consistency of local indices and local size.
251 */
252template<class Basis, class LocalView, class... Flags>
253Dune::TestSuite checkLocalView(const Basis& basis, const LocalView& localView, Flags... flags)
254{
255 Dune::TestSuite test(std::string("LocalView on ") + elementStr(localView.element(), basis.gridView()));
256
257 test.check(localView.size() <= localView.maxSize(), "localView.size() check")
258 << "localView.size() is " << localView.size() << " but localView.maxSize() is " << localView.maxSize();
259
260 // Count all local indices appearing in the tree.
261 std::vector<std::size_t> localIndices;
262 localIndices.resize(localView.size(), 0);
263 Dune::TypeTree::forEachLeafNode(localView.tree(), [&](const auto& node, auto&& treePath) {
264 test.check(node.size() == node.finiteElement().size())
265 << "Size of leaf node and finite element are different.";
266 for(std::size_t i=0; i<node.size(); ++i)
267 {
268 test.check(node.localIndex(i) < localView.size())
269 << "Local index exceeds localView.size().";
270 if (node.localIndex(i) < localView.size())
271 ++(localIndices[node.localIndex(i)]);
272 }
273 });
274
275 // Check if each local index appears exactly once.
276 for(std::size_t i=0; i<localView.size(); ++i)
277 {
278 if (localIndices[i])
279 test.check(localIndices[i]>=1)
280 << "Local index " << i << " did not appear";
281 test.check(localIndices[i]<=1)
282 << "Local index " << i << " appears multiple times";
283 }
284
285 // Check that all basis functions are non-zero.
286 if (not IsContained<AllowZeroBasisFunctions, Flags...>::value)
287 {
288 Dune::TypeTree::forEachLeafNode(localView.tree(), [&](const auto& node, auto&& treePath) {
289 test.subTest(checkNonZeroShapeFunctions(node.finiteElement()));
290 });
291 }
292
293 return test;
294}
295
296
297// Flag to enable a local continuity check for checking strong
298// continuity across an intersection within checkBasisContinuity().
299//
300// For each inside basis function this will compute the jump against
301// zero or the corresponding inside basis function. The latter is then
302// checked for being (up to a tolerance) zero on a set of quadrature points.
303struct EnableContinuityCheck
304{
305 std::size_t order_ = 5;
306 double tol_ = 1e-10;
307
308 template<class JumpEvaluator>
309 auto localJumpContinuityCheck(const JumpEvaluator& jumpEvaluator, std::size_t order, double tol) const
310 {
311 return [=](const auto& intersection, const auto& treePath, const auto& insideNode, const auto& outsideNode, const auto& insideToOutside) {
312 using Intersection = std::decay_t<decltype(intersection)>;
313 using Node = std::decay_t<decltype(insideNode)>;
314
315 std::vector<int> isContinuous(insideNode.size(), true);
316 const auto& quadRule = Dune::QuadratureRules<double, Intersection::mydimension>::rule(intersection.type(), order);
317
318 using Range = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
319 std::vector<std::vector<Range>> values;
320 std::vector<std::vector<Range>> neighborValues;
321
322 // Evaluate inside and outside basis functions.
323 values.resize(quadRule.size());
324 neighborValues.resize(quadRule.size());
325 for(std::size_t k=0; k<quadRule.size(); ++k)
326 {
327 auto pointInElement = intersection.geometryInInside().global(quadRule[k].position());
328 auto pointInNeighbor = intersection.geometryInOutside().global(quadRule[k].position());
329 insideNode.finiteElement().localBasis().evaluateFunction(pointInElement, values[k]);
330 outsideNode.finiteElement().localBasis().evaluateFunction(pointInNeighbor, neighborValues[k]);
331 }
332
333 // Check jump against outside basis function or zero.
334 for(std::size_t i=0; i<insideNode.size(); ++i)
335 {
336 for(std::size_t k=0; k<quadRule.size(); ++k)
337 {
338 auto jump = values[k][i];
339 if (insideToOutside[i].has_value())
340 jump -= neighborValues[k][insideToOutside[i].value()];
341 isContinuous[i] = isContinuous[i] and (jumpEvaluator(jump, intersection, quadRule[k].position()) < tol);
342 }
343 }
344 return isContinuous;
345 };
346 }
347
348 auto localContinuityCheck() const {
349 auto jumpNorm = [](auto&&jump, auto&& intersection, auto&& x) -> double {
350 return jump.infinity_norm();
351 };
352 return localJumpContinuityCheck(jumpNorm, order_, tol_);
353 }
354};
355
356// Flag to enable a local normal-continuity check for checking strong
357// continuity across an intersection within checkBasisContinuity().
358//
359// For each inside basis function this will compute the normal jump against
360// zero or the corresponding inside basis function. The latter is then
361// checked for being (up to a tolerance) zero on a set of quadrature points.
362struct EnableNormalContinuityCheck : public EnableContinuityCheck
363{
364 auto localContinuityCheck() const {
365 auto normalJump = [](auto&&jump, auto&& intersection, auto&& x) -> double {
366 return jump * intersection.unitOuterNormal(x);
367 };
368 return localJumpContinuityCheck(normalJump, order_, tol_);
369 }
370};
371
372// Flag to enable a local tangential-continuity check for checking continuity
373// of tangential parts of a vector-valued basis across an intersection
374// within checkBasisContinuity().
375//
376// For each inside basis function this will compute the tangential jump against
377// zero or the corresponding outside basis function. The jump is then
378// checked for being (up to a tolerance) zero on a set of quadrature points.
379struct EnableTangentialContinuityCheck : public EnableContinuityCheck
380{
381 auto localContinuityCheck() const {
382 auto tangentialJumpNorm = [](auto&&jump, auto&& intersection, auto&& x) -> double {
383 auto tangentialJump = jump - (jump * intersection.unitOuterNormal(x)) * intersection.unitOuterNormal(x);
384 return tangentialJump.two_norm();
385 };
386 return localJumpContinuityCheck(tangentialJumpNorm, order_, tol_);
387 }
388};
389
390// Flag to enable a center continuity check for checking continuity in the
391// center of an intersection within checkBasisContinuity().
392//
393// For each inside basis function this will compute the jump against
394// zero or the corresponding inside basis function. The latter is then
395// checked for being (up to a tolerance) zero in the center of mass
396// of the intersection.
397struct EnableCenterContinuityCheck : public EnableContinuityCheck
398{
399 template<class JumpEvaluator>
400 auto localJumpCenterContinuityCheck(const JumpEvaluator& jumpEvaluator, double tol) const
401 {
402 return [=](const auto& intersection, const auto& treePath, const auto& insideNode, const auto& outsideNode, const auto& insideToOutside) {
403 using Node = std::decay_t<decltype(insideNode)>;
404 using Range = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
405
406 std::vector<int> isContinuous(insideNode.size(), true);
407 std::vector<Range> insideValues;
408 std::vector<Range> outsideValues;
409
410 insideNode.finiteElement().localBasis().evaluateFunction(intersection.geometryInInside().center(), insideValues);
411 outsideNode.finiteElement().localBasis().evaluateFunction(intersection.geometryInOutside().center(), outsideValues);
412
413 auto centerLocal = intersection.geometry().local(intersection.geometry().center());
414
415 // Check jump against outside basis function or zero.
416 for(std::size_t i=0; i<insideNode.size(); ++i)
417 {
418 auto jump = insideValues[i];
419 if (insideToOutside[i].has_value())
420 jump -= outsideValues[insideToOutside[i].value()];
421 isContinuous[i] = isContinuous[i] and (jumpEvaluator(jump, intersection, centerLocal) < tol);
422 }
423 return isContinuous;
424 };
425 }
426
427 auto localContinuityCheck() const {
428 auto jumpNorm = [](auto&&jump, auto&& intersection, auto&& x) -> double {
429 return jump.infinity_norm();
430 };
431 return localJumpCenterContinuityCheck(jumpNorm, tol_);
432 }
433};
434
435
436/*
437 * Check if basis functions are continuous across faces.
438 * Continuity is checked by evaluation at a set of quadrature points
439 * from a quadrature rule of given order.
440 * If two basis functions (on neighboring elements) share the same
441 * global index, their values at the quadrature points (located on
442 * their intersection) should coincide up to the given tolerance.
443 *
444 * If a basis function only appears on one side of the intersection,
445 * it should be zero on the intersection.
446 */
447template<class Basis, class LocalCheck>
448Dune::TestSuite checkBasisContinuity(const Basis& basis, const LocalCheck& localCheck)
449{
450 Dune::TestSuite test("Global continuity check of basis functions");
451
452
453 auto localView = basis.localView();
454 auto neighborLocalView = basis.localView();
455
456 for (const auto& e : elements(basis.gridView()))
457 {
458 localView.bind(e);
459 for(const auto& intersection : intersections(basis.gridView(), e))
460 {
461 if (intersection.neighbor())
462 {
463 neighborLocalView.bind(intersection.outside());
464
465 Dune::TypeTree::forEachLeafNode(localView.tree(), [&](const auto& insideNode, auto&& treePath) {
466 const auto& outsideNode = Dune::TypeTree::child(neighborLocalView.tree(), treePath);
467
468 std::vector<std::optional<int>> insideToOutside;
469 insideToOutside.resize(insideNode.size());
470
471 // Map all inside DOFs to outside DOFs if possible
472 for(std::size_t i=0; i<insideNode.size(); ++i)
473 {
474 for(std::size_t j=0; j<outsideNode.size(); ++j)
475 {
476 if (localView.index(insideNode.localIndex(i)) == neighborLocalView.index(outsideNode.localIndex(j)))
477 {
478 // Basis function should only appear once in the neighbor element.
479 test.check(not insideToOutside[i].has_value())
480 << "Basis function " << localView.index(insideNode.localIndex(i))
481 << " appears twice in element " << elementStr(neighborLocalView.element(), basis.gridView());
482 insideToOutside[i] = j;
483 }
484 }
485 }
486
487 // Apply continuity check on given intersection with given inside/outside DOF node pair.
488 auto isContinuous = localCheck(intersection, treePath, insideNode, outsideNode, insideToOutside);
489
490 for(std::size_t i=0; i<insideNode.size(); ++i)
491 {
492 test.check(isContinuous[i])
493 << "Basis function " << localView.index(insideNode.localIndex(i))
494 << " is discontinuous across intersection of elements "
495 << elementStr(localView.element(), basis.gridView())
496 << " and " << elementStr(neighborLocalView.element(), basis.gridView());
497 }
498 });
499 }
500 }
501 }
502 return test;
503}
504
505template<class Basis, class... Flags>
506Dune::TestSuite checkConstBasis(const Basis& basis, Flags... flags)
507{
508 Dune::TestSuite test("const basis check");
509
510 using GridView = typename Basis::GridView;
511
512 // Check if basis models the GlobalBasis concept.
513 test.check(Dune::models<Dune::Functions::Concept::GlobalBasis<GridView>, Basis>(), "global basis concept check")
514 << "type passed to checkBasis() does not model the GlobalBasis concept";
515
516 // Perform all local tests.
517 auto localView = basis.localView();
518 for (const auto& e : elements(basis.gridView()))
519 {
520 localView.bind(e);
521 test.subTest(checkLocalView(basis, localView, flags...));
522 }
523
524 // Perform global index tests.
525 test.subTest(checkBasisIndices(basis));
526
527 // Perform continuity check.
528 // First capture flags in a tuple in order to iterate.
529 auto flagTuple = std::tie(flags...);
530 Dune::Hybrid::forEach(flagTuple, [&](auto&& flag) {
531 using Flag = std::decay_t<decltype(flag)>;
532 if constexpr (std::is_base_of_v<EnableContinuityCheck, Flag>)
533 test.subTest(checkBasisContinuity(basis, flag.localContinuityCheck()));
534 });
535
536 return test;
537}
538
539template<class Basis, class... Flags>
540Dune::TestSuite checkBasis(Basis& basis, Flags... flags)
541{
542 Dune::TestSuite test("basis check");
543
544 // Perform tests for a constant basis
545 test.subTest(checkConstBasis(basis,flags...));
546
547 // Check copy-construction / copy-assignable of a basis
548 {
549 Basis copy(basis);
550 test.subTest(checkConstBasis(copy,flags...));
551 copy = basis;
552 test.subTest(checkConstBasis(copy,flags...));
553 }
554
555 // Check update of gridView
556 auto gridView = basis.gridView();
557 basis.update(gridView);
558
559 return test;
560}
561
562
563
564
565#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_TEST_BASISTEST_HH
Definition: polynomial.hh:17
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Aug 14, 22:29, 2024)