preCICE
Loading...
Searching...
No Matches
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;
38 std::vector<double> position;
39 std::vector<double> value;
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*/
59 std::vector<VertexSpecification> vertices;
61 std::string meshName;
62};
63
65using ReferenceSpecification = std::vector<std::pair<int, std::vector<double>>>;
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 for (auto &vertex : meshSpec.vertices) {
74 if (vertex.rank == context.rank or vertex.rank == -1) {
75 if (vertex.position.size() == 3) // 3-dimensional
76 distributedMesh->createVertex(Eigen::Vector3d(vertex.position.data()));
77 else if (vertex.position.size() == 2) // 2-dimensional
78 distributedMesh->createVertex(Eigen::Vector2d(vertex.position.data()));
79
80 if (vertex.owner == context.rank)
81 distributedMesh->vertices().back().setOwner(true);
82 else
83 distributedMesh->vertices().back().setOwner(false);
84 }
85 }
86 addGlobalIndex(distributedMesh, globalIndexOffset);
87 // All tests use eight vertices
88 if (meshIsSmaller) {
89 distributedMesh->setGlobalNumberOfVertices(7);
90 } else {
91 distributedMesh->setGlobalNumberOfVertices(8);
92 }
93 return distributedMesh;
94}
95
96Eigen::VectorXd getDistributedData(const TestContext &context,
97 MeshSpecification const &meshSpec)
98{
99 Eigen::VectorXd d;
100
101 int i = 0;
102 for (auto &vertex : meshSpec.vertices) {
103 if (vertex.rank == context.rank or vertex.rank == -1) {
104 int valueDimension = vertex.value.size();
105 d.conservativeResize(i * valueDimension + valueDimension);
106 // Get data in every dimension
107 for (int dim = 0; dim < valueDimension; ++dim) {
108 d(i * valueDimension + dim) = vertex.value.at(dim);
109 }
110 i++;
111 }
112 }
113 return d;
114}
115
116void testDistributed(const TestContext &context,
118 MeshSpecification inMeshSpec,
119 MeshSpecification outMeshSpec,
120 ReferenceSpecification referenceSpec,
121 int inGlobalIndexOffset = 0,
122 bool meshIsSmaller = false)
123{
124 int valueDimension = inMeshSpec.vertices.at(0).value.size();
125
126 mesh::PtrMesh inMesh = getDistributedMesh(context, inMeshSpec, inGlobalIndexOffset);
127 Eigen::VectorXd inValues = getDistributedData(context, inMeshSpec);
128
129 mesh::PtrMesh outMesh = getDistributedMesh(context, outMeshSpec, 0, meshIsSmaller);
130 Eigen::VectorXd outValues = getDistributedData(context, outMeshSpec);
131 mapping.setMeshes(inMesh, outMesh);
132 BOOST_TEST(mapping.hasComputedMapping() == false);
133
134 mapping.computeMapping();
135 BOOST_TEST(mapping.hasComputedMapping() == true);
136 time::Sample inSample(valueDimension, inValues);
137 mapping.map(inSample, outValues);
138
139 int index = 0;
140 for (auto &referenceVertex : referenceSpec) {
141 if (referenceVertex.first == context.rank or referenceVertex.first == -1) {
142 for (int dim = 0; dim < valueDimension; ++dim) {
143 BOOST_TEST_INFO("Index of vertex: " << index << " - Dimension: " << dim);
144 BOOST_TEST(outValues(index * valueDimension + dim) == referenceVertex.second.at(dim));
145 }
146 ++index;
147 }
148 }
149 BOOST_TEST(outValues.size() == index * valueDimension);
150}
151
152constexpr int meshDims2D{2};
153
155PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
156BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1)
157{
158 PRECICE_TEST();
159 Multiquadrics fct(5.0);
160
161 std::vector<VertexSpecification> inVertexList{
162 {-1, 0, {0, 0}, {1}},
163 {-1, 0, {0, 1}, {2}},
164 {-1, 1, {1, 0}, {3}},
165 {-1, 1, {1, 1}, {4}},
166 {-1, 2, {2, 0}, {5}},
167 {-1, 2, {2, 1}, {6}},
168 {-1, 3, {3, 0}, {7}},
169 {-1, 3, {3, 1}, {8}}};
170
172 std::move(inVertexList),
174 "inMesh"};
175
176 std::vector<VertexSpecification> outVertexList{
177 {0, -1, {0, 0}, {0}},
178 {0, -1, {0, 1}, {0}},
179 {1, -1, {1, 0}, {0}},
180 {1, -1, {1, 1}, {0}},
181 {2, -1, {2, 0}, {0}},
182 {2, -1, {2, 1}, {0}},
183 {3, -1, {3, 0}, {0}},
184 {3, -1, {3, 1}, {0}}};
185
187 std::move(outVertexList),
189 "outMesh"};
190 ReferenceSpecification ref{// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
191 {0, {1}},
192 {0, {2}},
193 {1, {3}},
194 {1, {4}},
195 {2, {5}},
196 {2, {6}},
197 {3, {7}},
198 {3, {8}}};
199
201 testDistributed(context, mapping_on, in, out, ref);
203 testDistributed(context, mapping_sep, in, out, ref);
205 testDistributed(context, mapping_off, in, out, ref);
206}
207
208PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
209BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1Vector)
210{
211 PRECICE_TEST();
212 Multiquadrics fct(5.0);
213
214 std::vector<VertexSpecification> inVertexList{// Consistent mapping: The inMesh is communicated
215 {-1, 0, {0, 0}, {1, 4}},
216 {-1, 0, {0, 1}, {2, 5}},
217 {-1, 1, {1, 0}, {3, 6}},
218 {-1, 1, {1, 1}, {4, 7}},
219 {-1, 2, {2, 0}, {5, 8}},
220 {-1, 2, {2, 1}, {6, 9}},
221 {-1, 3, {3, 0}, {7, 10}},
222 {-1, 3, {3, 1}, {8, 11}}};
223 MeshSpecification in{// The outMesh is local, distributed among all ranks
224 std::move(inVertexList),
226 "inMesh"};
227
228 std::vector<VertexSpecification> outVertexList{// The outMesh is local, distributed among all ranks
229 {0, -1, {0, 0}, {0, 0}},
230 {0, -1, {0, 1}, {0, 0}},
231 {1, -1, {1, 0}, {0, 0}},
232 {1, -1, {1, 1}, {0, 0}},
233 {2, -1, {2, 0}, {0, 0}},
234 {2, -1, {2, 1}, {0, 0}},
235 {3, -1, {3, 0}, {0, 0}},
236 {3, -1, {3, 1}, {0, 0}}};
238 std::move(outVertexList),
240 "outMesh"};
241 ReferenceSpecification ref{// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
242 {0, {1, 4}},
243 {0, {2, 5}},
244 {1, {3, 6}},
245 {1, {4, 7}},
246 {2, {5, 8}},
247 {2, {6, 9}},
248 {3, {7, 10}},
249 {3, {8, 11}}};
250
252 testDistributed(context, mapping_on, in, out, ref);
254 testDistributed(context, mapping_sep, in, out, ref);
256 testDistributed(context, mapping_off, in, out, ref);
257}
258
260PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
261BOOST_AUTO_TEST_CASE(DistributedConsistent2DV2)
262{
263 PRECICE_TEST();
264 Multiquadrics fct(5.0);
265
266 std::vector<VertexSpecification> inVertexList{// Consistent mapping: The inMesh is communicated, rank 2 owns no vertices
267 {-1, 0, {0, 0}, {1}},
268 {-1, 0, {0, 1}, {2}},
269 {-1, 1, {1, 0}, {3}},
270 {-1, 1, {1, 1}, {4}},
271 {-1, 1, {2, 0}, {5}},
272 {-1, 3, {2, 1}, {6}},
273 {-1, 3, {3, 0}, {7}},
274 {-1, 3, {3, 1}, {8}}};
275
276 MeshSpecification in = {
277 std::move(inVertexList),
279 "inMesh"};
280 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
281 {0, -1, {0, 0}, {0}},
282 {0, -1, {0, 1}, {0}},
283 {0, -1, {1, 0}, {0}},
284 {2, -1, {1, 1}, {0}},
285 {2, -1, {2, 0}, {0}},
286 {2, -1, {2, 1}, {0}},
287 {3, -1, {3, 0}, {0}},
288 {3, -1, {3, 1}, {0}}};
289 MeshSpecification out = {
290 std::move(outVertexList),
292 "outMesh"};
293 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
294 // second rank (consistent with the outMesh) is empty, ...
295 {0, {1}},
296 {0, {2}},
297 {0, {3}},
298 {2, {4}},
299 {2, {5}},
300 {2, {6}},
301 {3, {7}},
302 {3, {8}}};
304 testDistributed(context, mapping_on, in, out, ref);
306 testDistributed(context, mapping_sep, in, out, ref);
308 testDistributed(context, mapping_off, in, out, ref);
309}
310
312PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
313BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3)
314{
315 PRECICE_TEST();
316 Multiquadrics fct(5.0);
317
318 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
319
320 std::vector<VertexSpecification> inVertexList{
321 // Rank 0 has part of the mesh, owns a subpart
322 {0, 0, {0, 0}, {1}},
323 {0, 0, {0, 1}, {2}},
324 {0, 0, {1, 0}, {3}},
325 {0, -1, {1, 1}, {4}},
326 {0, -1, {2, 0}, {5}},
327 {0, -1, {2, 1}, {6}},
328 // Rank 1 has no vertices
329 // Rank 2 has the entire mesh, but owns just 3 and 5.
330 {2, -1, {0, 0}, {1}},
331 {2, -1, {0, 1}, {2}},
332 {2, -1, {1, 0}, {3}},
333 {2, 2, {1, 1}, {4}},
334 {2, -1, {2, 0}, {5}},
335 {2, 2, {2, 1}, {6}},
336 {2, -1, {3, 0}, {7}},
337 {2, -1, {3, 1}, {8}},
338 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
339 {3, 3, {2, 0}, {5}},
340 {3, -1, {2, 1}, {6}},
341 {3, 3, {3, 0}, {7}},
342 {3, 3, {3, 1}, {8}},
343 };
345 std::move(inVertexList),
347 "inMesh"};
348 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
349 {0, -1, {0, 0}, {0}},
350 {0, -1, {0, 1}, {0}},
351 {0, -1, {1, 0}, {0}},
352 {2, -1, {1, 1}, {0}},
353 {2, -1, {2, 0}, {0}},
354 {2, -1, {2, 1}, {0}},
355 {3, -1, {3, 0}, {0}},
356 {3, -1, {3, 1}, {0}}};
358 std::move(outVertexList),
360 "outMesh"};
361 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
362 // second rank (consistent with the outMesh) is empty, ...
363 {0, {1}},
364 {0, {2}},
365 {0, {3}},
366 {2, {4}},
367 {2, {5}},
368 {2, {6}},
369 {3, {7}},
370 {3, {8}}};
371
373 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
375 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
377 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
378}
379
381PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
382BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3Vector)
383{
384 PRECICE_TEST();
385 Multiquadrics fct(5.0);
386
387 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
388
389 std::vector<VertexSpecification> inVertexList{
390 // Rank 0 has part of the mesh, owns a subpart
391 {0, 0, {0, 0}, {1, 4}},
392 {0, 0, {0, 1}, {2, 5}},
393 {0, 0, {1, 0}, {3, 6}},
394 {0, -1, {1, 1}, {4, 7}},
395 {0, -1, {2, 0}, {5, 8}},
396 {0, -1, {2, 1}, {6, 9}},
397 // Rank 1 has no vertices
398 // Rank 2 has the entire mesh, but owns just 3 and 5.
399 {2, -1, {0, 0}, {1, 4}},
400 {2, -1, {0, 1}, {2, 5}},
401 {2, -1, {1, 0}, {3, 6}},
402 {2, 2, {1, 1}, {4, 7}},
403 {2, -1, {2, 0}, {5, 8}},
404 {2, 2, {2, 1}, {6, 9}},
405 {2, -1, {3, 0}, {7, 10}},
406 {2, -1, {3, 1}, {8, 11}},
407 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
408 {3, 3, {2, 0}, {5, 8}},
409 {3, -1, {2, 1}, {6, 9}},
410 {3, 3, {3, 0}, {7, 10}},
411 {3, 3, {3, 1}, {8, 11}},
412 };
414 std::move(inVertexList),
416 "inMesh"};
417 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 1 is empty
418 {0, -1, {0, 0}, {0, 0}},
419 {0, -1, {0, 1}, {0, 0}},
420 {0, -1, {1, 0}, {0, 0}},
421 {2, -1, {1, 1}, {0, 0}},
422 {2, -1, {2, 0}, {0, 0}},
423 {2, -1, {2, 1}, {0, 0}},
424 {3, -1, {3, 0}, {0, 0}},
425 {3, -1, {3, 1}, {0, 0}}};
427 std::move(outVertexList),
429 "outMesh"};
430
431 ReferenceSpecification ref{// Tests for {0, 1, 2} on the first rank,
432 // second rank (consistent with the outMesh) is empty, ...
433 {0, {1, 4}},
434 {0, {2, 5}},
435 {0, {3, 6}},
436 {2, {4, 7}},
437 {2, {5, 8}},
438 {2, {6, 9}},
439 {3, {7, 10}},
440 {3, {8, 11}}};
442 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
444 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
446 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
447}
448
450PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
451BOOST_AUTO_TEST_CASE(DistributedConsistent2DV4)
452{
453 PRECICE_TEST();
455
456 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
457
458 std::vector<VertexSpecification> inVertexList{
459 // Rank 0 has no vertices
460 // Rank 1 has the entire mesh, owns a subpart
461 {1, 1, {0, 0}, {1.1}},
462 {1, 1, {0, 1}, {2.5}},
463 {1, 1, {1, 0}, {3}},
464 {1, 1, {1, 1}, {4}},
465 {1, -1, {2, 0}, {5}},
466 {1, -1, {2, 1}, {6}},
467 {1, -1, {3, 0}, {7}},
468 {1, -1, {3, 1}, {8}},
469 // Rank 2 has the entire mesh, owns a subpart
470 {2, -1, {0, 0}, {1.1}},
471 {2, -1, {0, 1}, {2.5}},
472 {2, -1, {1, 0}, {3}},
473 {2, -1, {1, 1}, {4}},
474 {2, 2, {2, 0}, {5}},
475 {2, 2, {2, 1}, {6}},
476 {2, 2, {3, 0}, {7}},
477 {2, 2, {3, 1}, {8}},
478 // Rank 3 has no vertices
479 };
481 std::move(inVertexList),
483 "inMesh"};
484 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
485 // not same order as input mesh and vertex (2,0) appears twice
486 {1, -1, {2, 0}, {0}},
487 {1, -1, {1, 0}, {0}},
488 {1, -1, {0, 1}, {0}},
489 {1, -1, {1, 1}, {0}},
490 {1, -1, {0, 0}, {0}},
491 {2, -1, {2, 0}, {0}},
492 {2, -1, {2, 1}, {0}},
493 {2, -1, {3, 0}, {0}},
494 {2, -1, {3, 1}, {0}}};
496 std::move(outVertexList),
498 "outMesh"};
499 ReferenceSpecification ref{{1, {5}},
500 {1, {3}},
501 {1, {2.5}},
502 {1, {4}},
503 {1, {1.1}},
504 {2, {5}},
505 {2, {6}},
506 {2, {7}},
507 {2, {8}}};
509 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
511 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
513 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
514}
515
516// same as 2DV4, but all ranks have vertices
517PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
518BOOST_AUTO_TEST_CASE(DistributedConsistent2DV5)
519{
520 PRECICE_TEST();
522
523 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
524 std::vector<VertexSpecification> inVertexList = {
525 // Every rank has the entire mesh and owns a subpart
526 {0, 0, {0, 0}, {1.1}},
527 {0, 0, {0, 1}, {2.5}},
528 {0, -1, {1, 0}, {3}},
529 {0, -1, {1, 1}, {4}},
530 {0, -1, {2, 0}, {5}},
531 {0, -1, {2, 1}, {6}},
532 {0, -1, {3, 0}, {7}},
533 {0, -1, {3, 1}, {8}},
534 {1, -1, {0, 0}, {1.1}},
535 {1, -1, {0, 1}, {2.5}},
536 {1, 1, {1, 0}, {3}},
537 {1, 1, {1, 1}, {4}},
538 {1, -1, {2, 0}, {5}},
539 {1, -1, {2, 1}, {6}},
540 {1, -1, {3, 0}, {7}},
541 {1, -1, {3, 1}, {8}},
542 {2, -1, {0, 0}, {1.1}},
543 {2, -1, {0, 1}, {2.5}},
544 {2, -1, {1, 0}, {3}},
545 {2, -1, {1, 1}, {4}},
546 {2, 2, {2, 0}, {5}},
547 {2, 2, {2, 1}, {6}},
548 {2, -1, {3, 0}, {7}},
549 {2, -1, {3, 1}, {8}},
550 {3, -1, {0, 0}, {1.1}},
551 {3, -1, {0, 1}, {2.5}},
552 {3, -1, {1, 0}, {3}},
553 {3, -1, {1, 1}, {4}},
554 {3, -1, {2, 0}, {5}},
555 {3, -1, {2, 1}, {6}},
556 {3, 3, {3, 0}, {7}},
557 {3, 3, {3, 1}, {8}},
558 };
560 std::move(inVertexList),
562 "inMesh"};
563 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
564 // not same order as input mesh and vertex (2,0) appears twice
565 {1, -1, {2, 0}, {0}},
566 {1, -1, {1, 0}, {0}},
567 {1, -1, {0, 1}, {0}},
568 {1, -1, {1, 1}, {0}},
569 {1, -1, {0, 0}, {0}},
570 {2, -1, {2, 0}, {0}},
571 {2, -1, {2, 1}, {0}},
572 {2, -1, {3, 0}, {0}},
573 {2, -1, {3, 1}, {0}}};
575 std::move(outVertexList),
577 "outMesh"};
578 ReferenceSpecification ref{{1, {5}},
579 {1, {3}},
580 {1, {2.5}},
581 {1, {4}},
582 {1, {1.1}},
583 {2, {5}},
584 {2, {6}},
585 {2, {7}},
586 {2, {8}}};
588 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
590 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
592 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
593}
594
596PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
597BOOST_AUTO_TEST_CASE(DistributedConsistent2DV6,
598 *boost::unit_test::tolerance(1e-7))
599{
600 PRECICE_TEST();
602 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
603
604 std::vector<VertexSpecification> inVertexList{
605 // Rank 0 has no vertices
606 // Rank 1 has the entire mesh, owns a subpart
607 {1, 1, {0, 0}, {1}},
608 {1, 1, {0, 1}, {2}},
609 {1, 1, {1, 0}, {3}},
610 {1, 1, {1, 1}, {4}},
611 {1, -1, {2, 0}, {5}},
612 {1, -1, {2, 1}, {6}},
613 {1, -1, {3, 0}, {7}},
614 {1, -1, {3, 1}, {8}},
615 // Rank 2 has the entire mesh, owns a subpart
616 {2, -1, {0, 0}, {1}},
617 {2, -1, {0, 1}, {2}},
618 {2, -1, {1, 0}, {3}},
619 {2, -1, {1, 1}, {4}},
620 {2, 2, {2, 0}, {5}},
621 {2, 2, {2, 1}, {6}},
622 {2, 2, {3, 0}, {7}},
623 {2, 2, {3, 1}, {8}},
624 // Rank 3 has no vertices
625 };
627 std::move(inVertexList),
629 "inMesh"};
630 std::vector<VertexSpecification> outVertexList{// The outMesh is local, rank 0 and 3 are empty
631 // not same order as input mesh and vertex (2,0) appears twice
632 {1, -1, {2, 0}, {0}},
633 {1, -1, {1, 0}, {0}},
634 {1, -1, {0, 1}, {0}},
635 {1, -1, {1, 1}, {0}},
636 {1, -1, {0, 0}, {0}},
637 {2, -1, {2, 0}, {0}},
638 {2, -1, {2, 1}, {0}},
639 {2, -1, {3, 0}, {0}},
640 {2, -1, {3, 1}, {0}}};
642 std::move(outVertexList),
644 "outMesh"};
645 ReferenceSpecification ref{{1, {5}},
646 {1, {3}},
647 {1, {2}},
648 {1, {4}},
649 {1, {1}},
650 {2, {5}},
651 {2, {6}},
652 {2, {7}},
653 {2, {8}}};
655 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
657 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
659 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
660}
661
663PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
664BOOST_AUTO_TEST_CASE(DistributedConservative2DV1)
665{
666 PRECICE_TEST();
667 Multiquadrics fct(5.0);
668
669 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
670 {0, -1, {0, 0}, {1}},
671 {0, -1, {0, 1}, {2}},
672 {1, -1, {1, 0}, {3}},
673 {1, -1, {1, 1}, {4}},
674 {2, -1, {2, 0}, {5}},
675 {2, -1, {2, 1}, {6}},
676 {3, -1, {3, 0}, {7}},
677 {3, -1, {3, 1}, {8}}};
679 std::move(inVertexList),
681 "inMesh"};
682 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed
683 {-1, 0, {0, 0}, {0}},
684 {-1, 0, {0, 1}, {0}},
685 {-1, 1, {1, 0}, {0}},
686 {-1, 1, {1, 1}, {0}},
687 {-1, 2, {2, 0}, {0}},
688 {-1, 2, {2, 1}, {0}},
689 {-1, 3, {3, 0}, {0}},
690 {-1, 3, {3, 1}, {0}}};
692 std::move(outVertexList),
694 "outMesh"};
695 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
696 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
697 {0, {1}},
698 {0, {2}},
699 {0, {0}},
700 {0, {0}},
701 {0, {0}},
702 {0, {0}},
703 {0, {0}},
704 {0, {0}},
705 {1, {0}},
706 {1, {0}},
707 {1, {3}},
708 {1, {4}},
709 {1, {0}},
710 {1, {0}},
711 {1, {0}},
712 {1, {0}},
713 {2, {0}},
714 {2, {0}},
715 {2, {0}},
716 {2, {0}},
717 {2, {5}},
718 {2, {6}},
719 {2, {0}},
720 {2, {0}},
721 {3, {0}},
722 {3, {0}},
723 {3, {0}},
724 {3, {0}},
725 {3, {0}},
726 {3, {0}},
727 {3, {7}},
728 {3, {8}}};
729
731 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
733 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
735 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
736}
737
739PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
740BOOST_AUTO_TEST_CASE(DistributedConservative2DV1Vector)
741{
742 PRECICE_TEST();
743 Multiquadrics fct(5.0);
744 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
745 {0, -1, {0, 0}, {1, 4}},
746 {0, -1, {0, 1}, {2, 5}},
747 {1, -1, {1, 0}, {3, 6}},
748 {1, -1, {1, 1}, {4, 7}},
749 {2, -1, {2, 0}, {5, 8}},
750 {2, -1, {2, 1}, {6, 9}},
751 {3, -1, {3, 0}, {7, 10}},
752 {3, -1, {3, 1}, {8, 11}}};
754 std::move(inVertexList),
756 "inMesh"};
757 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed
758 {-1, 0, {0, 0}, {0, 0}},
759 {-1, 0, {0, 1}, {0, 0}},
760 {-1, 1, {1, 0}, {0, 0}},
761 {-1, 1, {1, 1}, {0, 0}},
762 {-1, 2, {2, 0}, {0, 0}},
763 {-1, 2, {2, 1}, {0, 0}},
764 {-1, 3, {3, 0}, {0, 0}},
765 {-1, 3, {3, 1}, {0, 0}}};
767 std::move(outVertexList),
769 "outMesh"};
770 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
771 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
772 {0, {1, 4}},
773 {0, {2, 5}},
774 {0, {0, 0}},
775 {0, {0, 0}},
776 {0, {0, 0}},
777 {0, {0, 0}},
778 {0, {0, 0}},
779 {0, {0, 0}},
780 {1, {0, 0}},
781 {1, {0, 0}},
782 {1, {3, 6}},
783 {1, {4, 7}},
784 {1, {0, 0}},
785 {1, {0, 0}},
786 {1, {0, 0}},
787 {1, {0, 0}},
788 {2, {0, 0}},
789 {2, {0, 0}},
790 {2, {0, 0}},
791 {2, {0, 0}},
792 {2, {5, 8}},
793 {2, {6, 9}},
794 {2, {0, 0}},
795 {2, {0, 0}},
796 {3, {0, 0}},
797 {3, {0, 0}},
798 {3, {0, 0}},
799 {3, {0, 0}},
800 {3, {0, 0}},
801 {3, {0, 0}},
802 {3, {7, 10}},
803 {3, {8, 11}}};
805 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
807 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
809 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
810}
811
813PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
814BOOST_AUTO_TEST_CASE(DistributedConservative2DV2)
815{
816 PRECICE_TEST();
817 Multiquadrics fct(5.0);
818
819 std::vector<int> globalIndexOffsets = {0, 0, 4, 6};
820
821 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local but rank 0 has no vertices
822 {1, -1, {0, 0}, {1}},
823 {1, -1, {0, 1}, {2}},
824 {1, -1, {1, 0}, {3}},
825 {1, -1, {1, 1}, {4}},
826 {2, -1, {2, 0}, {5}},
827 {2, -1, {2, 1}, {6}},
828 {3, -1, {3, 0}, {7}},
829 {3, -1, {3, 1}, {8}}};
831 std::move(inVertexList),
833 "inMesh"};
834 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 owns no vertex
835 {-1, 1, {0, 0}, {0}},
836 {-1, 1, {0, 1}, {0}},
837 {-1, 1, {1, 0}, {0}},
838 {-1, 1, {1, 1}, {0}},
839 {-1, 2, {2, 0}, {0}},
840 {-1, 2, {2, 1}, {0}},
841 {-1, 3, {3, 0}, {0}},
842 {-1, 3, {3, 1}, {0}}};
844 std::move(outVertexList),
846 "outMesh"};
847 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
848 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
849 {0, {0}},
850 {0, {0}},
851 {0, {0}},
852 {0, {0}},
853 {0, {0}},
854 {0, {0}},
855 {0, {0}},
856 {0, {0}},
857 {1, {1}},
858 {1, {2}},
859 {1, {3}},
860 {1, {4}},
861 {1, {0}},
862 {1, {0}},
863 {1, {0}},
864 {1, {0}},
865 {2, {0}},
866 {2, {0}},
867 {2, {0}},
868 {2, {0}},
869 {2, {5}},
870 {2, {6}},
871 {2, {0}},
872 {2, {0}},
873 {3, {0}},
874 {3, {0}},
875 {3, {0}},
876 {3, {0}},
877 {3, {0}},
878 {3, {0}},
879 {3, {7}},
880 {3, {8}}};
881
883 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
885 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
887 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
888}
889
891PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
892BOOST_AUTO_TEST_CASE(DistributedConservative2DV3)
893{
894 PRECICE_TEST();
895 Multiquadrics fct(2.0);
896 std::vector<int> globalIndexOffsets = {0, 0, 3, 5};
897
898 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local but rank 0 has no vertices
899 {1, -1, {0, 0}, {1}},
900 {1, -1, {1, 0}, {3}},
901 {1, -1, {1, 1}, {4}},
902 {2, -1, {2, 0}, {5}},
903 {2, -1, {2, 1}, {6}},
904 {3, -1, {3, 0}, {7}},
905 {3, -1, {3, 1}, {8}}}; // Sum of all vertices is 34
907 std::move(inVertexList),
909 "inMesh"};
910 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 owns no vertex
911 {-1, 1, {0, 0}, {0}},
912 {-1, 1, {0, 1}, {0}},
913 {-1, 1, {1, 0}, {0}},
914 {-1, 1, {1, 1}, {0}},
915 {-1, 2, {2, 0}, {0}},
916 {-1, 2, {2, 1}, {0}},
917 {-1, 3, {3, 0}, {0}},
918 {-1, 3, {3, 1}, {0}}};
920 std::move(outVertexList),
922 "outMesh"};
923 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
924 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
925 {0, {0}},
926 {0, {0}},
927 {0, {0}},
928 {0, {0}},
929 {0, {0}},
930 {0, {0}},
931 {0, {0}},
932 {0, {0}},
933 {1, {1}},
934 {1, {0}},
935 {1, {3}},
936 {1, {4}},
937 {1, {0}},
938 {1, {0}},
939 {1, {0}},
940 {1, {0}},
941 {2, {0}},
942 {2, {0}},
943 {2, {0}},
944 {2, {0}},
945 {2, {5}},
946 {2, {6}},
947 {2, {0}},
948 {2, {0}},
949 {3, {0}},
950 {3, {0}},
951 {3, {0}},
952 {3, {0}},
953 {3, {0}},
954 {3, {0}},
955 {3, {7}},
956 {3, {8}}};
957 // Sum of reference is also 34
959 testDistributed(context, mapping_on, in, out, ref, globalIndexOffsets.at(context.rank));
961 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank));
963 testDistributed(context, mapping_off, in, out, ref, globalIndexOffsets.at(context.rank));
964}
965
967PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
968BOOST_AUTO_TEST_CASE(DistributedConservative2DV4,
969 *boost::unit_test::tolerance(1e-6))
970{
971 PRECICE_TEST();
973 std::vector<int> globalIndexOffsets = {0, 2, 4, 6};
974
975 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
976 {0, -1, {0, 0}, {1}},
977 {0, -1, {0, 1}, {2}},
978 {1, -1, {1, 0}, {3}},
979 {1, -1, {1, 1}, {4}},
980 {2, -1, {2, 0}, {5}},
981 {2, -1, {2, 1}, {6}},
982 {3, -1, {3, 0}, {7}},
983 {3, -1, {3, 1}, {8}}}; // Sum is 36
985 std::move(inVertexList),
987 "inMesh"};
988 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed, rank 0 has no vertex at all
989 {-1, 1, {0, 1}, {0}},
990 {-1, 1, {1, 0}, {0}},
991 {-1, 1, {1, 1}, {0}},
992 {-1, 2, {2, 0}, {0}},
993 {-1, 2, {2, 1}, {0}},
994 {-1, 3, {3, 0}, {0}},
995 {-1, 3, {3, 1}, {0}}};
997 std::move(outVertexList),
999 "outMesh"};
1000 ReferenceSpecification ref{// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
1001 // {2, 3, 4, 3, 0, 0, 0, 0} on the second, ...
1002 {0, {0}},
1003 {0, {0}},
1004 {0, {0}},
1005 {0, {0}},
1006 {0, {0}},
1007 {0, {0}},
1008 {0, {0}},
1009 {1, {3.1307987056295525}},
1010 {1, {4.0734031184906971}},
1011 {1, {3.0620533966233006}},
1012 {1, {0}},
1013 {1, {0}},
1014 {1, {0}},
1015 {1, {0}},
1016 {2, {0}},
1017 {2, {0}},
1018 {2, {0}},
1019 {2, {4.4582564956060873}},
1020 {2, {5.8784343572772633}},
1021 {2, {0}},
1022 {2, {0}},
1023 {3, {0}},
1024 {3, {0}},
1025 {3, {0}},
1026 {3, {0}},
1027 {3, {0}},
1028 {3, {7.4683403859032156}},
1029 {3, {7.9287135404698859}}}; // Sum is ~36
1031 testDistributed(context, mapping_sep, in, out, ref, globalIndexOffsets.at(context.rank), true);
1032 // Polynomial == OFF won't reach the desired accuracy
1033}
1034
1036PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1037BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5)
1038{
1039 PRECICE_TEST();
1040 Multiquadrics fct(5.0);
1041 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
1042 {0, -1, {0, 0}, {1}},
1043 {0, -1, {0, 1}, {2}},
1044 {1, -1, {1, 0}, {3}},
1045 {1, -1, {1, 1}, {4}},
1046 {2, -1, {2, 0}, {5}},
1047 {2, -1, {2, 1}, {6}},
1048 {3, -1, {3, 0}, {7}},
1049 {3, -1, {3, 1}, {8}}};
1051 std::move(inVertexList),
1052 meshDims2D,
1053 "inMesh"};
1054 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed and non-contigous
1055 {-1, 0, {0, 0}, {0}},
1056 {-1, 1, {0, 1}, {0}},
1057 {-1, 1, {1, 0}, {0}},
1058 {-1, 0, {1, 1}, {0}},
1059 {-1, 2, {2, 0}, {0}},
1060 {-1, 2, {2, 1}, {0}},
1061 {-1, 3, {3, 0}, {0}},
1062 {-1, 3, {3, 1}, {0}}};
1064 std::move(outVertexList),
1065 meshDims2D,
1066 "outMesh"};
1067 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
1068 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
1069 {0, {1}},
1070 {0, {0}},
1071 {0, {0}},
1072 {0, {4}},
1073 {0, {0}},
1074 {0, {0}},
1075 {0, {0}},
1076 {0, {0}},
1077 {1, {0}},
1078 {1, {2}},
1079 {1, {3}},
1080 {1, {0}},
1081 {1, {0}},
1082 {1, {0}},
1083 {1, {0}},
1084 {1, {0}},
1085 {2, {0}},
1086 {2, {0}},
1087 {2, {0}},
1088 {2, {0}},
1089 {2, {5}},
1090 {2, {6}},
1091 {2, {0}},
1092 {2, {0}},
1093 {3, {0}},
1094 {3, {0}},
1095 {3, {0}},
1096 {3, {0}},
1097 {3, {0}},
1098 {3, {0}},
1099 {3, {7}},
1100 {3, {8}}};
1101
1103 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
1105 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
1107 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
1108}
1109
1111PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1112BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5Vector)
1113{
1114 PRECICE_TEST();
1115 Multiquadrics fct(5.0);
1116
1117 std::vector<VertexSpecification> inVertexList{// Conservative mapping: The inMesh is local
1118 {0, -1, {0, 0}, {1, 4}},
1119 {0, -1, {0, 1}, {2, 5}},
1120 {1, -1, {1, 0}, {3, 6}},
1121 {1, -1, {1, 1}, {4, 7}},
1122 {2, -1, {2, 0}, {5, 8}},
1123 {2, -1, {2, 1}, {6, 9}},
1124 {3, -1, {3, 0}, {7, 10}},
1125 {3, -1, {3, 1}, {8, 11}}};
1127 std::move(inVertexList),
1128 meshDims2D,
1129 "inMesh"};
1130 std::vector<VertexSpecification> outVertexList{// The outMesh is distributed and non-contigous
1131 {-1, 0, {0, 0}, {0, 0}},
1132 {-1, 1, {0, 1}, {0, 0}},
1133 {-1, 1, {1, 0}, {0, 0}},
1134 {-1, 0, {1, 1}, {0, 0}},
1135 {-1, 2, {2, 0}, {0, 0}},
1136 {-1, 2, {2, 1}, {0, 0}},
1137 {-1, 3, {3, 0}, {0, 0}},
1138 {-1, 3, {3, 1}, {0, 0}}};
1140 std::move(outVertexList),
1141 meshDims2D,
1142 "outMesh"};
1143 ReferenceSpecification ref{// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
1144 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
1145 {0, {1, 4}},
1146 {0, {0, 0}},
1147 {0, {0, 0}},
1148 {0, {4, 7}},
1149 {0, {0, 0}},
1150 {0, {0, 0}},
1151 {0, {0, 0}},
1152 {0, {0, 0}},
1153 {1, {0, 0}},
1154 {1, {2, 5}},
1155 {1, {3, 6}},
1156 {1, {0, 0}},
1157 {1, {0, 0}},
1158 {1, {0, 0}},
1159 {1, {0, 0}},
1160 {1, {0, 0}},
1161 {2, {0, 0}},
1162 {2, {0, 0}},
1163 {2, {0, 0}},
1164 {2, {0, 0}},
1165 {2, {5, 8}},
1166 {2, {6, 9}},
1167 {2, {0, 0}},
1168 {2, {0, 0}},
1169 {3, {0, 0}},
1170 {3, {0, 0}},
1171 {3, {0, 0}},
1172 {3, {0, 0}},
1173 {3, {0, 0}},
1174 {3, {0, 0}},
1175 {3, {7, 10}},
1176 {3, {8, 11}}};
1178 testDistributed(context, mapping_on, in, out, ref, context.rank * 2);
1180 testDistributed(context, mapping_sep, in, out, ref, context.rank * 2);
1182 testDistributed(context, mapping_off, in, out, ref, context.rank * 2);
1183}
1184
1185void testTagging(const TestContext &context,
1186 MeshSpecification inMeshSpec,
1187 MeshSpecification outMeshSpec,
1188 MeshSpecification shouldTagFirstRound,
1189 MeshSpecification shouldTagSecondRound,
1190 bool consistent)
1191{
1192 int meshDimension = inMeshSpec.meshDimension;
1193
1194 mesh::PtrMesh inMesh = getDistributedMesh(context, inMeshSpec);
1195 Eigen::VectorXd inValues = getDistributedData(context, inMeshSpec);
1196
1197 mesh::PtrMesh outMesh = getDistributedMesh(context, outMeshSpec);
1198 Eigen::VectorXd outValues = getDistributedData(context, outMeshSpec);
1199
1200 Gaussian fct(4.5); // Support radius approx. 1
1203 inMesh->computeBoundingBox();
1204 outMesh->computeBoundingBox();
1205
1206 mapping.setMeshes(inMesh, outMesh);
1207 mapping.tagMeshFirstRound();
1208
1209 for (const auto &v : inMesh->vertices()) {
1210 auto pos = std::find_if(shouldTagFirstRound.vertices.begin(), shouldTagFirstRound.vertices.end(),
1211 [meshDimension, &v](const VertexSpecification &spec) {
1212 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1213 });
1214 bool found = pos != shouldTagFirstRound.vertices.end();
1215 BOOST_TEST(found >= v.isTagged(),
1216 "FirstRound: Vertex " << v << " is tagged, but should not be.");
1217 BOOST_TEST(found <= v.isTagged(),
1218 "FirstRound: Vertex " << v << " is not tagged, but should be.");
1219 }
1220
1221 mapping.tagMeshSecondRound();
1222
1223 for (const auto &v : inMesh->vertices()) {
1224 auto posFirst = std::find_if(shouldTagFirstRound.vertices.begin(), shouldTagFirstRound.vertices.end(),
1225 [meshDimension, &v](const VertexSpecification &spec) {
1226 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1227 });
1228 bool foundFirst = posFirst != shouldTagFirstRound.vertices.end();
1229 auto posSecond = std::find_if(shouldTagSecondRound.vertices.begin(), shouldTagSecondRound.vertices.end(),
1230 [meshDimension, &v](const VertexSpecification &spec) {
1231 return std::equal(spec.position.data(), spec.position.data() + meshDimension, v.getCoords().data());
1232 });
1233 bool foundSecond = posSecond != shouldTagSecondRound.vertices.end();
1234 BOOST_TEST(foundFirst <= v.isTagged(), "SecondRound: Vertex " << v
1235 << " is not tagged, but should be from the first round.");
1236 BOOST_TEST(foundSecond <= v.isTagged(), "SecondRound: Vertex " << v
1237 << " is not tagged, but should be.");
1238 BOOST_TEST((foundSecond or foundFirst) >= v.isTagged(), "SecondRound: Vertex " << v
1239 << " is tagged, but should not be.");
1240 }
1241}
1242
1243PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm())
1244BOOST_AUTO_TEST_CASE(testTagFirstRound)
1245{
1246 PRECICE_TEST();
1247 // *
1248 // + <-- owned
1249 //* * x * *
1250 // *
1251 // *
1252 std::vector<VertexSpecification> outVertexList{
1253 {0, -1, {0, 0}, {0}}};
1254 MeshSpecification outMeshSpec{
1255 std::move(outVertexList),
1256 meshDims2D,
1257 "outMesh"};
1258 std::vector<VertexSpecification> inVertexList{
1259 {0, -1, {-1, 0}, {1}}, // inside
1260 {0, -1, {-2, 0}, {1}}, // outside
1261 {0, 0, {1, 0}, {1}}, // inside, owner
1262 {0, -1, {2, 0}, {1}}, // outside
1263 {0, -1, {0, -1}, {1}}, // inside
1264 {0, -1, {0, -2}, {1}}, // outside
1265 {0, -1, {0, 1}, {1}}, // inside
1266 {0, -1, {0, 2}, {1}} // outside
1267 };
1268 MeshSpecification inMeshSpec{
1269 std::move(inVertexList),
1270 meshDims2D,
1271 "inMesh"};
1272 std::vector<VertexSpecification> firstRoundVertices = {
1273 {0, -1, {-1, 0}, {1}},
1274 {0, -1, {1, 0}, {1}},
1275 {0, -1, {0, -1}, {1}},
1276 {0, -1, {0, 1}, {1}}};
1277
1278 MeshSpecification shouldTagFirstRound = {
1279 firstRoundVertices,
1280 static_cast<int>(firstRoundVertices.at(0).position.size()),
1281 ""};
1282 std::vector<VertexSpecification> secondRoundVertices = {
1283 {0, -1, {2, 0}, {1}}};
1284 MeshSpecification shouldTagSecondRound = {
1285 secondRoundVertices,
1286 static_cast<int>(secondRoundVertices.at(0).position.size()),
1287 ""};
1288 testTagging(context, inMeshSpec, outMeshSpec, shouldTagFirstRound, shouldTagSecondRound, true);
1289 // For conservative just swap meshes.
1290 testTagging(context, outMeshSpec, inMeshSpec, shouldTagFirstRound, shouldTagSecondRound, false);
1291}
1292
1293BOOST_AUTO_TEST_SUITE_END() // Parallel
1294
1296
1297#undef doLocalCode
1298#define doLocalCode(Type, function, polynomial) \
1299 { \
1300 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap2D(Mapping::CONSISTENT, 2, function, {{false, false, false}}, polynomial); \
1301 perform2DTestConsistentMapping(consistentMap2D); \
1302 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap2DVector(Mapping::CONSISTENT, 2, function, {{false, false, false}}, polynomial); \
1303 perform2DTestConsistentMappingVector(consistentMap2DVector); \
1304 RadialBasisFctMapping<RadialBasisFctSolver<Type>> consistentMap3D(Mapping::CONSISTENT, 3, function, {{false, false, false}}, polynomial); \
1305 perform3DTestConsistentMapping(consistentMap3D); \
1306 RadialBasisFctMapping<RadialBasisFctSolver<Type>> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, function, {{false, false, false}}, polynomial); \
1307 perform2DTestScaledConsistentMapping(scaledConsistentMap2D); \
1308 RadialBasisFctMapping<RadialBasisFctSolver<Type>> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, function, {{false, false, false}}, polynomial); \
1309 perform3DTestScaledConsistentMapping(scaledConsistentMap3D); \
1310 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap2D(Mapping::CONSERVATIVE, 2, function, {{false, false, false}}, polynomial); \
1311 perform2DTestConservativeMapping(conservativeMap2D); \
1312 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap2DVector(Mapping::CONSERVATIVE, 2, function, {{false, false, false}}, polynomial); \
1313 perform2DTestConservativeMappingVector(conservativeMap2DVector); \
1314 RadialBasisFctMapping<RadialBasisFctSolver<Type>> conservativeMap3D(Mapping::CONSERVATIVE, 3, function, {{false, false, false}}, polynomial); \
1315 perform3DTestConservativeMapping(conservativeMap3D); \
1316 }
1317
1318PRECICE_TEST_SETUP(1_rank)
1326
1327PRECICE_TEST_SETUP(1_rank)
1328BOOST_AUTO_TEST_CASE(MapMultiquadrics)
1329{
1330 PRECICE_TEST();
1331 Multiquadrics fct(1e-3);
1334}
1335
1336PRECICE_TEST_SETUP(1_rank)
1337BOOST_AUTO_TEST_CASE(MapInverseMultiquadrics)
1338{
1339 PRECICE_TEST();
1340 InverseMultiquadrics fct(1e-3);
1342}
1343
1344PRECICE_TEST_SETUP(1_rank)
1352
1353PRECICE_TEST_SETUP(1_rank)
1355{
1356 PRECICE_TEST();
1357 Gaussian fct(1.0);
1359}
1360
1361PRECICE_TEST_SETUP(1_rank)
1362BOOST_AUTO_TEST_CASE(MapCompactThinPlateSplinesC2)
1363{
1364 PRECICE_TEST();
1365 double supportRadius = 1.2;
1366 CompactThinPlateSplinesC2 fct(supportRadius);
1368}
1369
1370PRECICE_TEST_SETUP(1_rank)
1371BOOST_AUTO_TEST_CASE(MapCompactPolynomialC0)
1372{
1373 PRECICE_TEST();
1374 double supportRadius = 1.2;
1375 CompactPolynomialC0 fct(supportRadius);
1377}
1378
1379PRECICE_TEST_SETUP(1_rank)
1380BOOST_AUTO_TEST_CASE(MapCompactPolynomialC2)
1381{
1382 PRECICE_TEST();
1383 double supportRadius = 1.2;
1384 CompactPolynomialC2 fct(supportRadius);
1386}
1387
1388PRECICE_TEST_SETUP(1_rank)
1389BOOST_AUTO_TEST_CASE(MapCompactPolynomialC4)
1390{
1391 PRECICE_TEST();
1392 double supportRadius = 1.2;
1393 CompactPolynomialC4 fct(supportRadius);
1395}
1396
1397PRECICE_TEST_SETUP(1_rank)
1398BOOST_AUTO_TEST_CASE(MapCompactPolynomialC6)
1399{
1400 PRECICE_TEST();
1401 double supportRadius = 1.2;
1402 CompactPolynomialC6 fct(supportRadius);
1404}
1405
1406PRECICE_TEST_SETUP(1_rank)
1407BOOST_AUTO_TEST_CASE(MapCompactPolynomialC8)
1408{
1409 PRECICE_TEST();
1410 double supportRadius = 1.2;
1411 CompactPolynomialC8 fct(supportRadius);
1413}
1414#undef doLocalCode
1415
1417{
1418 using Eigen::Vector2d;
1419 int dimensions = 2;
1420
1421 bool xDead = false;
1422 bool yDead = true;
1423 bool zDead = false;
1424
1425 ThinPlateSplines fct;
1427 {{xDead, yDead, zDead}}, polynomial);
1428
1429 // Create mesh to map from
1430 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1431 inMesh->createVertex(Vector2d(0.0, 1.0));
1432 inMesh->createVertex(Vector2d(1.0, 1.0));
1433 inMesh->createVertex(Vector2d(2.0, 1.0));
1434 inMesh->createVertex(Vector2d(3.0, 1.0));
1435 addGlobalIndex(inMesh);
1436 inMesh->setGlobalNumberOfVertices(inMesh->nVertices());
1437
1438 Eigen::VectorXd values(4);
1439 values << 1.0, 2.0, 2.0, 1.0;
1440
1441 // Create mesh to map to
1442 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1443 outMesh->createVertex(Vector2d(0, 1.));
1444 outMesh->createVertex(Vector2d(3, 1.));
1445 outMesh->createVertex(Vector2d(1.3, 1.));
1446 outMesh->createVertex(Vector2d(5, 1.));
1447 addGlobalIndex(outMesh);
1448 outMesh->setGlobalNumberOfVertices(outMesh->nVertices());
1449
1450 // Setup mapping with mapping coordinates and geometry used
1451 mapping.setMeshes(inMesh, outMesh);
1452 BOOST_TEST(mapping.hasComputedMapping() == false);
1453
1454 mapping.computeMapping();
1455 Eigen::VectorXd outValues(4);
1456 time::Sample inSample(1, values);
1457 mapping.map(inSample, outValues);
1458 BOOST_TEST(mapping.hasComputedMapping() == true);
1459 if (constraint == Mapping::CONSISTENT) {
1460 if (polynomial == Polynomial::OFF) {
1461 BOOST_TEST(testing::equals(outValues(2), 2.0522549299731567, 1e-7));
1462 } else if (polynomial == Polynomial::SEPARATE) {
1463 BOOST_TEST(testing::equals(outValues(2), 2.0896514371485777, 1e-7));
1464 } else {
1465 BOOST_TEST(testing::equals(outValues(2), 2.1180354377884774, 1e-7));
1466 }
1467 } else {
1468 if (polynomial == Polynomial::OFF) {
1469 BOOST_TEST(testing::equals(outValues(1), 1.8471144693068295, 1e-7));
1470 } else if (polynomial == Polynomial::SEPARATE) {
1471 BOOST_TEST(testing::equals(outValues(1), 1.8236736422730249, 1e-7));
1472 } else {
1473 BOOST_TEST(testing::equals(outValues(1), 1.7587181970483183, 1e-7));
1474 }
1475 }
1476}
1477
1479{
1480 using Eigen::Vector3d;
1481 int dimensions = 3;
1482
1483 ThinPlateSplines fct;
1484 bool xDead = false;
1485 bool yDead = true;
1486 bool zDead = false;
1488 Mapping mapping(constraint, dimensions, fct, {{xDead, yDead, zDead}}, polynomial);
1489
1490 // Create mesh to map from
1491 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1492 inMesh->createVertex(Vector3d(0.0, 3.0, 0.0));
1493 inMesh->createVertex(Vector3d(1.0, 3.0, 0.0));
1494 inMesh->createVertex(Vector3d(0.0, 3.0, 1.0));
1495 inMesh->createVertex(Vector3d(1.0, 3.0, 1.0));
1496 addGlobalIndex(inMesh);
1497 inMesh->setGlobalNumberOfVertices(inMesh->nVertices());
1498
1499 Eigen::VectorXd values(4);
1500 values << 1.0, 2.0, 3.0, 4.0;
1501
1502 // Create mesh to map to
1503 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1504 outMesh->createVertex(Vector3d(0.0, 2.9, 0.0));
1505 outMesh->createVertex(Vector3d(0.8, 2.9, 0.1));
1506 outMesh->createVertex(Vector3d(0.1, 2.9, 0.9));
1507 outMesh->createVertex(Vector3d(1.1, 2.9, 1.1));
1508 addGlobalIndex(outMesh);
1509 outMesh->setGlobalNumberOfVertices(outMesh->nVertices());
1510
1511 // Setup mapping with mapping coordinates and geometry used
1512 mapping.setMeshes(inMesh, outMesh);
1513 BOOST_TEST(mapping.hasComputedMapping() == false);
1514
1515 mapping.computeMapping();
1516 Eigen::VectorXd outValues(4);
1517 time::Sample inSample(1, values);
1518 mapping.map(inSample, outValues);
1519 BOOST_TEST(mapping.hasComputedMapping() == true);
1520
1521 if (constraint == Mapping::CONSISTENT) {
1522 if (polynomial == Polynomial::OFF) {
1523 const double tolerance = 1e-7;
1524 Eigen::VectorXd expected(4);
1525 expected << 1.0, -0.454450524334, 0.99146426249, 6.98958304876;
1526 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1527 } else {
1528 Eigen::VectorXd expected(4);
1529 expected << 1.0, 2.0, 2.9, 4.3;
1530 BOOST_TEST(testing::equals(outValues, expected));
1531 }
1532 } else {
1533 if (polynomial == Polynomial::OFF) {
1534 const double tolerance = 1e-6;
1535 Eigen::VectorXd expected(4);
1536 expected << 1.17251596926, 4.10368825944, 3.56931954192, 3.40160932341;
1537 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1538 } else if (polynomial == Polynomial::ON) {
1539 const double tolerance = 1e-6;
1540 Eigen::VectorXd expected(4);
1541 expected << 0.856701171969, 2.38947124326, 3.34078733786, 3.41304024691;
1542 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1543 } else {
1544 const double tolerance = 1e-6;
1545 Eigen::VectorXd expected(4);
1546 expected << 0.380480856704, 2.83529451713, 3.73088270249, 3.05334192368;
1547 BOOST_TEST(testing::equals(outValues, expected, tolerance));
1548 }
1549 }
1550}
1551
1552PRECICE_TEST_SETUP(1_rank)
1560
1561PRECICE_TEST_SETUP(1_rank)
1569
1570PRECICE_TEST_SETUP(1_rank)
1578
1579PRECICE_TEST_SETUP(1_rank)
1587
1588BOOST_AUTO_TEST_SUITE_END() // Serial
1589
1591
1592PRECICE_TEST_SETUP(1_rank)
1593BOOST_AUTO_TEST_CASE(inverseTriangularMatrix)
1594{
1595 PRECICE_TEST();
1596 // Define the size of the matrix
1597 const int inSize = 112;
1598
1599 // Generate a random matrix A
1600 Eigen::MatrixXd A = Eigen::MatrixXd::Random(inSize, inSize);
1601
1602 // Create an SPD matrix M
1603 Eigen::MatrixXd M = A * A.transpose();
1604
1605 // Test our custom inversion
1606 Eigen::MatrixXd M_inv_custom = utils::invertLowerTriangularBlockwise(M);
1607
1608 // Second, using Eigen's built-in method
1609 Eigen::MatrixXd M_inv_builtin = Eigen::MatrixXd::Identity(inSize, inSize);
1610 M.triangularView<Eigen::Lower>().solveInPlace(M_inv_builtin);
1611
1612 // Now compare the two inverses and compute the differences
1613 Eigen::MatrixXd diff = M_inv_custom - M_inv_builtin;
1614
1615 double max_abs_diff = diff.maxCoeff();
1616 double min_abs_diff = diff.minCoeff();
1617
1618 // Set a tolerance
1619 double tolerance = 1e-12;
1620
1621 // Use BOOST_CHECK_SMALL to check that the maximum absolute difference is within the tolerance
1622 BOOST_CHECK_SMALL(max_abs_diff, tolerance);
1623 BOOST_CHECK_SMALL(min_abs_diff, tolerance);
1624}
1625
1626BOOST_AUTO_TEST_SUITE_END() // Helper
1627BOOST_AUTO_TEST_SUITE_END() // RadialBasisFunctionMapping
BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithSubsteps)
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
void addGlobalIndex(mesh::PtrMesh &mesh, int offset=0)
std::vector< std::pair< int, std::vector< double > > > ReferenceSpecification
Contains which values are expected on which rank: rank -> vector of data.
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)
#define doLocalCode(Type, function, polynomial)
#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
Utility class for managing MPI operations.
Definition Parallel.hpp:24
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
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:278
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:54
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:64
void computeBoundingBox()
Computes the boundingBox for the vertices.
Definition Mesh.cpp:266
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:104
Rank rank
the rank of the current participant
contains data mapping from points to meshes.
Polynomial
How to handle the polynomial?
provides Mesh, Data and primitives.
std::shared_ptr< Mesh > PtrMesh
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.
std::vector< VertexSpecification > vertices
Holds rank, owner, position and value of a single vertex.