preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
RadialBasisFctMappingTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <memory>
4#include <ostream>
5#include <string>
6#include <utility>
7#include <vector>
8#include "logging/Logger.hpp"
9#include "mapping/Mapping.hpp"
15#include "mesh/Data.hpp"
16#include "mesh/Mesh.hpp"
18#include "mesh/Utils.hpp"
19#include "mesh/Vertex.hpp"
21#include "testing/Testing.hpp"
22
23using namespace precice;
24using namespace precice::mesh;
25using namespace precice::mapping;
26using namespace precice::testing;
28
29BOOST_AUTO_TEST_SUITE(MappingTests)
30BOOST_AUTO_TEST_SUITE(RadialBasisFunctionMapping)
31
33
34
36 int rank;
37 int owner;
40};
41
42/*
43MeshSpecification format:
44{ {rank, owner rank, {x, y, z}, {v}}, ... }
45
46also see struct VertexSpecification.
47
48- -1 on rank means all ranks
49- -1 on owner rank means no rank
50- x, y, z is position of vertex, z is optional, 2D mesh will be created then
51- v is the value of the respective vertex. Only 1D supported at this time.
52
53ReferenceSpecification format:
54{ {rank, {v}, ... }
55- -1 on rank means all ranks
56- v is the expected value of n-th vertex on that particular rank
57*/
63
66
68 MeshSpecification &meshSpec,
69 int globalIndexOffset = 0,
70 bool meshIsSmaller = false)
71{
72 mesh::PtrMesh distributedMesh(new mesh::Mesh(meshSpec.meshName, meshSpec.meshDimension, testing::nextMeshID()));
73 int i = 0;
74 for (auto &vertex : meshSpec.vertices) {
75 if (vertex.rank == context.rank or vertex.rank == -1) {
76 if (vertex.position.size() == 3) // 3-dimensional
77 distributedMesh->createVertex(Eigen::Vector3d(vertex.position.data()));
78 else if (vertex.position.size() == 2) // 2-dimensional
79 distributedMesh->createVertex(Eigen::Vector2d(vertex.position.data()));
80
81 if (vertex.owner == context.rank)
82 distributedMesh->vertices().back().setOwner(true);
83 else
84 distributedMesh->vertices().back().setOwner(false);
85
86 i++;
87 }
88 }
89 addGlobalIndex(distributedMesh, globalIndexOffset);
90 // All tests use eight vertices
91 if (meshIsSmaller) {
92 distributedMesh->setGlobalNumberOfVertices(7);
93 } else {
94 distributedMesh->setGlobalNumberOfVertices(8);
95 }
96 return distributedMesh;
97}
98
99Eigen::VectorXd getDistributedData(const TestContext &context,
100 MeshSpecification const &meshSpec)
101{
102 Eigen::VectorXd d;
103
104 int i = 0;
105 for (auto &vertex : meshSpec.vertices) {
106 if (vertex.rank == context.rank or vertex.rank == -1) {
107 int valueDimension = vertex.value.size();
108 d.conservativeResize(i * valueDimension + valueDimension);
109 // Get data in every dimension
110 for (int dim = 0; dim < valueDimension; ++dim) {
111 d(i * valueDimension + dim) = vertex.value.at(dim);
112 }
113 i++;
114 }
115 }
116 return d;
117}
118
119void testDistributed(const TestContext &context,
120 Mapping &mapping,
121 MeshSpecification inMeshSpec,
122 MeshSpecification outMeshSpec,
123 ReferenceSpecification referenceSpec,
124 int inGlobalIndexOffset = 0,
125 bool meshIsSmaller = false)
126{
127 int valueDimension = inMeshSpec.vertices.at(0).value.size();
128
129 mesh::PtrMesh inMesh = getDistributedMesh(context, inMeshSpec, inGlobalIndexOffset);
130 Eigen::VectorXd inValues = getDistributedData(context, inMeshSpec);
131
132 mesh::PtrMesh outMesh = getDistributedMesh(context, outMeshSpec, 0, meshIsSmaller);
133 Eigen::VectorXd outValues = getDistributedData(context, outMeshSpec);
134 mapping.setMeshes(inMesh, outMesh);
135 BOOST_TEST(mapping.hasComputedMapping() == false);
136
137 mapping.computeMapping();
138 BOOST_TEST(mapping.hasComputedMapping() == true);
139 time::Sample inSample(valueDimension, inValues);
140 mapping.map(inSample, outValues);
141
142 int index = 0;
143 for (auto &referenceVertex : referenceSpec) {
144 if (referenceVertex.first == context.rank or referenceVertex.first == -1) {
145 for (int dim = 0; dim < valueDimension; ++dim) {
146 BOOST_TEST_INFO("Index of vertex: " << index << " - Dimension: " << dim);
147 BOOST_TEST(outValues(index * valueDimension + dim) == referenceVertex.second.at(dim));
148 }
149 ++index;
150 }
151 }
152 BOOST_TEST(outValues.size() == index * valueDimension);
153}
154
155constexpr int meshDims2D{2};
156
158PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
159BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1)
160{
161 PRECICE_TEST();
162 Multiquadrics fct(5.0);
163
165 {-1, 0, {0, 0}, {1}},
166 {-1, 0, {0, 1}, {2}},
167 {-1, 1, {1, 0}, {3}},
168 {-1, 1, {1, 1}, {4}},
169 {-1, 2, {2, 0}, {5}},
170 {-1, 2, {2, 1}, {6}},
171 {-1, 3, {3, 0}, {7}},
172 {-1, 3, {3, 1}, {8}}};
173
175 std::move(inVertexList),
177 "inMesh"};
178
180 {0, -1, {0, 0}, {0}},
181 {0, -1, {0, 1}, {0}},
182 {1, -1, {1, 0}, {0}},
183 {1, -1, {1, 1}, {0}},
184 {2, -1, {2, 0}, {0}},
185 {2, -1, {2, 1}, {0}},
186 {3, -1, {3, 0}, {0}},
187 {3, -1, {3, 1}, {0}}};
188
190 std::move(outVertexList),
192 "outMesh"};
193 ReferenceSpecification ref{// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
194 {0, {1}},
195 {0, {2}},
196 {1, {3}},
197 {1, {4}},
198 {2, {5}},
199 {2, {6}},
200 {3, {7}},
201 {3, {8}}};
202
203 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
204 testDistributed(context, mapping_on, in, out, ref);
205 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
206 testDistributed(context, mapping_sep, in, out, ref);
207 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
208 testDistributed(context, mapping_off, in, out, ref);
209}
210
211PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
212BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1Vector)
213{
214 PRECICE_TEST();
215 Multiquadrics fct(5.0);
216
217 std::vector<VertexSpecification> inVertexList{// Consistent mapping: The inMesh is communicated
218 {-1, 0, {0, 0}, {1, 4}},
219 {-1, 0, {0, 1}, {2, 5}},
220 {-1, 1, {1, 0}, {3, 6}},
221 {-1, 1, {1, 1}, {4, 7}},
222 {-1, 2, {2, 0}, {5, 8}},
223 {-1, 2, {2, 1}, {6, 9}},
224 {-1, 3, {3, 0}, {7, 10}},
225 {-1, 3, {3, 1}, {8, 11}}};
226 MeshSpecification in{// The outMesh is local, distributed among all ranks
227 std::move(inVertexList),
229 "inMesh"};
230
231 std::vector<VertexSpecification> outVertexList{// The outMesh is local, distributed among all ranks
232 {0, -1, {0, 0}, {0, 0}},
233 {0, -1, {0, 1}, {0, 0}},
234 {1, -1, {1, 0}, {0, 0}},
235 {1, -1, {1, 1}, {0, 0}},
236 {2, -1, {2, 0}, {0, 0}},
237 {2, -1, {2, 1}, {0, 0}},
238 {3, -1, {3, 0}, {0, 0}},
239 {3, -1, {3, 1}, {0, 0}}};
241 std::move(outVertexList),
243 "outMesh"};
244 ReferenceSpecification ref{// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
245 {0, {1, 4}},
246 {0, {2, 5}},
247 {1, {3, 6}},
248 {1, {4, 7}},
249 {2, {5, 8}},
250 {2, {6, 9}},
251 {3, {7, 10}},
252 {3, {8, 11}}};
253
254 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
255 testDistributed(context, mapping_on, in, out, ref);
256 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
257 testDistributed(context, mapping_sep, in, out, ref);
258 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
259 testDistributed(context, mapping_off, in, out, ref);
260}
261
263PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
264BOOST_AUTO_TEST_CASE(DistributedConsistent2DV2)
265{
266 PRECICE_TEST();
267 Multiquadrics fct(5.0);
268
269 std::vector<VertexSpecification> inVertexList{// Consistent mapping: The inMesh is communicated, rank 2 owns no vertices
270 {-1, 0, {0, 0}, {1}},
271 {-1, 0, {0, 1}, {2}},
272 {-1, 1, {1, 0}, {3}},
273 {-1, 1, {1, 1}, {4}},
274 {-1, 1, {2, 0}, {5}},
275 {-1, 3, {2, 1}, {6}},
276 {-1, 3, {3, 0}, {7}},
277 {-1, 3, {3, 1}, {8}}};
278
279 MeshSpecification in = {
280 std::move(inVertexList),
282 "inMesh"};
283 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
284 {0, -1, {0, 0}, {0}},
285 {0, -1, {0, 1}, {0}},
286 {0, -1, {1, 0}, {0}},
287 {2, -1, {1, 1}, {0}},
288 {2, -1, {2, 0}, {0}},
289 {2, -1, {2, 1}, {0}},
290 {3, -1, {3, 0}, {0}},
291 {3, -1, {3, 1}, {0}}};
293 std::move(outVertexList),
295 "outMesh"};
296 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
297 // second rank (consistent with the outMesh) is empty, ...
298 {0, {1}},
299 {0, {2}},
300 {0, {3}},
301 {2, {4}},
302 {2, {5}},
303 {2, {6}},
304 {3, {7}},
305 {3, {8}}};
306 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
307 testDistributed(context, mapping_on, in, out, ref);
308 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
309 testDistributed(context, mapping_sep, in, out, ref);
310 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
311 testDistributed(context, mapping_off, in, out, ref);
312}
313
315PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
316BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3)
317{
318 PRECICE_TEST();
319 Multiquadrics fct(5.0);
320
321 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
322
324 // Rank 0 has part of the mesh, owns a subpart
325 {0, 0, {0, 0}, {1}},
326 {0, 0, {0, 1}, {2}},
327 {0, 0, {1, 0}, {3}},
328 {0, -1, {1, 1}, {4}},
329 {0, -1, {2, 0}, {5}},
330 {0, -1, {2, 1}, {6}},
331 // Rank 1 has no vertices
332 // Rank 2 has the entire mesh, but owns just 3 and 5.
333 {2, -1, {0, 0}, {1}},
334 {2, -1, {0, 1}, {2}},
335 {2, -1, {1, 0}, {3}},
336 {2, 2, {1, 1}, {4}},
337 {2, -1, {2, 0}, {5}},
338 {2, 2, {2, 1}, {6}},
339 {2, -1, {3, 0}, {7}},
340 {2, -1, {3, 1}, {8}},
341 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
342 {3, 3, {2, 0}, {5}},
343 {3, -1, {2, 1}, {6}},
344 {3, 3, {3, 0}, {7}},
345 {3, 3, {3, 1}, {8}},
346 };
348 std::move(inVertexList),
350 "inMesh"};
351 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
352 {0, -1, {0, 0}, {0}},
353 {0, -1, {0, 1}, {0}},
354 {0, -1, {1, 0}, {0}},
355 {2, -1, {1, 1}, {0}},
356 {2, -1, {2, 0}, {0}},
357 {2, -1, {2, 1}, {0}},
358 {3, -1, {3, 0}, {0}},
359 {3, -1, {3, 1}, {0}}};
361 std::move(outVertexList),
363 "outMesh"};
364 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
365 // second rank (consistent with the outMesh) is empty, ...
366 {0, {1}},
367 {0, {2}},
368 {0, {3}},
369 {2, {4}},
370 {2, {5}},
371 {2, {6}},
372 {3, {7}},
373 {3, {8}}};
374
375 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
376 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
377 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
378 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
379 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
380 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
381}
382
384PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
385BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3Vector)
386{
387 PRECICE_TEST();
388 Multiquadrics fct(5.0);
389
390 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
391
393 // Rank 0 has part of the mesh, owns a subpart
394 {0, 0, {0, 0}, {1, 4}},
395 {0, 0, {0, 1}, {2, 5}},
396 {0, 0, {1, 0}, {3, 6}},
397 {0, -1, {1, 1}, {4, 7}},
398 {0, -1, {2, 0}, {5, 8}},
399 {0, -1, {2, 1}, {6, 9}},
400 // Rank 1 has no vertices
401 // Rank 2 has the entire mesh, but owns just 3 and 5.
402 {2, -1, {0, 0}, {1, 4}},
403 {2, -1, {0, 1}, {2, 5}},
404 {2, -1, {1, 0}, {3, 6}},
405 {2, 2, {1, 1}, {4, 7}},
406 {2, -1, {2, 0}, {5, 8}},
407 {2, 2, {2, 1}, {6, 9}},
408 {2, -1, {3, 0}, {7, 10}},
409 {2, -1, {3, 1}, {8, 11}},
410 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
411 {3, 3, {2, 0}, {5, 8}},
412 {3, -1, {2, 1}, {6, 9}},
413 {3, 3, {3, 0}, {7, 10}},
414 {3, 3, {3, 1}, {8, 11}},
415 };
417 std::move(inVertexList),
419 "inMesh"};
420 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
421 {0, -1, {0, 0}, {0, 0}},
422 {0, -1, {0, 1}, {0, 0}},
423 {0, -1, {1, 0}, {0, 0}},
424 {2, -1, {1, 1}, {0, 0}},
425 {2, -1, {2, 0}, {0, 0}},
426 {2, -1, {2, 1}, {0, 0}},
427 {3, -1, {3, 0}, {0, 0}},
428 {3, -1, {3, 1}, {0, 0}}};
430 std::move(outVertexList),
432 "outMesh"};
433
434 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
435 // second rank (consistent with the outMesh) is empty, ...
436 {0, {1, 4}},
437 {0, {2, 5}},
438 {0, {3, 6}},
439 {2, {4, 7}},
440 {2, {5, 8}},
441 {2, {6, 9}},
442 {3, {7, 10}},
443 {3, {8, 11}}};
444 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
445 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
446 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
447 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
448 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
449 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
450}
451
453PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
454BOOST_AUTO_TEST_CASE(DistributedConsistent2DV4)
455{
456 PRECICE_TEST();
458
459 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
460
462 // Rank 0 has no vertices
463 // Rank 1 has the entire mesh, owns a subpart
464 {1, 1, {0, 0}, {1.1}},
465 {1, 1, {0, 1}, {2.5}},
466 {1, 1, {1, 0}, {3}},
467 {1, 1, {1, 1}, {4}},
468 {1, -1, {2, 0}, {5}},
469 {1, -1, {2, 1}, {6}},
470 {1, -1, {3, 0}, {7}},
471 {1, -1, {3, 1}, {8}},
472 // Rank 2 has the entire mesh, owns a subpart
473 {2, -1, {0, 0}, {1.1}},
474 {2, -1, {0, 1}, {2.5}},
475 {2, -1, {1, 0}, {3}},
476 {2, -1, {1, 1}, {4}},
477 {2, 2, {2, 0}, {5}},
478 {2, 2, {2, 1}, {6}},
479 {2, 2, {3, 0}, {7}},
480 {2, 2, {3, 1}, {8}},
481 // Rank 3 has no vertices
482 };
484 std::move(inVertexList),
486 "inMesh"};
487 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
488 // not same order as input mesh and vertex (2,0) appears twice
489 {1, -1, {2, 0}, {0}},
490 {1, -1, {1, 0}, {0}},
491 {1, -1, {0, 1}, {0}},
492 {1, -1, {1, 1}, {0}},
493 {1, -1, {0, 0}, {0}},
494 {2, -1, {2, 0}, {0}},
495 {2, -1, {2, 1}, {0}},
496 {2, -1, {3, 0}, {0}},
497 {2, -1, {3, 1}, {0}}};
499 std::move(outVertexList),
501 "outMesh"};
503 {1, {3}},
504 {1, {2.5}},
505 {1, {4}},
506 {1, {1.1}},
507 {2, {5}},
508 {2, {6}},
509 {2, {7}},
510 {2, {8}}};
511 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
512 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
513 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
514 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
515 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
516 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
517}
518
519// same as 2DV4, but all ranks have vertices
520PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
521BOOST_AUTO_TEST_CASE(DistributedConsistent2DV5)
522{
523 PRECICE_TEST();
525
526 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
527 std::vector<VertexSpecification> inVertexList = {
528 // Every rank has the entire mesh and owns a subpart
529 {0, 0, {0, 0}, {1.1}},
530 {0, 0, {0, 1}, {2.5}},
531 {0, -1, {1, 0}, {3}},
532 {0, -1, {1, 1}, {4}},
533 {0, -1, {2, 0}, {5}},
534 {0, -1, {2, 1}, {6}},
535 {0, -1, {3, 0}, {7}},
536 {0, -1, {3, 1}, {8}},
537 {1, -1, {0, 0}, {1.1}},
538 {1, -1, {0, 1}, {2.5}},
539 {1, 1, {1, 0}, {3}},
540 {1, 1, {1, 1}, {4}},
541 {1, -1, {2, 0}, {5}},
542 {1, -1, {2, 1}, {6}},
543 {1, -1, {3, 0}, {7}},
544 {1, -1, {3, 1}, {8}},
545 {2, -1, {0, 0}, {1.1}},
546 {2, -1, {0, 1}, {2.5}},
547 {2, -1, {1, 0}, {3}},
548 {2, -1, {1, 1}, {4}},
549 {2, 2, {2, 0}, {5}},
550 {2, 2, {2, 1}, {6}},
551 {2, -1, {3, 0}, {7}},
552 {2, -1, {3, 1}, {8}},
553 {3, -1, {0, 0}, {1.1}},
554 {3, -1, {0, 1}, {2.5}},
555 {3, -1, {1, 0}, {3}},
556 {3, -1, {1, 1}, {4}},
557 {3, -1, {2, 0}, {5}},
558 {3, -1, {2, 1}, {6}},
559 {3, 3, {3, 0}, {7}},
560 {3, 3, {3, 1}, {8}},
561 };
563 std::move(inVertexList),
565 "inMesh"};
566 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
567 // not same order as input mesh and vertex (2,0) appears twice
568 {1, -1, {2, 0}, {0}},
569 {1, -1, {1, 0}, {0}},
570 {1, -1, {0, 1}, {0}},
571 {1, -1, {1, 1}, {0}},
572 {1, -1, {0, 0}, {0}},
573 {2, -1, {2, 0}, {0}},
574 {2, -1, {2, 1}, {0}},
575 {2, -1, {3, 0}, {0}},
576 {2, -1, {3, 1}, {0}}};
578 std::move(outVertexList),
580 "outMesh"};
582 {1, {3}},
583 {1, {2.5}},
584 {1, {4}},
585 {1, {1.1}},
586 {2, {5}},
587 {2, {6}},
588 {2, {7}},
589 {2, {8}}};
590 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
591 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
592 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
593 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
594 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
595 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
596}
597
599PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
600BOOST_AUTO_TEST_CASE(DistributedConsistent2DV6,
601 *boost::unit_test::tolerance(1e-7))
602{
603 PRECICE_TEST();
605 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
606
608 // Rank 0 has no vertices
609 // Rank 1 has the entire mesh, owns a subpart
610 {1, 1, {0, 0}, {1}},
611 {1, 1, {0, 1}, {2}},
612 {1, 1, {1, 0}, {3}},
613 {1, 1, {1, 1}, {4}},
614 {1, -1, {2, 0}, {5}},
615 {1, -1, {2, 1}, {6}},
616 {1, -1, {3, 0}, {7}},
617 {1, -1, {3, 1}, {8}},
618 // Rank 2 has the entire mesh, owns a subpart
619 {2, -1, {0, 0}, {1}},
620 {2, -1, {0, 1}, {2}},
621 {2, -1, {1, 0}, {3}},
622 {2, -1, {1, 1}, {4}},
623 {2, 2, {2, 0}, {5}},
624 {2, 2, {2, 1}, {6}},
625 {2, 2, {3, 0}, {7}},
626 {2, 2, {3, 1}, {8}},
627 // Rank 3 has no vertices
628 };
630 std::move(inVertexList),
632 "inMesh"};
633 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
634 // not same order as input mesh and vertex (2,0) appears twice
635 {1, -1, {2, 0}, {0}},
636 {1, -1, {1, 0}, {0}},
637 {1, -1, {0, 1}, {0}},
638 {1, -1, {1, 1}, {0}},
639 {1, -1, {0, 0}, {0}},
640 {2, -1, {2, 0}, {0}},
641 {2, -1, {2, 1}, {0}},
642 {2, -1, {3, 0}, {0}},
643 {2, -1, {3, 1}, {0}}};
645 std::move(outVertexList),
647 "outMesh"};
649 {1, {3}},
650 {1, {2}},
651 {1, {4}},
652 {1, {1}},
653 {2, {5}},
654 {2, {6}},
655 {2, {7}},
656 {2, {8}}};
657 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_on(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::ON);
658 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
659 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_sep(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
660 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
661 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_off(Mapping::CONSISTENT, 2, fct, {{false, false, false}}, Polynomial::OFF);
662 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
663}
664
666PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
667BOOST_AUTO_TEST_CASE(DistributedConservative2DV1)
668{
669 PRECICE_TEST();
670 Multiquadrics fct(5.0);
671
672 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
673 {0, -1, {0, 0}, {1}},
674 {0, -1, {0, 1}, {2}},
675 {1, -1, {1, 0}, {3}},
676 {1, -1, {1, 1}, {4}},
677 {2, -1, {2, 0}, {5}},
678 {2, -1, {2, 1}, {6}},
679 {3, -1, {3, 0}, {7}},
680 {3, -1, {3, 1}, {8}}};
682 std::move(inVertexList),
684 "inMesh"};
685 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed
686 {-1, 0, {0, 0}, {0}},
687 {-1, 0, {0, 1}, {0}},
688 {-1, 1, {1, 0}, {0}},
689 {-1, 1, {1, 1}, {0}},
690 {-1, 2, {2, 0}, {0}},
691 {-1, 2, {2, 1}, {0}},
692 {-1, 3, {3, 0}, {0}},
693 {-1, 3, {3, 1}, {0}}};
695 std::move(outVertexList),
697 "outMesh"};
698 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
699 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
700 {0, {1}},
701 {0, {2}},
702 {0, {0}},
703 {0, {0}},
704 {0, {0}},
705 {0, {0}},
706 {0, {0}},
707 {0, {0}},
708 {1, {0}},
709 {1, {0}},
710 {1, {3}},
711 {1, {4}},
712 {1, {0}},
713 {1, {0}},
714 {1, {0}},
715 {1, {0}},
716 {2, {0}},
717 {2, {0}},
718 {2, {0}},
719 {2, {0}},
720 {2, {5}},
721 {2, {6}},
722 {2, {0}},
723 {2, {0}},
724 {3, {0}},
725 {3, {0}},
726 {3, {0}},
727 {3, {0}},
728 {3, {0}},
729 {3, {0}},
730 {3, {7}},
731 {3, {8}}};
732
733 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
734 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
735 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
736 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
737 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
738 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
739}
740
742PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
743BOOST_AUTO_TEST_CASE(DistributedConservative2DV1Vector)
744{
745 PRECICE_TEST();
746 Multiquadrics fct(5.0);
747 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
748 {0, -1, {0, 0}, {1, 4}},
749 {0, -1, {0, 1}, {2, 5}},
750 {1, -1, {1, 0}, {3, 6}},
751 {1, -1, {1, 1}, {4, 7}},
752 {2, -1, {2, 0}, {5, 8}},
753 {2, -1, {2, 1}, {6, 9}},
754 {3, -1, {3, 0}, {7, 10}},
755 {3, -1, {3, 1}, {8, 11}}};
757 std::move(inVertexList),
759 "inMesh"};
760 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed
761 {-1, 0, {0, 0}, {0, 0}},
762 {-1, 0, {0, 1}, {0, 0}},
763 {-1, 1, {1, 0}, {0, 0}},
764 {-1, 1, {1, 1}, {0, 0}},
765 {-1, 2, {2, 0}, {0, 0}},
766 {-1, 2, {2, 1}, {0, 0}},
767 {-1, 3, {3, 0}, {0, 0}},
768 {-1, 3, {3, 1}, {0, 0}}};
770 std::move(outVertexList),
772 "outMesh"};
773 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
774 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
775 {0, {1, 4}},
776 {0, {2, 5}},
777 {0, {0, 0}},
778 {0, {0, 0}},
779 {0, {0, 0}},
780 {0, {0, 0}},
781 {0, {0, 0}},
782 {0, {0, 0}},
783 {1, {0, 0}},
784 {1, {0, 0}},
785 {1, {3, 6}},
786 {1, {4, 7}},
787 {1, {0, 0}},
788 {1, {0, 0}},
789 {1, {0, 0}},
790 {1, {0, 0}},
791 {2, {0, 0}},
792 {2, {0, 0}},
793 {2, {0, 0}},
794 {2, {0, 0}},
795 {2, {5, 8}},
796 {2, {6, 9}},
797 {2, {0, 0}},
798 {2, {0, 0}},
799 {3, {0, 0}},
800 {3, {0, 0}},
801 {3, {0, 0}},
802 {3, {0, 0}},
803 {3, {0, 0}},
804 {3, {0, 0}},
805 {3, {7, 10}},
806 {3, {8, 11}}};
807 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
808 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
809 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
810 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
811 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
812 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
813}
814
816PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
817BOOST_AUTO_TEST_CASE(DistributedConservative2DV2)
818{
819 PRECICE_TEST();
820 Multiquadrics fct(5.0);
821
822 std::vector<int> globalIndexOffsets = {0, 0, 4, 6};
823
824 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local but rank 0 has no vertices
825 {1, -1, {0, 0}, {1}},
826 {1, -1, {0, 1}, {2}},
827 {1, -1, {1, 0}, {3}},
828 {1, -1, {1, 1}, {4}},
829 {2, -1, {2, 0}, {5}},
830 {2, -1, {2, 1}, {6}},
831 {3, -1, {3, 0}, {7}},
832 {3, -1, {3, 1}, {8}}};
834 std::move(inVertexList),
836 "inMesh"};
837 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 owns no vertex
838 {-1, 1, {0, 0}, {0}},
839 {-1, 1, {0, 1}, {0}},
840 {-1, 1, {1, 0}, {0}},
841 {-1, 1, {1, 1}, {0}},
842 {-1, 2, {2, 0}, {0}},
843 {-1, 2, {2, 1}, {0}},
844 {-1, 3, {3, 0}, {0}},
845 {-1, 3, {3, 1}, {0}}};
847 std::move(outVertexList),
849 "outMesh"};
850 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
851 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
852 {0, {0}},
853 {0, {0}},
854 {0, {0}},
855 {0, {0}},
856 {0, {0}},
857 {0, {0}},
858 {0, {0}},
859 {0, {0}},
860 {1, {1}},
861 {1, {2}},
862 {1, {3}},
863 {1, {4}},
864 {1, {0}},
865 {1, {0}},
866 {1, {0}},
867 {1, {0}},
868 {2, {0}},
869 {2, {0}},
870 {2, {0}},
871 {2, {0}},
872 {2, {5}},
873 {2, {6}},
874 {2, {0}},
875 {2, {0}},
876 {3, {0}},
877 {3, {0}},
878 {3, {0}},
879 {3, {0}},
880 {3, {0}},
881 {3, {0}},
882 {3, {7}},
883 {3, {8}}};
884
885 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
886 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
887 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
888 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
889 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
890 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
891}
892
894PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
895BOOST_AUTO_TEST_CASE(DistributedConservative2DV3)
896{
897 PRECICE_TEST();
898 Multiquadrics fct(2.0);
899 std::vector<int> globalIndexOffsets = {0, 0, 3, 5};
900
901 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local but rank 0 has no vertices
902 {1, -1, {0, 0}, {1}},
903 {1, -1, {1, 0}, {3}},
904 {1, -1, {1, 1}, {4}},
905 {2, -1, {2, 0}, {5}},
906 {2, -1, {2, 1}, {6}},
907 {3, -1, {3, 0}, {7}},
908 {3, -1, {3, 1}, {8}}}; // Sum of all vertices is 34
910 std::move(inVertexList),
912 "inMesh"};
913 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 owns no vertex
914 {-1, 1, {0, 0}, {0}},
915 {-1, 1, {0, 1}, {0}},
916 {-1, 1, {1, 0}, {0}},
917 {-1, 1, {1, 1}, {0}},
918 {-1, 2, {2, 0}, {0}},
919 {-1, 2, {2, 1}, {0}},
920 {-1, 3, {3, 0}, {0}},
921 {-1, 3, {3, 1}, {0}}};
923 std::move(outVertexList),
925 "outMesh"};
926 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
927 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
928 {0, {0}},
929 {0, {0}},
930 {0, {0}},
931 {0, {0}},
932 {0, {0}},
933 {0, {0}},
934 {0, {0}},
935 {0, {0}},
936 {1, {1}},
937 {1, {0}},
938 {1, {3}},
939 {1, {4}},
940 {1, {0}},
941 {1, {0}},
942 {1, {0}},
943 {1, {0}},
944 {2, {0}},
945 {2, {0}},
946 {2, {0}},
947 {2, {0}},
948 {2, {5}},
949 {2, {6}},
950 {2, {0}},
951 {2, {0}},
952 {3, {0}},
953 {3, {0}},
954 {3, {0}},
955 {3, {0}},
956 {3, {0}},
957 {3, {0}},
958 {3, {7}},
959 {3, {8}}};
960 // Sum of reference is also 34
961 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
962 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
963 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
964 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
965 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
966 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
967}
968
970PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
971BOOST_AUTO_TEST_CASE(DistributedConservative2DV4,
972 *boost::unit_test::tolerance(1e-6))
973{
974 PRECICE_TEST();
976 std::vector<int> globalIndexOffsets = {0, 2, 4, 6};
977
978 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
979 {0, -1, {0, 0}, {1}},
980 {0, -1, {0, 1}, {2}},
981 {1, -1, {1, 0}, {3}},
982 {1, -1, {1, 1}, {4}},
983 {2, -1, {2, 0}, {5}},
984 {2, -1, {2, 1}, {6}},
985 {3, -1, {3, 0}, {7}},
986 {3, -1, {3, 1}, {8}}}; // Sum is 36
988 std::move(inVertexList),
990 "inMesh"};
991 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 has no vertex at all
992 {-1, 1, {0, 1}, {0}},
993 {-1, 1, {1, 0}, {0}},
994 {-1, 1, {1, 1}, {0}},
995 {-1, 2, {2, 0}, {0}},
996 {-1, 2, {2, 1}, {0}},
997 {-1, 3, {3, 0}, {0}},
998 {-1, 3, {3, 1}, {0}}};
1000 std::move(outVertexList),
1001 meshDims2D,
1002 "outMesh"};
1003 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
1004 // {2, 3, 4, 3, 0, 0, 0, 0} on the second, ...
1005 {0, {0}},
1006 {0, {0}},
1007 {0, {0}},
1008 {0, {0}},
1009 {0, {0}},
1010 {0, {0}},
1011 {0, {0}},
1012 {1, {3.1307987056295525}},
1013 {1, {4.0734031184906971}},
1014 {1, {3.0620533966233006}},
1015 {1, {0}},
1016 {1, {0}},
1017 {1, {0}},
1018 {1, {0}},
1019 {2, {0}},
1020 {2, {0}},
1021 {2, {0}},
1022 {2, {4.4582564956060873}},
1023 {2, {5.8784343572772633}},
1024 {2, {0}},
1025 {2, {0}},
1026 {3, {0}},
1027 {3, {0}},
1028 {3, {0}},
1029 {3, {0}},
1030 {3, {0}},
1031 {3, {7.4683403859032156}},
1032 {3, {7.9287135404698859}}}; // Sum is ~36
1033 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
1034 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank), true);
1035 // Polynomial == OFF won't reach the desired accuracy
1036}
1037
1039PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1040BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5)
1041{
1042 PRECICE_TEST();
1043 Multiquadrics fct(5.0);
1044 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
1045 {0, -1, {0, 0}, {1}},
1046 {0, -1, {0, 1}, {2}},
1047 {1, -1, {1, 0}, {3}},
1048 {1, -1, {1, 1}, {4}},
1049 {2, -1, {2, 0}, {5}},
1050 {2, -1, {2, 1}, {6}},
1051 {3, -1, {3, 0}, {7}},
1052 {3, -1, {3, 1}, {8}}};
1054 std::move(inVertexList),
1055 meshDims2D,
1056 "inMesh"};
1057 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed and non-contigous
1058 {-1, 0, {0, 0}, {0}},
1059 {-1, 1, {0, 1}, {0}},
1060 {-1, 1, {1, 0}, {0}},
1061 {-1, 0, {1, 1}, {0}},
1062 {-1, 2, {2, 0}, {0}},
1063 {-1, 2, {2, 1}, {0}},
1064 {-1, 3, {3, 0}, {0}},
1065 {-1, 3, {3, 1}, {0}}};
1067 std::move(outVertexList),
1068 meshDims2D,
1069 "outMesh"};
1070 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
1071 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
1072 {0, {1}},
1073 {0, {0}},
1074 {0, {0}},
1075 {0, {4}},
1076 {0, {0}},
1077 {0, {0}},
1078 {0, {0}},
1079 {0, {0}},
1080 {1, {0}},
1081 {1, {2}},
1082 {1, {3}},
1083 {1, {0}},
1084 {1, {0}},
1085 {1, {0}},
1086 {1, {0}},
1087 {1, {0}},
1088 {2, {0}},
1089 {2, {0}},
1090 {2, {0}},
1091 {2, {0}},
1092 {2, {5}},
1093 {2, {6}},
1094 {2, {0}},
1095 {2, {0}},
1096 {3, {0}},
1097 {3, {0}},
1098 {3, {0}},
1099 {3, {0}},
1100 {3, {0}},
1101 {3, {0}},
1102 {3, {7}},
1103 {3, {8}}};
1104
1105 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
1106 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
1107 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
1108 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
1109 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
1110 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
1111}
1112
1114PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1115BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5Vector)
1116{
1117 PRECICE_TEST();
1118 Multiquadrics fct(5.0);
1119
1120 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
1121 {0, -1, {0, 0}, {1, 4}},
1122 {0, -1, {0, 1}, {2, 5}},
1123 {1, -1, {1, 0}, {3, 6}},
1124 {1, -1, {1, 1}, {4, 7}},
1125 {2, -1, {2, 0}, {5, 8}},
1126 {2, -1, {2, 1}, {6, 9}},
1127 {3, -1, {3, 0}, {7, 10}},
1128 {3, -1, {3, 1}, {8, 11}}};
1130 std::move(inVertexList),
1131 meshDims2D,
1132 "inMesh"};
1133 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed and non-contigous
1134 {-1, 0, {0, 0}, {0, 0}},
1135 {-1, 1, {0, 1}, {0, 0}},
1136 {-1, 1, {1, 0}, {0, 0}},
1137 {-1, 0, {1, 1}, {0, 0}},
1138 {-1, 2, {2, 0}, {0, 0}},
1139 {-1, 2, {2, 1}, {0, 0}},
1140 {-1, 3, {3, 0}, {0, 0}},
1141 {-1, 3, {3, 1}, {0, 0}}};
1143 std::move(outVertexList),
1144 meshDims2D,
1145 "outMesh"};
1146 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
1147 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
1148 {0, {1, 4}},
1149 {0, {0, 0}},
1150 {0, {0, 0}},
1151 {0, {4, 7}},
1152 {0, {0, 0}},
1153 {0, {0, 0}},
1154 {0, {0, 0}},
1155 {0, {0, 0}},
1156 {1, {0, 0}},
1157 {1, {2, 5}},
1158 {1, {3, 6}},
1159 {1, {0, 0}},
1160 {1, {0, 0}},
1161 {1, {0, 0}},
1162 {1, {0, 0}},
1163 {1, {0, 0}},
1164 {2, {0, 0}},
1165 {2, {0, 0}},
1166 {2, {0, 0}},
1167 {2, {0, 0}},
1168 {2, {5, 8}},
1169 {2, {6, 9}},
1170 {2, {0, 0}},
1171 {2, {0, 0}},
1172 {3, {0, 0}},
1173 {3, {0, 0}},
1174 {3, {0, 0}},
1175 {3, {0, 0}},
1176 {3, {0, 0}},
1177 {3, {0, 0}},
1178 {3, {7, 10}},
1179 {3, {8, 11}}};
1180 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_on(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::ON);
1181 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
1182 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_sep(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
1183 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
1184 RadialBasisFctMapping<RadialBasisFctSolver<Multiquadrics>> mapping_off(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}}, Polynomial::OFF);
1185 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
1186}
1187
1188void testTagging(const TestContext &context,
1189 MeshSpecification inMeshSpec,
1190 MeshSpecification outMeshSpec,
1191 MeshSpecification shouldTagFirstRound,
1192 MeshSpecification shouldTagSecondRound,
1193 bool consistent)
1194{
1195 int meshDimension = inMeshSpec.meshDimension;
1196
1197 mesh::PtrMesh inMesh = getDistributedMesh(context, inMeshSpec);
1198 Eigen::VectorXd inValues = getDistributedData(context, inMeshSpec);
1199
1200 mesh::PtrMesh outMesh = getDistributedMesh(context, outMeshSpec);
1201 Eigen::VectorXd outValues = getDistributedData(context, outMeshSpec);
1202
1203 Gaussian fct(4.5); // Support radius approx. 1
1205 RadialBasisFctMapping<RadialBasisFctSolver<Gaussian>> mapping(constr, 2, fct, {{false, false, false}}, Polynomial::SEPARATE);
1206 inMesh->computeBoundingBox();
1207 outMesh->computeBoundingBox();
1208
1209 mapping.setMeshes(inMesh, outMesh);
1210 mapping.tagMeshFirstRound();
1211
1212 for (const auto &v : inMesh->vertices()) {
1213 auto pos = std::find_if(shouldTagFirstRound.vertices.begin(), shouldTagFirstRound.vertices.end(),
1214 [meshDimension, &v](const VertexSpecification &spec) {
1215 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1216 });
1217 bool found = pos != shouldTagFirstRound.vertices.end();
1218 BOOST_TEST(found >= v.isTagged(),
1219 "FirstRound: Vertex " << v << " is tagged, but should not be.");
1220 BOOST_TEST(found <= v.isTagged(),
1221 "FirstRound: Vertex " << v << " is not tagged, but should be.");
1222 }
1223
1224 mapping.tagMeshSecondRound();
1225
1226 for (const auto &v : inMesh->vertices()) {
1227 auto posFirst = std::find_if(shouldTagFirstRound.vertices.begin(), shouldTagFirstRound.vertices.end(),
1228 [meshDimension, &v](const VertexSpecification &spec) {
1229 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1230 });
1231 bool foundFirst = posFirst != shouldTagFirstRound.vertices.end();
1232 auto posSecond = std::find_if(shouldTagSecondRound.vertices.begin(), shouldTagSecondRound.vertices.end(),
1233 [meshDimension, &v](const VertexSpecification &spec) {
1234 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1235 });
1236 bool foundSecond = posSecond != shouldTagSecondRound.vertices.end();
1237 BOOST_TEST(foundFirst <= v.isTagged(), "SecondRound: Vertex " << v
1238 << " is not tagged, but should be from the first round.");
1239 BOOST_TEST(foundSecond <= v.isTagged(), "SecondRound: Vertex " << v
1240 << " is not tagged, but should be.");
1241 BOOST_TEST((foundSecond or foundFirst) >= v.isTagged(), "SecondRound: Vertex " << v
1242 << " is tagged, but should not be.");
1243 }
1244}
1245
1246PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1247BOOST_AUTO_TEST_CASE(testTagFirstRound)
1248{
1249 PRECICE_TEST();
1250 // *
1251 // + <-- owned
1252 //* * x * *
1253 // *
1254 // *
1256 {0, -1, {0, 0}, {0}}};
1257 MeshSpecification outMeshSpec{
1258 std::move(outVertexList),
1259 meshDims2D,
1260 "outMesh"};
1262 {0, -1, {-1, 0}, {1}}, // inside
1263 {0, -1, {-2, 0}, {1}}, // outside
1264 {0, 0, {1, 0}, {1}}, // inside, owner
1265 {0, -1, {2, 0}, {1}}, // outside
1266 {0, -1, {0, -1}, {1}}, // inside
1267 {0, -1, {0, -2}, {1}}, // outside
1268 {0, -1, {0, 1}, {1}}, // inside
1269 {0, -1, {0, 2}, {1}} // outside
1270 };
1271 MeshSpecification inMeshSpec{
1272 std::move(inVertexList),
1273 meshDims2D,
1274 "inMesh"};
1275 std::vector<VertexSpecification> firstRoundVertices = {
1276 {0, -1, {-1, 0}, {1}},
1277 {0, -1, {1, 0}, {1}},
1278 {0, -1, {0, -1}, {1}},
1279 {0, -1, {0, 1}, {1}}};
1280
1281 MeshSpecification shouldTagFirstRound = {
1282 firstRoundVertices,
1283 static_cast<int>(firstRoundVertices.at(0).position.size()),
1284 ""};
1285 std::vector<VertexSpecification> secondRoundVertices = {
1286 {0, -1, {2, 0}, {1}}};
1287 MeshSpecification shouldTagSecondRound = {
1288 secondRoundVertices,
1289 static_cast<int>(secondRoundVertices.at(0).position.size()),
1290 ""};
1291 testTagging(context, inMeshSpec, outMeshSpec, shouldTagFirstRound, shouldTagSecondRound, true);
1292 // For conservative just swap meshes.
1293 testTagging(context, outMeshSpec, inMeshSpec, shouldTagFirstRound, shouldTagSecondRound, false);
1294}
1295
1296BOOST_AUTO_TEST_SUITE_END() // Parallel
1297
1299
1300#undef doLocalCode
1301#define doLocalCode(Type, function, polynomial) \
1302 { \
1303 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap2D(Mapping::CONSISTENT, 2, function, {{false, false, false}}, polynomial); \
1304 perform2DTestConsistentMapping(consistentMap2D); \
1305 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap2DVector(Mapping::CONSISTENT, 2, function, {{false, false, false}}, polynomial); \
1306 perform2DTestConsistentMappingVector(consistentMap2DVector); \
1307 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap3D(Mapping::CONSISTENT, 3, function, {{false, false, false}}, polynomial); \
1308 perform3DTestConsistentMapping(consistentMap3D); \
1309 RadialBasisFctMapping<RadialBasisFctSolver<Type>> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, function, {{false, false, false}}, polynomial); \
1310 perform2DTestScaledConsistentMapping(scaledConsistentMap2D); \
1311 RadialBasisFctMapping<RadialBasisFctSolver<Type>> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, function, {{false, false, false}}, polynomial); \
1312 perform3DTestScaledConsistentMapping(scaledConsistentMap3D); \
1313 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap2D(Mapping::CONSERVATIVE, 2, function, {{false, false, false}}, polynomial); \
1314 perform2DTestConservativeMapping(conservativeMap2D); \
1315 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap2DVector(Mapping::CONSERVATIVE, 2, function, {{false, false, false}}, polynomial); \
1316 perform2DTestConservativeMappingVector(conservativeMap2DVector); \
1317 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap3D(Mapping::CONSERVATIVE, 3, function, {{false, false, false}}, polynomial); \
1318 perform3DTestConservativeMapping(conservativeMap3D); \
1319 }
1320
1321PRECICE_TEST_SETUP(1_rank)
1322BOOST_AUTO_TEST_CASE(MapThinPlateSplines)
1323{
1324 PRECICE_TEST();
1325 ThinPlateSplines fct;
1326 doLocalCode(ThinPlateSplines, fct, Polynomial::ON);
1327 doLocalCode(ThinPlateSplines, fct, Polynomial::SEPARATE);
1328}
1329
1330PRECICE_TEST_SETUP(1_rank)
1331BOOST_AUTO_TEST_CASE(MapMultiquadrics)
1332{
1333 PRECICE_TEST();
1334 Multiquadrics fct(1e-3);
1335 doLocalCode(Multiquadrics, fct, Polynomial::ON);
1336 doLocalCode(Multiquadrics, fct, Polynomial::SEPARATE);
1337}
1338
1339PRECICE_TEST_SETUP(1_rank)
1340BOOST_AUTO_TEST_CASE(MapInverseMultiquadrics)
1341{
1342 PRECICE_TEST();
1343 InverseMultiquadrics fct(1e-3);
1344 doLocalCode(InverseMultiquadrics, fct, Polynomial::SEPARATE);
1345}
1346
1347PRECICE_TEST_SETUP(1_rank)
1348BOOST_AUTO_TEST_CASE(MapVolumeSplines)
1349{
1350 PRECICE_TEST();
1351 VolumeSplines fct;
1352 doLocalCode(VolumeSplines, fct, Polynomial::ON);
1353 doLocalCode(VolumeSplines, fct, Polynomial::SEPARATE);
1354}
1355
1356PRECICE_TEST_SETUP(1_rank)
1358{
1359 PRECICE_TEST();
1360 Gaussian fct(1.0);
1361 doLocalCode(Gaussian, fct, Polynomial::SEPARATE);
1362}
1363
1364PRECICE_TEST_SETUP(1_rank)
1365BOOST_AUTO_TEST_CASE(MapCompactThinPlateSplinesC2)
1366{
1367 PRECICE_TEST();
1368 double supportRadius = 1.2;
1369 CompactThinPlateSplinesC2 fct(supportRadius);
1370 doLocalCode(CompactThinPlateSplinesC2, fct, Polynomial::SEPARATE);
1371}
1372
1373PRECICE_TEST_SETUP(1_rank)
1374BOOST_AUTO_TEST_CASE(MapCompactPolynomialC0)
1375{
1376 PRECICE_TEST();
1377 double supportRadius = 1.2;
1378 CompactPolynomialC0 fct(supportRadius);
1379 doLocalCode(CompactPolynomialC0, fct, Polynomial::SEPARATE);
1380}
1381
1382PRECICE_TEST_SETUP(1_rank)
1383BOOST_AUTO_TEST_CASE(MapCompactPolynomialC2)
1384{
1385 PRECICE_TEST();
1386 double supportRadius = 1.2;
1387 CompactPolynomialC2 fct(supportRadius);
1388 doLocalCode(CompactPolynomialC2, fct, Polynomial::SEPARATE);
1389}
1390
1391PRECICE_TEST_SETUP(1_rank)
1392BOOST_AUTO_TEST_CASE(MapCompactPolynomialC4)
1393{
1394 PRECICE_TEST();
1395 double supportRadius = 1.2;
1396 CompactPolynomialC4 fct(supportRadius);
1397 doLocalCode(CompactPolynomialC4, fct, Polynomial::SEPARATE);
1398}
1399
1400PRECICE_TEST_SETUP(1_rank)
1401BOOST_AUTO_TEST_CASE(MapCompactPolynomialC6)
1402{
1403 PRECICE_TEST();
1404 double supportRadius = 1.2;
1405 CompactPolynomialC6 fct(supportRadius);
1406 doLocalCode(CompactPolynomialC6, fct, Polynomial::SEPARATE);
1407}
1408
1409PRECICE_TEST_SETUP(1_rank)
1410BOOST_AUTO_TEST_CASE(MapCompactPolynomialC8)
1411{
1412 PRECICE_TEST();
1413 double supportRadius = 1.2;
1414 CompactPolynomialC8 fct(supportRadius);
1415 doLocalCode(CompactPolynomialC8, fct, Polynomial::SEPARATE);
1416}
1417#undef doLocalCode
1418
1420{
1421 using Eigen::Vector2d;
1422 int dimensions = 2;
1423
1424 bool xDead = false;
1425 bool yDead = true;
1426 bool zDead = false;
1427
1428 ThinPlateSplines fct;
1429 RadialBasisFctMapping<RadialBasisFctSolver<ThinPlateSplines>> mapping(constraint, dimensions, fct,
1430 {{xDead, yDead, zDead}}, polynomial);
1431
1432 // Create mesh to map from
1433 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1434 inMesh->createVertex(Vector2d(0.0, 1.0));
1435 inMesh->createVertex(Vector2d(1.0, 1.0));
1436 inMesh->createVertex(Vector2d(2.0, 1.0));
1437 inMesh->createVertex(Vector2d(3.0, 1.0));
1438 addGlobalIndex(inMesh);
1439 inMesh->setGlobalNumberOfVertices(inMesh->nVertices());
1440
1441 Eigen::VectorXd values(4);
1442 values << 1.0, 2.0, 2.0, 1.0;
1443
1444 // Create mesh to map to
1445 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1446 outMesh->createVertex(Vector2d(0, 1.));
1447 outMesh->createVertex(Vector2d(3, 1.));
1448 outMesh->createVertex(Vector2d(1.3, 1.));
1449 outMesh->createVertex(Vector2d(5, 1.));
1450 addGlobalIndex(outMesh);
1451 outMesh->setGlobalNumberOfVertices(outMesh->nVertices());
1452
1453 // Setup mapping with mapping coordinates and geometry used
1454 mapping.setMeshes(inMesh, outMesh);
1455 BOOST_TEST(mapping.hasComputedMapping() == false);
1456
1457 mapping.computeMapping();
1458 Eigen::VectorXd outValues(4);
1459 time::Sample inSample(1, values);
1460 mapping.map(inSample, outValues);
1461 BOOST_TEST(mapping.hasComputedMapping() == true);
1462 if (constraint == Mapping::CONSISTENT) {
1463 if (polynomial == Polynomial::OFF) {
1464 BOOST_TEST(testing::equals(outValues(2), 2.0522549299731567, 1e-7));
1465 } else if (polynomial == Polynomial::SEPARATE) {
1466 BOOST_TEST(testing::equals(outValues(2), 2.0896514371485777, 1e-7));
1467 } else {
1468 BOOST_TEST(testing::equals(outValues(2), 2.1180354377884774, 1e-7));
1469 }
1470 } else {
1471 if (polynomial == Polynomial::OFF) {
1472 BOOST_TEST(testing::equals(outValues(1), 1.8471144693068295, 1e-7));
1473 } else if (polynomial == Polynomial::SEPARATE) {
1474 BOOST_TEST(testing::equals(outValues(1), 1.8236736422730249, 1e-7));
1475 } else {
1476 BOOST_TEST(testing::equals(outValues(1), 1.7587181970483183, 1e-7));
1477 }
1478 }
1479}
1480
1482{
1483 using Eigen::Vector3d;
1484 int dimensions = 3;
1485
1486 ThinPlateSplines fct;
1487 bool xDead = false;
1488 bool yDead = true;
1489 bool zDead = false;
1491 Mapping mapping(constraint, dimensions, fct, {{xDead, yDead, zDead}}, polynomial);
1492
1493 // Create mesh to map from
1494 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1495 inMesh->createVertex(Vector3d(0.0, 3.0, 0.0));
1496 inMesh->createVertex(Vector3d(1.0, 3.0, 0.0));
1497 inMesh->createVertex(Vector3d(0.0, 3.0, 1.0));
1498 inMesh->createVertex(Vector3d(1.0, 3.0, 1.0));
1499 addGlobalIndex(inMesh);
1500 inMesh->setGlobalNumberOfVertices(inMesh->nVertices());
1501
1502 Eigen::VectorXd values(4);
1503 values << 1.0, 2.0, 3.0, 4.0;
1504
1505 // Create mesh to map to
1506 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1507 outMesh->createVertex(Vector3d(0.0, 2.9, 0.0));
1508 outMesh->createVertex(Vector3d(0.8, 2.9, 0.1));
1509 outMesh->createVertex(Vector3d(0.1, 2.9, 0.9));
1510 outMesh->createVertex(Vector3d(1.1, 2.9, 1.1));
1511 addGlobalIndex(outMesh);
1512 outMesh->setGlobalNumberOfVertices(outMesh->nVertices());
1513
1514 // Setup mapping with mapping coordinates and geometry used
1515 mapping.setMeshes(inMesh, outMesh);
1516 BOOST_TEST(mapping.hasComputedMapping() == false);
1517
1518 mapping.computeMapping();
1519 Eigen::VectorXd outValues(4);
1520 time::Sample inSample(1, values);
1521 mapping.map(inSample, outValues);
1522 BOOST_TEST(mapping.hasComputedMapping() == true);
1523
1524 if (constraint == Mapping::CONSISTENT) {
1525 if (polynomial == Polynomial::OFF) {
1526 const double tolerance = 1e-7;
1527 Eigen::VectorXd expected(4);
1528 expected << 1.0, -0.454450524334, 0.99146426249, 6.98958304876;
1529 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1530 } else {
1531 Eigen::VectorXd expected(4);
1532 expected << 1.0, 2.0, 2.9, 4.3;
1533 BOOST_TEST(testing::equals(outValues, expected));
1534 }
1535 } else {
1536 if (polynomial == Polynomial::OFF) {
1537 const double tolerance = 1e-6;
1538 Eigen::VectorXd expected(4);
1539 expected << 1.17251596926, 4.10368825944, 3.56931954192, 3.40160932341;
1540 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1541 } else if (polynomial == Polynomial::ON) {
1542 const double tolerance = 1e-6;
1543 Eigen::VectorXd expected(4);
1544 expected << 0.856701171969, 2.38947124326, 3.34078733786, 3.41304024691;
1545 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1546 } else {
1547 const double tolerance = 1e-6;
1548 Eigen::VectorXd expected(4);
1549 expected << 0.380480856704, 2.83529451713, 3.73088270249, 3.05334192368;
1550 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1551 }
1552 }
1553}
1554
1555PRECICE_TEST_SETUP(1_rank)
1556BOOST_AUTO_TEST_CASE(DeadAxis2Consistent)
1557{
1558 PRECICE_TEST();
1559 testDeadAxis2d(Polynomial::ON, Mapping::CONSISTENT);
1560 testDeadAxis2d(Polynomial::OFF, Mapping::CONSISTENT);
1561 testDeadAxis2d(Polynomial::SEPARATE, Mapping::CONSISTENT);
1562}
1563
1564PRECICE_TEST_SETUP(1_rank)
1565BOOST_AUTO_TEST_CASE(DeadAxis3DConsistent)
1566{
1567 PRECICE_TEST();
1568 testDeadAxis3d(Polynomial::ON, Mapping::CONSISTENT);
1569 testDeadAxis3d(Polynomial::OFF, Mapping::CONSISTENT);
1570 testDeadAxis3d(Polynomial::SEPARATE, Mapping::CONSISTENT);
1571}
1572
1573PRECICE_TEST_SETUP(1_rank)
1574BOOST_AUTO_TEST_CASE(DeadAxis2Conservative)
1575{
1576 PRECICE_TEST();
1577 testDeadAxis2d(Polynomial::ON, Mapping::CONSERVATIVE);
1578 testDeadAxis2d(Polynomial::OFF, Mapping::CONSERVATIVE);
1579 testDeadAxis2d(Polynomial::SEPARATE, Mapping::CONSERVATIVE);
1580}
1581
1582PRECICE_TEST_SETUP(1_rank)
1583BOOST_AUTO_TEST_CASE(DeadAxis3DConervative)
1584{
1585 PRECICE_TEST();
1586 testDeadAxis3d(Polynomial::ON, Mapping::CONSERVATIVE);
1587 testDeadAxis3d(Polynomial::OFF, Mapping::CONSERVATIVE);
1588 testDeadAxis3d(Polynomial::SEPARATE, Mapping::CONSERVATIVE);
1589}
1590
1591BOOST_AUTO_TEST_SUITE_END() // Serial
1592
1594
1595PRECICE_TEST_SETUP(1_rank)
1596BOOST_AUTO_TEST_CASE(inverseTriangularMatrix)
1597{
1598 PRECICE_TEST();
1599 // Define the size of the matrix
1600 const int inSize = 112;
1601
1602 // Generate a random matrix A
1603 Eigen::MatrixXd A = Eigen::MatrixXd::Random(inSize, inSize);
1604
1605 // Create an SPD matrix M
1606 Eigen::MatrixXd M = A * A.transpose();
1607
1608 // Test our custom inversion
1609 Eigen::MatrixXd M_inv_custom = utils::invertLowerTriangularBlockwise(M);
1610
1611 // Second, using Eigen's built-in method
1612 Eigen::MatrixXd M_inv_builtin = Eigen::MatrixXd::Identity(inSize, inSize);
1613 M.triangularView<Eigen::Lower>().solveInPlace(M_inv_builtin);
1614
1615 // Now compare the two inverses and compute the differences
1616 Eigen::MatrixXd diff = M_inv_custom - M_inv_builtin;
1617
1618 double max_abs_diff = diff.maxCoeff();
1619 double min_abs_diff = diff.minCoeff();
1620
1621 // Set a tolerance
1622 double tolerance = 1e-12;
1623
1624 // Use BOOST_CHECK_SMALL to check that the maximum absolute difference is within the tolerance
1625 BOOST_CHECK_SMALL(max_abs_diff, tolerance);
1626 BOOST_CHECK_SMALL(min_abs_diff, tolerance);
1627}
1628
1629BOOST_AUTO_TEST_SUITE_END() // Helper
1630BOOST_AUTO_TEST_SUITE_END() // RadialBasisFunctionMapping
std::ostream & out
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
void addGlobalIndex(mesh::PtrMesh &mesh, int offset=0)
void testTagging(const TestContext &context, MeshSpecification inMeshSpec, MeshSpecification outMeshSpec, MeshSpecification shouldTagFirstRound, MeshSpecification shouldTagSecondRound, bool consistent)
void testDeadAxis3d(Polynomial polynomial, Mapping::Constraint constraint)
void testDistributed(const TestContext &context, Mapping &mapping, MeshSpecification inMeshSpec, MeshSpecification outMeshSpec, ReferenceSpecification referenceSpec, int inGlobalIndexOffset=0, bool meshIsSmaller=false)
void testDeadAxis2d(Polynomial polynomial, Mapping::Constraint constraint)
constexpr int meshDims2D
Eigen::VectorXd getDistributedData(const TestContext &context, MeshSpecification const &meshSpec)
mesh::PtrMesh getDistributedMesh(const TestContext &context, MeshSpecification &meshSpec, int globalIndexOffset=0, bool meshIsSmaller=false)
BOOST_AUTO_TEST_CASE(MapThinPlateSplines)
#define doLocalCode(Type, function, polynomial)
unsigned int index
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
T at(T... args)
T back(T... args)
T begin(T... args)
Wendland radial basis function with compact support.
Wendland radial basis function with compact support.
Wendland radial basis function with compact support.
Wendland radial basis function with compact support.
Wendland radial basis function with compact support.
Radial basis function with compact support.
Radial basis function with global and compact support.
Radial basis function with global support.
Abstract base class for mapping of data from one mesh to another.
Definition Mapping.hpp:16
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
void setMeshes(const mesh::PtrMesh &input, const mesh::PtrMesh &output)
Sets input and output meshes carrying data to be mapped.
Definition Mapping.cpp:28
virtual void computeMapping()=0
Computes the mapping coefficients from the in- and output mesh.
bool hasComputedMapping() const
Returns true, if the mapping has been computed.
Definition Mapping.cpp:253
void map(int inputDataID, int outputDataID)
Definition Mapping.cpp:127
Radial basis function with global support.
Mapping with radial basis functions.
Radial basis function with global support.
Radial basis function with global support.
Container and creator for meshes.
Definition Mesh.hpp:38
void setGlobalNumberOfVertices(int num)
Definition Mesh.hpp:270
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:55
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
void computeBoundingBox()
Computes the boundingBox for the vertices.
Definition Mesh.cpp:244
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
void setOwner(bool owner)
Definition Vertex.cpp:27
Rank rank
the rank of the current participant
T end(T... args)
T find_if(T... args)
contains data mapping from points to meshes.
Polynomial
How to handle the polynomial?
provides Mesh, Data and primitives.
contains the testing framework.
Definition helper.hpp:9
boost::test_tools::predicate_result equals(const std::vector< float > &VectorA, const std::vector< float > &VectorB, float tolerance)
equals to be used in tests. Compares two std::vectors using a given tolerance. Prints both operands o...
Definition Testing.cpp:93
Eigen::MatrixXd invertLowerTriangularBlockwise(const Eigen::MatrixXd &L)
Implements an iterative block scheme to determine the inverse of a lower triangular matrix which is m...
Main namespace of the precice library.
T ref(T... args)
T size(T... args)
std::vector< VertexSpecification > vertices
Holds rank, owner, position and value of a single vertex.