preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PetRadialBasisFctMappingTest.cpp
Go to the documentation of this file.
1#ifndef PRECICE_NO_PETSC
2
3#include <Eigen/Core>
4#include <algorithm>
5#include <map>
6#include <memory>
7#include <ostream>
8#include <set>
9#include <string>
10#include <utility>
11#include <vector>
12#include "logging/Logger.hpp"
13#include "mapping/Mapping.hpp"
17#include "mesh/Data.hpp"
18#include "mesh/Mesh.hpp"
20#include "mesh/Utils.hpp"
21#include "mesh/Vertex.hpp"
22#include "precice/impl/versions.hpp"
24#include "testing/Testing.hpp"
26#include "utils/Petsc.hpp"
27
28using namespace precice;
29using namespace precice::mesh;
30using namespace precice::mapping;
31using namespace precice::testing;
33
34BOOST_AUTO_TEST_SUITE(MappingTests)
35BOOST_AUTO_TEST_SUITE(PetRadialBasisFunctionMapping)
36
37void addGlobalIndex(mesh::PtrMesh &mesh, int offset = 0)
38{
39 for (mesh::Vertex &v : mesh->vertices()) {
40 v.setGlobalIndex(v.getID() + offset);
41 }
42}
43
45{
46 auto inputIntegral = mesh::integrateSurface(inMesh, inData->values());
47 auto outputIntegral = mesh::integrateSurface(outMesh, outData->values());
48
49 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
50 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
51 }
52}
53
55
56
58 int rank;
59 int owner;
62
63 const Eigen::Map<const Eigen::VectorXd> asEigen() const
64 {
65 return Eigen::Map<const Eigen::VectorXd>{position.data(), static_cast<Eigen::Index>(position.size())};
66 }
67};
68
69/*
70MeshSpecification format:
71{ {rank, owner rank, {x, y, z}, {v}}, ... }
72
73also see struct VertexSpecification.
74
75- -1 on rank means all ranks
76- -1 on owner rank means no rank
77- x, y, z is position of vertex, z is optional, 2D mesh will be created then
78- v is the value of the respective vertex. Only 1D supported at this time.
79
80ReferenceSpecification format:
81{ {rank, {v}, ... }
82- -1 on rank means all ranks
83- v is the expected value of n-th vertex on that particular rank
84*/
86
89
90void getDistributedMesh(const TestContext &context,
91 MeshSpecification const &vertices,
92 mesh::PtrMesh &mesh,
93 mesh::PtrData &data,
94 int globalIndexOffset = 0)
95{
96 Eigen::VectorXd d;
97
98 int i = 0;
99 for (auto &vertex : vertices) {
100 if (vertex.rank == context.rank or vertex.rank == -1) {
101 if (vertex.position.size() == 3) // 3-dimensional
102 mesh->createVertex(Eigen::Vector3d(vertex.position.data()));
103 else if (vertex.position.size() == 2) // 2-dimensional
104 mesh->createVertex(Eigen::Vector2d(vertex.position.data()));
105
106 int valueDimension = vertex.value.size();
107
108 if (vertex.owner == context.rank)
109 mesh->vertices().back().setOwner(true);
110 else
111 mesh->vertices().back().setOwner(false);
112
113 d.conservativeResize(i * valueDimension + valueDimension);
114 // Get data in every dimension
115 for (int dim = 0; dim < valueDimension; ++dim) {
116 d(i * valueDimension + dim) = vertex.value.at(dim);
117 }
118 i++;
119 }
120 }
121 addGlobalIndex(mesh, globalIndexOffset);
122 mesh->allocateDataValues();
123 data->values() = d;
124}
125
126void testDistributed(const TestContext &context,
127 Mapping &mapping,
128 MeshSpecification inMeshSpec,
129 MeshSpecification outMeshSpec,
130 ReferenceSpecification referenceSpec,
131 int inGlobalIndexOffset = 0)
132{
133 int meshDimension = inMeshSpec.at(0).position.size();
134 int valueDimension = inMeshSpec.at(0).value.size();
135
136 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", meshDimension, testing::nextMeshID()));
137 mesh::PtrData inData = inMesh->createData("InData", valueDimension, 0_dataID);
138 int inDataID = inData->getID();
139
140 getDistributedMesh(context, inMeshSpec, inMesh, inData, inGlobalIndexOffset);
141
142 mesh::PtrMesh outMesh(new mesh::Mesh("outMesh", meshDimension, testing::nextMeshID()));
143 mesh::PtrData outData = outMesh->createData("OutData", valueDimension, 1_dataID);
144 int outDataID = outData->getID();
145
146 getDistributedMesh(context, outMeshSpec, outMesh, outData);
147
148 mapping.setMeshes(inMesh, outMesh);
149 BOOST_TEST(mapping.hasComputedMapping() == false);
150
151 mapping.computeMapping();
152 BOOST_TEST(mapping.hasComputedMapping() == true);
153 Eigen::VectorXd guess;
154 // Rerun the map and see if a different initial guess influences the result
155 // Run 2 uses last solution as initial guess
156 // Run 3 uses a zero initial guess
157 // Run 4 uses a random initial guess
158 for (auto run : {1, 2, 3, 4}) {
159 if (run == 3) {
160 guess.setZero();
161 }
162 if (run == 4) {
163 guess.setRandom();
164 }
165 mapping.map(inDataID, outDataID, guess);
166
167 int index = 0;
168 for (auto &referenceVertex : referenceSpec) {
169 if (referenceVertex.first == context.rank or referenceVertex.first == -1) {
170 for (int dim = 0; dim < valueDimension; ++dim) {
171 BOOST_TEST_INFO("Index of vertex: " << index << " - Dimension: " << dim);
172 BOOST_TEST(outData->values()(index * valueDimension + dim) == referenceVertex.second.at(dim));
173 }
174 ++index;
175 }
176 }
177 BOOST_TEST(outData->values().size() == index * valueDimension);
178 }
179}
180
182PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
183BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1)
184{
185 PRECICE_TEST();
186 Gaussian fct(5.0);
187 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
188
189 testDistributed(context, mapping,
190 {// Consistent mapping: The inMesh is communicated
191 {-1, 0, {0, 0}, {1}},
192 {-1, 0, {0, 1}, {2}},
193 {-1, 1, {1, 0}, {3}},
194 {-1, 1, {1, 1}, {4}},
195 {-1, 2, {2, 0}, {5}},
196 {-1, 2, {2, 1}, {6}},
197 {-1, 3, {3, 0}, {7}},
198 {-1, 3, {3, 1}, {8}}},
199 {// The outMesh is local, distributed among all ranks
200 {0, -1, {0, 0}, {0}},
201 {0, -1, {0, 1}, {0}},
202 {1, -1, {1, 0}, {0}},
203 {1, -1, {1, 1}, {0}},
204 {2, -1, {2, 0}, {0}},
205 {2, -1, {2, 1}, {0}},
206 {3, -1, {3, 0}, {0}},
207 {3, -1, {3, 1}, {0}}},
208 {// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
209 {0, {1}},
210 {0, {2}},
211 {1, {3}},
212 {1, {4}},
213 {2, {5}},
214 {2, {6}},
215 {3, {7}},
216 {3, {8}}});
217}
218
219PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
220BOOST_AUTO_TEST_CASE(DistributedConsistent2DV1Vector)
221{
222 PRECICE_TEST();
223 Gaussian fct(5.0);
224 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
225
226 testDistributed(context, mapping,
227 {// Consistent mapping: The inMesh is communicated
228 {-1, 0, {0, 0}, {1, 4}},
229 {-1, 0, {0, 1}, {2, 5}},
230 {-1, 1, {1, 0}, {3, 6}},
231 {-1, 1, {1, 1}, {4, 7}},
232 {-1, 2, {2, 0}, {5, 8}},
233 {-1, 2, {2, 1}, {6, 9}},
234 {-1, 3, {3, 0}, {7, 10}},
235 {-1, 3, {3, 1}, {8, 11}}},
236 {// The outMesh is local, distributed among all ranks
237 {0, -1, {0, 0}, {0, 0}},
238 {0, -1, {0, 1}, {0, 0}},
239 {1, -1, {1, 0}, {0, 0}},
240 {1, -1, {1, 1}, {0, 0}},
241 {2, -1, {2, 0}, {0, 0}},
242 {2, -1, {2, 1}, {0, 0}},
243 {3, -1, {3, 0}, {0, 0}},
244 {3, -1, {3, 1}, {0, 0}}},
245 {// Tests for {0, 1} on the first rank, {1, 2} on the second, ...
246 {0, {1, 4}},
247 {0, {2, 5}},
248 {1, {3, 6}},
249 {1, {4, 7}},
250 {2, {5, 8}},
251 {2, {6, 9}},
252 {3, {7, 10}},
253 {3, {8, 11}}});
254}
255
257PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
258BOOST_AUTO_TEST_CASE(DistributedConsistent2DV2)
259{
260 PRECICE_TEST();
261 Gaussian fct(5.0);
262 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
263
264 testDistributed(context, mapping,
265 {// Consistent mapping: The inMesh is communicated, rank 2 owns no vertices
266 {-1, 0, {0, 0}, {1}},
267 {-1, 0, {0, 1}, {2}},
268 {-1, 1, {1, 0}, {3}},
269 {-1, 1, {1, 1}, {4}},
270 {-1, 1, {2, 0}, {5}},
271 {-1, 3, {2, 1}, {6}},
272 {-1, 3, {3, 0}, {7}},
273 {-1, 3, {3, 1}, {8}}},
274 {// The outMesh is local, rank 1 is empty
275 {0, -1, {0, 0}, {0}},
276 {0, -1, {0, 1}, {0}},
277 {0, -1, {1, 0}, {0}},
278 {2, -1, {1, 1}, {0}},
279 {2, -1, {2, 0}, {0}},
280 {2, -1, {2, 1}, {0}},
281 {3, -1, {3, 0}, {0}},
282 {3, -1, {3, 1}, {0}}},
283 {// Tests for {0, 1, 2} on the first rank,
284 // second rank (consistent with the outMesh) is empty, ...
285 {0, {1}},
286 {0, {2}},
287 {0, {3}},
288 {2, {4}},
289 {2, {5}},
290 {2, {6}},
291 {3, {7}},
292 {3, {8}}});
293}
294
296PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
297BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3)
298{
299 PRECICE_TEST();
300 Gaussian fct(5.0);
301 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
302
303 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
304
305 testDistributed(context, mapping,
306 {
307 // Rank 0 has part of the mesh, owns a subpart
308 {0, 0, {0, 0}, {1}},
309 {0, 0, {0, 1}, {2}},
310 {0, 0, {1, 0}, {3}},
311 {0, -1, {1, 1}, {4}},
312 {0, -1, {2, 0}, {5}},
313 {0, -1, {2, 1}, {6}},
314 // Rank 1 has no vertices
315 // Rank 2 has the entire mesh, but owns just 3 and 5.
316 {2, -1, {0, 0}, {1}},
317 {2, -1, {0, 1}, {2}},
318 {2, -1, {1, 0}, {3}},
319 {2, 2, {1, 1}, {4}},
320 {2, -1, {2, 0}, {5}},
321 {2, 2, {2, 1}, {6}},
322 {2, -1, {3, 0}, {7}},
323 {2, -1, {3, 1}, {8}},
324 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
325 {3, 3, {2, 0}, {5}},
326 {3, -1, {2, 1}, {6}},
327 {3, 3, {3, 0}, {7}},
328 {3, 3, {3, 1}, {8}},
329 },
330 {// The outMesh is local, rank 1 is empty
331 {0, -1, {0, 0}, {0}},
332 {0, -1, {0, 1}, {0}},
333 {0, -1, {1, 0}, {0}},
334 {2, -1, {1, 1}, {0}},
335 {2, -1, {2, 0}, {0}},
336 {2, -1, {2, 1}, {0}},
337 {3, -1, {3, 0}, {0}},
338 {3, -1, {3, 1}, {0}}},
339 {// Tests for {0, 1, 2} on the first rank,
340 // second rank (consistent with the outMesh) is empty, ...
341 {0, {1}},
342 {0, {2}},
343 {0, {3}},
344 {2, {4}},
345 {2, {5}},
346 {2, {6}},
347 {3, {7}},
348 {3, {8}}},
349 globalIndexOffsets.at(context.rank));
350}
351
353PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
354BOOST_AUTO_TEST_CASE(DistributedConsistent2DV3Vector)
355{
356 PRECICE_TEST();
357 Gaussian fct(5.0);
358 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
359
360 std::vector<int> globalIndexOffsets = {0, 0, 0, 4};
361
362 testDistributed(context, mapping,
363 {
364 // Rank 0 has part of the mesh, owns a subpart
365 {0, 0, {0, 0}, {1, 4}},
366 {0, 0, {0, 1}, {2, 5}},
367 {0, 0, {1, 0}, {3, 6}},
368 {0, -1, {1, 1}, {4, 7}},
369 {0, -1, {2, 0}, {5, 8}},
370 {0, -1, {2, 1}, {6, 9}},
371 // Rank 1 has no vertices
372 // Rank 2 has the entire mesh, but owns just 3 and 5.
373 {2, -1, {0, 0}, {1, 4}},
374 {2, -1, {0, 1}, {2, 5}},
375 {2, -1, {1, 0}, {3, 6}},
376 {2, 2, {1, 1}, {4, 7}},
377 {2, -1, {2, 0}, {5, 8}},
378 {2, 2, {2, 1}, {6, 9}},
379 {2, -1, {3, 0}, {7, 10}},
380 {2, -1, {3, 1}, {8, 11}},
381 // Rank 3 has the last 4 vertices, owns 4, 6 and 7
382 {3, 3, {2, 0}, {5, 8}},
383 {3, -1, {2, 1}, {6, 9}},
384 {3, 3, {3, 0}, {7, 10}},
385 {3, 3, {3, 1}, {8, 11}},
386 },
387 {// The outMesh is local, rank 1 is empty
388 {0, -1, {0, 0}, {0, 0}},
389 {0, -1, {0, 1}, {0, 0}},
390 {0, -1, {1, 0}, {0, 0}},
391 {2, -1, {1, 1}, {0, 0}},
392 {2, -1, {2, 0}, {0, 0}},
393 {2, -1, {2, 1}, {0, 0}},
394 {3, -1, {3, 0}, {0, 0}},
395 {3, -1, {3, 1}, {0, 0}}},
396 {// Tests for {0, 1, 2} on the first rank,
397 // second rank (consistent with the outMesh) is empty, ...
398 {0, {1, 4}},
399 {0, {2, 5}},
400 {0, {3, 6}},
401 {2, {4, 7}},
402 {2, {5, 8}},
403 {2, {6, 9}},
404 {3, {7, 10}},
405 {3, {8, 11}}},
406 globalIndexOffsets.at(context.rank));
407}
408
410PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
411BOOST_AUTO_TEST_CASE(DistributedConsistent2DV4)
412{
413 PRECICE_TEST();
415 PetRadialBasisFctMapping<ThinPlateSplines> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
416
417 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
418
419 testDistributed(context, mapping,
420 {
421 // Rank 0 has no vertices
422 // Rank 1 has the entire mesh, owns a subpart
423 {1, 1, {0, 0}, {1.1}},
424 {1, 1, {0, 1}, {2.5}},
425 {1, 1, {1, 0}, {3}},
426 {1, 1, {1, 1}, {4}},
427 {1, -1, {2, 0}, {5}},
428 {1, -1, {2, 1}, {6}},
429 {1, -1, {3, 0}, {7}},
430 {1, -1, {3, 1}, {8}},
431 // Rank 2 has the entire mesh, owns a subpart
432 {2, -1, {0, 0}, {1.1}},
433 {2, -1, {0, 1}, {2.5}},
434 {2, -1, {1, 0}, {3}},
435 {2, -1, {1, 1}, {4}},
436 {2, 2, {2, 0}, {5}},
437 {2, 2, {2, 1}, {6}},
438 {2, 2, {3, 0}, {7}},
439 {2, 2, {3, 1}, {8}},
440 // Rank 3 has no vertices
441 },
442 {// The outMesh is local, rank 0 and 3 are empty
443 // not same order as input mesh and vertex (2,0) appears twice
444 {1, -1, {2, 0}, {0}},
445 {1, -1, {1, 0}, {0}},
446 {1, -1, {0, 1}, {0}},
447 {1, -1, {1, 1}, {0}},
448 {1, -1, {0, 0}, {0}},
449 {2, -1, {2, 0}, {0}},
450 {2, -1, {2, 1}, {0}},
451 {2, -1, {3, 0}, {0}},
452 {2, -1, {3, 1}, {0}}},
453 {{1, {5}},
454 {1, {3}},
455 {1, {2.5}},
456 {1, {4}},
457 {1, {1.1}},
458 {2, {5}},
459 {2, {6}},
460 {2, {7}},
461 {2, {8}}},
462 globalIndexOffsets.at(context.rank));
463}
464
465// same as 2DV4, but all ranks have vertices
466PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
467BOOST_AUTO_TEST_CASE(DistributedConsistent2DV5)
468{
469 PRECICE_TEST();
471 PetRadialBasisFctMapping<ThinPlateSplines> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
472
473 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
474
475 testDistributed(context, mapping,
476 {
477 // Every rank has the entire mesh and owns a subpart
478 {0, 0, {0, 0}, {1.1}},
479 {0, 0, {0, 1}, {2.5}},
480 {0, -1, {1, 0}, {3}},
481 {0, -1, {1, 1}, {4}},
482 {0, -1, {2, 0}, {5}},
483 {0, -1, {2, 1}, {6}},
484 {0, -1, {3, 0}, {7}},
485 {0, -1, {3, 1}, {8}},
486 {1, -1, {0, 0}, {1.1}},
487 {1, -1, {0, 1}, {2.5}},
488 {1, 1, {1, 0}, {3}},
489 {1, 1, {1, 1}, {4}},
490 {1, -1, {2, 0}, {5}},
491 {1, -1, {2, 1}, {6}},
492 {1, -1, {3, 0}, {7}},
493 {1, -1, {3, 1}, {8}},
494 {2, -1, {0, 0}, {1.1}},
495 {2, -1, {0, 1}, {2.5}},
496 {2, -1, {1, 0}, {3}},
497 {2, -1, {1, 1}, {4}},
498 {2, 2, {2, 0}, {5}},
499 {2, 2, {2, 1}, {6}},
500 {2, -1, {3, 0}, {7}},
501 {2, -1, {3, 1}, {8}},
502 {3, -1, {0, 0}, {1.1}},
503 {3, -1, {0, 1}, {2.5}},
504 {3, -1, {1, 0}, {3}},
505 {3, -1, {1, 1}, {4}},
506 {3, -1, {2, 0}, {5}},
507 {3, -1, {2, 1}, {6}},
508 {3, 3, {3, 0}, {7}},
509 {3, 3, {3, 1}, {8}},
510 },
511 {// The outMesh is local, rank 0 and 3 are empty
512 // not same order as input mesh and vertex (2,0) appears twice
513 {1, -1, {2, 0}, {0}},
514 {1, -1, {1, 0}, {0}},
515 {1, -1, {0, 1}, {0}},
516 {1, -1, {1, 1}, {0}},
517 {1, -1, {0, 0}, {0}},
518 {2, -1, {2, 0}, {0}},
519 {2, -1, {2, 1}, {0}},
520 {2, -1, {3, 0}, {0}},
521 {2, -1, {3, 1}, {0}}},
522 {{1, {5}},
523 {1, {3}},
524 {1, {2.5}},
525 {1, {4}},
526 {1, {1.1}},
527 {2, {5}},
528 {2, {6}},
529 {2, {7}},
530 {2, {8}}},
531 globalIndexOffsets.at(context.rank));
532}
533
535PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
536BOOST_AUTO_TEST_CASE(DistributedConsistent2DV6,
537 *boost::unit_test::tolerance(1e-7))
538{
539 PRECICE_TEST();
541 PetRadialBasisFctMapping<ThinPlateSplines> mapping(Mapping::CONSISTENT, 2, fct, {{false, false, false}});
542
543 std::vector<int> globalIndexOffsets = {0, 0, 0, 0};
544
545 testDistributed(context, mapping,
546 {
547 // Rank 0 has no vertices
548 // Rank 1 has the entire mesh, owns a subpart
549 {1, 1, {0, 0}, {1}},
550 {1, 1, {0, 1}, {2}},
551 {1, 1, {1, 0}, {3}},
552 {1, 1, {1, 1}, {4}},
553 {1, -1, {2, 0}, {5}},
554 {1, -1, {2, 1}, {6}},
555 {1, -1, {3, 0}, {7}},
556 {1, -1, {3, 1}, {8}},
557 // Rank 2 has the entire mesh, owns a subpart
558 {2, -1, {0, 0}, {1}},
559 {2, -1, {0, 1}, {2}},
560 {2, -1, {1, 0}, {3}},
561 {2, -1, {1, 1}, {4}},
562 {2, 2, {2, 0}, {5}},
563 {2, 2, {2, 1}, {6}},
564 {2, 2, {3, 0}, {7}},
565 {2, 2, {3, 1}, {8}},
566 // Rank 3 has no vertices
567 },
568 {// The outMesh is local, rank 0 and 3 are empty
569 // not same order as input mesh and vertex (2,0) appears twice
570 {1, -1, {2, 0}, {0}},
571 {1, -1, {1, 0}, {0}},
572 {1, -1, {0, 1}, {0}},
573 {1, -1, {1, 1}, {0}},
574 {1, -1, {0, 0}, {0}},
575 {2, -1, {2, 0}, {0}},
576 {2, -1, {2, 1}, {0}},
577 {2, -1, {3, 0}, {0}},
578 {2, -1, {3, 1}, {0}}},
579 {{1, {5}},
580 {1, {3}},
581 {1, {2}},
582 {1, {4}},
583 {1, {1}},
584 {2, {5}},
585 {2, {6}},
586 {2, {7}},
587 {2, {8}}},
588 globalIndexOffsets.at(context.rank));
589}
590
592PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
593BOOST_AUTO_TEST_CASE(DistributedConservative2DV1)
594{
595 PRECICE_TEST();
596 Gaussian fct(5.0);
597 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
598
599 testDistributed(context, mapping,
600 {// Conservative mapping: The inMesh is local
601 {0, -1, {0, 0}, {1}},
602 {0, -1, {0, 1}, {2}},
603 {1, -1, {1, 0}, {3}},
604 {1, -1, {1, 1}, {4}},
605 {2, -1, {2, 0}, {5}},
606 {2, -1, {2, 1}, {6}},
607 {3, -1, {3, 0}, {7}},
608 {3, -1, {3, 1}, {8}}},
609 {// The outMesh is distributed
610 {-1, 0, {0, 0}, {0}},
611 {-1, 0, {0, 1}, {0}},
612 {-1, 1, {1, 0}, {0}},
613 {-1, 1, {1, 1}, {0}},
614 {-1, 2, {2, 0}, {0}},
615 {-1, 2, {2, 1}, {0}},
616 {-1, 3, {3, 0}, {0}},
617 {-1, 3, {3, 1}, {0}}},
618 {// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
619 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
620 {0, {1}},
621 {0, {2}},
622 {0, {0}},
623 {0, {0}},
624 {0, {0}},
625 {0, {0}},
626 {0, {0}},
627 {0, {0}},
628 {1, {0}},
629 {1, {0}},
630 {1, {3}},
631 {1, {4}},
632 {1, {0}},
633 {1, {0}},
634 {1, {0}},
635 {1, {0}},
636 {2, {0}},
637 {2, {0}},
638 {2, {0}},
639 {2, {0}},
640 {2, {5}},
641 {2, {6}},
642 {2, {0}},
643 {2, {0}},
644 {3, {0}},
645 {3, {0}},
646 {3, {0}},
647 {3, {0}},
648 {3, {0}},
649 {3, {0}},
650 {3, {7}},
651 {3, {8}}},
652 context.rank * 2);
653}
654
656PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
657BOOST_AUTO_TEST_CASE(DistributedConservative2DV1Vector)
658{
659 PRECICE_TEST();
660 Gaussian fct(5.0);
661 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
662
663 testDistributed(context, mapping,
664 {// Conservative mapping: The inMesh is local
665 {0, -1, {0, 0}, {1, 4}},
666 {0, -1, {0, 1}, {2, 5}},
667 {1, -1, {1, 0}, {3, 6}},
668 {1, -1, {1, 1}, {4, 7}},
669 {2, -1, {2, 0}, {5, 8}},
670 {2, -1, {2, 1}, {6, 9}},
671 {3, -1, {3, 0}, {7, 10}},
672 {3, -1, {3, 1}, {8, 11}}},
673 {// The outMesh is distributed
674 {-1, 0, {0, 0}, {0, 0}},
675 {-1, 0, {0, 1}, {0, 0}},
676 {-1, 1, {1, 0}, {0, 0}},
677 {-1, 1, {1, 1}, {0, 0}},
678 {-1, 2, {2, 0}, {0, 0}},
679 {-1, 2, {2, 1}, {0, 0}},
680 {-1, 3, {3, 0}, {0, 0}},
681 {-1, 3, {3, 1}, {0, 0}}},
682 {// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
683 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
684 {0, {1, 4}},
685 {0, {2, 5}},
686 {0, {0, 0}},
687 {0, {0, 0}},
688 {0, {0, 0}},
689 {0, {0, 0}},
690 {0, {0, 0}},
691 {0, {0, 0}},
692 {1, {0, 0}},
693 {1, {0, 0}},
694 {1, {3, 6}},
695 {1, {4, 7}},
696 {1, {0, 0}},
697 {1, {0, 0}},
698 {1, {0, 0}},
699 {1, {0, 0}},
700 {2, {0, 0}},
701 {2, {0, 0}},
702 {2, {0, 0}},
703 {2, {0, 0}},
704 {2, {5, 8}},
705 {2, {6, 9}},
706 {2, {0, 0}},
707 {2, {0, 0}},
708 {3, {0, 0}},
709 {3, {0, 0}},
710 {3, {0, 0}},
711 {3, {0, 0}},
712 {3, {0, 0}},
713 {3, {0, 0}},
714 {3, {7, 10}},
715 {3, {8, 11}}},
716 context.rank * 2);
717}
718
720PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
721BOOST_AUTO_TEST_CASE(DistributedConservative2DV2)
722{
723 PRECICE_TEST();
724 Gaussian fct(5.0);
725 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
726
727 std::vector<int> globalIndexOffsets = {0, 0, 4, 6};
728
729 testDistributed(context, mapping,
730 {// Conservative mapping: The inMesh is local but rank 0 has no vertices
731 {1, -1, {0, 0}, {1}},
732 {1, -1, {0, 1}, {2}},
733 {1, -1, {1, 0}, {3}},
734 {1, -1, {1, 1}, {4}},
735 {2, -1, {2, 0}, {5}},
736 {2, -1, {2, 1}, {6}},
737 {3, -1, {3, 0}, {7}},
738 {3, -1, {3, 1}, {8}}},
739 {// The outMesh is distributed, rank 0 owns no vertex
740 {-1, 1, {0, 0}, {0}},
741 {-1, 1, {0, 1}, {0}},
742 {-1, 1, {1, 0}, {0}},
743 {-1, 1, {1, 1}, {0}},
744 {-1, 2, {2, 0}, {0}},
745 {-1, 2, {2, 1}, {0}},
746 {-1, 3, {3, 0}, {0}},
747 {-1, 3, {3, 1}, {0}}},
748 {// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
749 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
750 {0, {0}},
751 {0, {0}},
752 {0, {0}},
753 {0, {0}},
754 {0, {0}},
755 {0, {0}},
756 {0, {0}},
757 {0, {0}},
758 {1, {1}},
759 {1, {2}},
760 {1, {3}},
761 {1, {4}},
762 {1, {0}},
763 {1, {0}},
764 {1, {0}},
765 {1, {0}},
766 {2, {0}},
767 {2, {0}},
768 {2, {0}},
769 {2, {0}},
770 {2, {5}},
771 {2, {6}},
772 {2, {0}},
773 {2, {0}},
774 {3, {0}},
775 {3, {0}},
776 {3, {0}},
777 {3, {0}},
778 {3, {0}},
779 {3, {0}},
780 {3, {7}},
781 {3, {8}}},
782 globalIndexOffsets.at(context.rank));
783}
784
786PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
787BOOST_AUTO_TEST_CASE(DistributedConservative2DV3)
788{
789 PRECICE_TEST();
790 Gaussian fct(2.0);
791 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
792
793 std::vector<int> globalIndexOffsets = {0, 0, 3, 5};
794
795 testDistributed(context, mapping,
796 {// Conservative mapping: The inMesh is local but rank 0 has no vertices
797 {1, -1, {0, 0}, {1}},
798 {1, -1, {1, 0}, {3}},
799 {1, -1, {1, 1}, {4}},
800 {2, -1, {2, 0}, {5}},
801 {2, -1, {2, 1}, {6}},
802 {3, -1, {3, 0}, {7}},
803 {3, -1, {3, 1}, {8}}}, // Sum of all vertices is 34
804 { // The outMesh is distributed, rank 0 owns no vertex
805 {-1, 1, {0, 0}, {0}},
806 {-1, 1, {0, 1}, {0}},
807 {-1, 1, {1, 0}, {0}},
808 {-1, 1, {1, 1}, {0}},
809 {-1, 2, {2, 0}, {0}},
810 {-1, 2, {2, 1}, {0}},
811 {-1, 3, {3, 0}, {0}},
812 {-1, 3, {3, 1}, {0}}},
813 {// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
814 // {1, 2, 2, 3, 0, 0, 0, 0} on the second, ...
815 {0, {0}},
816 {0, {0}},
817 {0, {0}},
818 {0, {0}},
819 {0, {0}},
820 {0, {0}},
821 {0, {0}},
822 {0, {0}},
823 {1, {1}},
824 {1, {0}},
825 {1, {3}},
826 {1, {4}},
827 {1, {0}},
828 {1, {0}},
829 {1, {0}},
830 {1, {0}},
831 {2, {0}},
832 {2, {0}},
833 {2, {0}},
834 {2, {0}},
835 {2, {5}},
836 {2, {6}},
837 {2, {0}},
838 {2, {0}},
839 {3, {0}},
840 {3, {0}},
841 {3, {0}},
842 {3, {0}},
843 {3, {0}},
844 {3, {0}},
845 {3, {7}},
846 {3, {8}}}, // Sum of reference is also 34
847 globalIndexOffsets.at(context.rank));
848}
849
851PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
852BOOST_AUTO_TEST_CASE(DistributedConservative2DV4,
853 *boost::unit_test::tolerance(1e-6))
854{
855 PRECICE_TEST();
856 Gaussian fct(4.0);
857 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
858
859 std::vector<int> globalIndexOffsets = {0, 2, 4, 6};
860
861 testDistributed(context, mapping,
862 {// Conservative mapping: The inMesh is local
863 {0, -1, {0, 0}, {1}},
864 {0, -1, {0, 1}, {2}},
865 {1, -1, {1, 0}, {3}},
866 {1, -1, {1, 1}, {4}},
867 {2, -1, {2, 0}, {5}},
868 {2, -1, {2, 1}, {6}},
869 {3, -1, {3, 0}, {7}},
870 {3, -1, {3, 1}, {8}}}, // Sum is 36
871 { // The outMesh is distributed, rank 0 has no vertex at all
872 {-1, 1, {0, 1}, {0}},
873 {-1, 1, {1, 0}, {0}},
874 {-1, 1, {1, 1}, {0}},
875 {-1, 2, {2, 0}, {0}},
876 {-1, 2, {2, 1}, {0}},
877 {-1, 3, {3, 0}, {0}},
878 {-1, 3, {3, 1}, {0}}},
879 {// Tests for {0, 0, 0, 0, 0, 0, 0, 0} on the first rank,
880 // {2, 3, 4, 3, 0, 0, 0, 0} on the second, ...
881 {0, {0}},
882 {0, {0}},
883 {0, {0}},
884 {0, {0}},
885 {0, {0}},
886 {0, {0}},
887 {0, {0}},
888 {1, {2.4285714526861519}},
889 {1, {3.61905}},
890 {1, {4.14286}},
891 {1, {0}},
892 {1, {0}},
893 {1, {0}},
894 {1, {0}},
895 {2, {0}},
896 {2, {0}},
897 {2, {0}},
898 {2, {5.333333295}},
899 {2, {5.85714}},
900 {2, {0}},
901 {2, {0}},
902 {3, {0}},
903 {3, {0}},
904 {3, {0}},
905 {3, {0}},
906 {3, {0}},
907 {3, {7.047619}},
908 {3, {7.571428}}}, // Sum is ~36
909 globalIndexOffsets.at(context.rank));
910}
911
913PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
914BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5)
915{
916 PRECICE_TEST();
917 Gaussian fct(5.0);
918 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
919
920 testDistributed(context, mapping,
921 {// Conservative mapping: The inMesh is local
922 {0, -1, {0, 0}, {1}},
923 {0, -1, {0, 1}, {2}},
924 {1, -1, {1, 0}, {3}},
925 {1, -1, {1, 1}, {4}},
926 {2, -1, {2, 0}, {5}},
927 {2, -1, {2, 1}, {6}},
928 {3, -1, {3, 0}, {7}},
929 {3, -1, {3, 1}, {8}}},
930 {// The outMesh is distributed and non-contigous
931 {-1, 0, {0, 0}, {0}},
932 {-1, 1, {0, 1}, {0}},
933 {-1, 1, {1, 0}, {0}},
934 {-1, 0, {1, 1}, {0}},
935 {-1, 2, {2, 0}, {0}},
936 {-1, 2, {2, 1}, {0}},
937 {-1, 3, {3, 0}, {0}},
938 {-1, 3, {3, 1}, {0}}},
939 {// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
940 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
941 {0, {1}},
942 {0, {0}},
943 {0, {0}},
944 {0, {4}},
945 {0, {0}},
946 {0, {0}},
947 {0, {0}},
948 {0, {0}},
949 {1, {0}},
950 {1, {2}},
951 {1, {3}},
952 {1, {0}},
953 {1, {0}},
954 {1, {0}},
955 {1, {0}},
956 {1, {0}},
957 {2, {0}},
958 {2, {0}},
959 {2, {0}},
960 {2, {0}},
961 {2, {5}},
962 {2, {6}},
963 {2, {0}},
964 {2, {0}},
965 {3, {0}},
966 {3, {0}},
967 {3, {0}},
968 {3, {0}},
969 {3, {0}},
970 {3, {0}},
971 {3, {7}},
972 {3, {8}}},
973 context.rank * 2);
974}
975
977PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
978BOOST_AUTO_TEST_CASE(testDistributedConservative2DV5Vector)
979{
980 PRECICE_TEST();
981 Gaussian fct(5.0);
982 PetRadialBasisFctMapping<Gaussian> mapping(Mapping::CONSERVATIVE, 2, fct, {{false, false, false}});
983
984 testDistributed(context, mapping,
985 {// Conservative mapping: The inMesh is local
986 {0, -1, {0, 0}, {1, 4}},
987 {0, -1, {0, 1}, {2, 5}},
988 {1, -1, {1, 0}, {3, 6}},
989 {1, -1, {1, 1}, {4, 7}},
990 {2, -1, {2, 0}, {5, 8}},
991 {2, -1, {2, 1}, {6, 9}},
992 {3, -1, {3, 0}, {7, 10}},
993 {3, -1, {3, 1}, {8, 11}}},
994 {// The outMesh is distributed and non-contigous
995 {-1, 0, {0, 0}, {0, 0}},
996 {-1, 1, {0, 1}, {0, 0}},
997 {-1, 1, {1, 0}, {0, 0}},
998 {-1, 0, {1, 1}, {0, 0}},
999 {-1, 2, {2, 0}, {0, 0}},
1000 {-1, 2, {2, 1}, {0, 0}},
1001 {-1, 3, {3, 0}, {0, 0}},
1002 {-1, 3, {3, 1}, {0, 0}}},
1003 {// Tests for {0, 1, 0, 0, 0, 0, 0, 0} on the first rank,
1004 // {0, 0, 2, 3, 0, 0, 0, 0} on the second, ...
1005 {0, {1, 4}},
1006 {0, {0, 0}},
1007 {0, {0, 0}},
1008 {0, {4, 7}},
1009 {0, {0, 0}},
1010 {0, {0, 0}},
1011 {0, {0, 0}},
1012 {0, {0, 0}},
1013 {1, {0, 0}},
1014 {1, {2, 5}},
1015 {1, {3, 6}},
1016 {1, {0, 0}},
1017 {1, {0, 0}},
1018 {1, {0, 0}},
1019 {1, {0, 0}},
1020 {1, {0, 0}},
1021 {2, {0, 0}},
1022 {2, {0, 0}},
1023 {2, {0, 0}},
1024 {2, {0, 0}},
1025 {2, {5, 8}},
1026 {2, {6, 9}},
1027 {2, {0, 0}},
1028 {2, {0, 0}},
1029 {3, {0, 0}},
1030 {3, {0, 0}},
1031 {3, {0, 0}},
1032 {3, {0, 0}},
1033 {3, {0, 0}},
1034 {3, {0, 0}},
1035 {3, {7, 10}},
1036 {3, {8, 11}}},
1037 context.rank * 2);
1038}
1039
1040void testTagging(const TestContext &context,
1041 MeshSpecification inMeshSpec,
1042 MeshSpecification outMeshSpec,
1043 MeshSpecification shouldTagFirstRound,
1044 MeshSpecification shouldTagSecondRound,
1045 bool consistent)
1046{
1047 int meshDimension = inMeshSpec.at(0).position.size();
1048 int valueDimension = inMeshSpec.at(0).value.size();
1049
1050 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", meshDimension, testing::nextMeshID()));
1051 mesh::PtrData inData = inMesh->createData("InData", valueDimension, 0_dataID);
1052 getDistributedMesh(context, inMeshSpec, inMesh, inData);
1053
1054 mesh::PtrMesh outMesh(new mesh::Mesh("outMesh", meshDimension, testing::nextMeshID()));
1055 mesh::PtrData outData = outMesh->createData("OutData", valueDimension, 1_dataID);
1056 getDistributedMesh(context, outMeshSpec, outMesh, outData);
1057 BOOST_TEST_MESSAGE("Mesh sizes in: " << inMesh->nVertices() << " out: " << outMesh->nVertices());
1058
1059 Gaussian fct(4.5); // Support radius approx. 1
1060 BOOST_TEST_MESSAGE("Basis function has support radius " << fct.getSupportRadius());
1061 BOOST_TEST(fct.getSupportRadius() > 1.0);
1062 BOOST_TEST(fct.hasCompactSupport());
1063
1065 PetRadialBasisFctMapping<Gaussian> mapping(constr, 2, fct, {{false, false, false}});
1066 inMesh->computeBoundingBox();
1067 outMesh->computeBoundingBox();
1068
1069 mapping.setMeshes(inMesh, outMesh);
1070 mapping.tagMeshFirstRound();
1071
1072 const auto &taggedMesh = consistent ? inMesh : outMesh;
1073
1074 // Expected set of tagged elements for first round
1076 for (const auto &vspec : shouldTagFirstRound) {
1077 expectedFirst.emplace(vspec.asEigen());
1078 }
1079
1080 for (const auto &v : taggedMesh->vertices()) {
1081 bool found = expectedFirst.count(v.getCoords()) != 0;
1082 BOOST_TEST((!found || v.isTagged()),
1083 "FirstRound: Vertex " << v << " is tagged, but should not be.");
1084 BOOST_TEST((found || !v.isTagged()),
1085 "FirstRound: Vertex " << v << " is not tagged, but should be.");
1086 }
1087
1088 // Expected set of tagged elements for second round
1090 expectedFirst.begin(), expectedFirst.end());
1091 for (const auto &vspec : shouldTagSecondRound) {
1092 expectedSecond.emplace(vspec.asEigen());
1093 }
1094
1095 mapping.tagMeshSecondRound();
1096
1097 for (const auto &v : taggedMesh->vertices()) {
1098 bool found = expectedSecond.count(v.getCoords()) != 0;
1099 BOOST_TEST((!found || v.isTagged()),
1100 "SecondRound: Vertex " << v << " is tagged, but should not be.");
1101 BOOST_TEST((found || !v.isTagged()),
1102 "SecondRound: Vertex " << v << " is not tagged, but should be.");
1103 }
1104}
1105
1106PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
1107BOOST_AUTO_TEST_CASE(TaggingConsistent)
1108{
1109 PRECICE_TEST();
1110 // *
1111 // + <-- owned
1112 //* * x * *
1113 // *
1114 // *
1115 MeshSpecification outMeshSpec = {
1116 {0, -1, {0, 0}, {0}}};
1117 MeshSpecification inMeshSpec = {
1118 {0, -1, {-1, 0}, {1}}, // inside
1119 {0, -1, {-2, 0}, {1}}, // outside
1120 {0, 0, {1, 0}, {1}}, // inside, owner
1121 {0, -1, {2, 0}, {1}}, // outside
1122 {0, -1, {0, -1}, {1}}, // inside
1123 {0, -1, {0, -2}, {1}}, // outside
1124 {0, -1, {0, 1}, {1}}, // inside
1125 {0, -1, {0, 2}, {1}} // outside
1126 };
1127 MeshSpecification shouldTagFirstRound = {
1128 {0, -1, {-1, 0}, {1}},
1129 {0, -1, {1, 0}, {1}},
1130 {0, -1, {0, -1}, {1}},
1131 {0, -1, {0, 1}, {1}}};
1132 MeshSpecification shouldTagSecondRound = {
1133 {0, -1, {2, 0}, {1}}};
1134 testTagging(context, inMeshSpec, outMeshSpec, shouldTagFirstRound, shouldTagSecondRound, true);
1135}
1136
1137PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::PETSc)
1138BOOST_AUTO_TEST_CASE(TaggingConservative)
1139{
1140 PRECICE_TEST();
1141 // *
1142 // + <-- owned
1143 //* * x * *
1144 // *
1145 // *
1146 MeshSpecification outMeshSpec = {
1147 {0, -1, {0, 0}, {0}}};
1148 MeshSpecification inMeshSpec = {
1149 {0, -1, {-1, 0}, {1}}, // inside
1150 {0, -1, {-2, 0}, {1}}, // outside
1151 {0, 0, {1, 0}, {1}}, // inside, owner
1152 {0, -1, {2, 0}, {1}}, // outside
1153 {0, -1, {0, -1}, {1}}, // inside
1154 {0, -1, {0, -2}, {1}}, // outside
1155 {0, -1, {0, 1}, {1}}, // inside
1156 {0, -1, {0, 2}, {1}} // outside
1157 };
1158 MeshSpecification shouldTagFirstRound = {
1159 {0, -1, {0, 0}, {1}}};
1160 MeshSpecification shouldTagSecondRound = {
1161 {0, -1, {0, 0}, {1}}};
1162 testTagging(context, inMeshSpec, outMeshSpec, shouldTagFirstRound, shouldTagSecondRound, false);
1163}
1164
1165BOOST_AUTO_TEST_SUITE_END() // Parallel
1166
1168
1170{
1171 int dimensions = 2;
1172 using Eigen::Vector2d;
1173
1174 // Create mesh to map from
1175 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1176 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1177 int inDataID = inData->getID();
1178 inMesh->createVertex(Vector2d(0.0, 0.0));
1179 inMesh->createVertex(Vector2d(1.0, 0.0));
1180 inMesh->createVertex(Vector2d(1.0, 1.0));
1181 inMesh->createVertex(Vector2d(0.0, 1.0));
1182 inMesh->allocateDataValues();
1183 addGlobalIndex(inMesh);
1184
1185 auto &values = inData->values();
1186 values << 1.0, 2.0, 2.0, 1.0;
1187
1188 // Create mesh to map to
1189 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1190 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1191 int outDataID = outData->getID();
1192 mesh::Vertex &vertex = outMesh->createVertex(Vector2d(0, 0));
1193 outMesh->allocateDataValues();
1194 addGlobalIndex(outMesh);
1195
1196 // Setup mapping with mapping coordinates and geometry used
1197 mapping.setMeshes(inMesh, outMesh);
1198 BOOST_TEST(mapping.hasComputedMapping() == false);
1199
1200 vertex.setCoords(Vector2d(0.0, 0.0));
1201 mapping.computeMapping();
1202 Eigen::VectorXd guess;
1203 mapping.map(inDataID, outDataID, guess);
1204 double value = outData->values()(0);
1205 BOOST_TEST(mapping.hasComputedMapping() == true);
1206 BOOST_TEST(value == 1.0);
1207 BOOST_TEST(guess.size() > 0);
1208 guess.resize(0);
1209
1210 vertex.setCoords(Vector2d(0.0, 0.5));
1211 mapping.computeMapping();
1212 mapping.map(inDataID, outDataID, guess);
1213 value = outData->values()(0);
1214 BOOST_TEST(mapping.hasComputedMapping() == true);
1215 BOOST_TEST(value == 1.0);
1216 BOOST_TEST(guess.size() > 0);
1217 guess.resize(0);
1218
1219 vertex.setCoords(Vector2d(0.0, 1.0));
1220 mapping.computeMapping();
1221 mapping.map(inDataID, outDataID, guess);
1222 value = outData->values()(0);
1223 BOOST_TEST(mapping.hasComputedMapping() == true);
1224 BOOST_TEST(value == 1.0);
1225 BOOST_TEST(guess.size() > 0);
1226 guess.resize(0);
1227
1228 vertex.setCoords(Vector2d(1.0, 0.0));
1229 mapping.computeMapping();
1230 mapping.map(inDataID, outDataID, guess);
1231 value = outData->values()(0);
1232 BOOST_TEST(mapping.hasComputedMapping() == true);
1233 BOOST_TEST(value == 2.0);
1234 BOOST_TEST(guess.size() > 0);
1235 guess.resize(0);
1236
1237 vertex.setCoords(Vector2d(1.0, 0.5));
1238 mapping.computeMapping();
1239 mapping.map(inDataID, outDataID, guess);
1240 value = outData->values()(0);
1241 BOOST_TEST(mapping.hasComputedMapping() == true);
1242 BOOST_TEST(value == 2.0);
1243 BOOST_TEST(guess.size() > 0);
1244 guess.resize(0);
1245
1246 vertex.setCoords(Vector2d(1.0, 1.0));
1247 mapping.computeMapping();
1248 mapping.map(inDataID, outDataID, guess);
1249 value = outData->values()(0);
1250 BOOST_TEST(mapping.hasComputedMapping() == true);
1251 BOOST_TEST(value == 2.0);
1252 BOOST_TEST(guess.size() > 0);
1253 guess.resize(0);
1254
1255 vertex.setCoords(Vector2d(0.5, 0.0));
1256 mapping.computeMapping();
1257 mapping.map(inDataID, outDataID, guess);
1258 value = outData->values()(0);
1259 BOOST_TEST(mapping.hasComputedMapping() == true);
1260 BOOST_TEST(value == 1.5);
1261 BOOST_TEST(guess.size() > 0);
1262 guess.resize(0);
1263
1264 vertex.setCoords(Vector2d(0.5, 0.5));
1265 mapping.computeMapping();
1266 mapping.map(inDataID, outDataID, guess);
1267 value = outData->values()(0);
1268 BOOST_TEST(mapping.hasComputedMapping() == true);
1269 BOOST_TEST(value == 1.5);
1270 BOOST_TEST(guess.size() > 0);
1271 guess.resize(0);
1272
1273 vertex.setCoords(Vector2d(0.5, 1.0));
1274 mapping.computeMapping();
1275 mapping.map(inDataID, outDataID, guess);
1276 value = outData->values()(0);
1277 BOOST_TEST(mapping.hasComputedMapping() == true);
1278 BOOST_TEST(value == 1.5);
1279 BOOST_TEST(guess.size() > 0);
1280}
1281
1283{
1284 int dimensions = 2;
1285 using Eigen::Vector2d;
1286
1287 // Create mesh to map from
1288 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1289 mesh::PtrData inData = inMesh->createData("InData", 2, 0_dataID);
1290 int inDataID = inData->getID();
1291 inMesh->createVertex(Vector2d(0.0, 0.0));
1292 inMesh->createVertex(Vector2d(1.0, 0.0));
1293 inMesh->createVertex(Vector2d(1.0, 1.0));
1294 inMesh->createVertex(Vector2d(0.0, 1.0));
1295 inMesh->allocateDataValues();
1296 addGlobalIndex(inMesh);
1297
1298 auto &values = inData->values();
1299 values << 1.0, 4.0, 2.0, 5.0, 2.0, 5.0, 1.0, 4.0;
1300
1301 // Create mesh to map to
1302 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1303 mesh::PtrData outData = outMesh->createData("OutData", 2, 1_dataID);
1304 int outDataID = outData->getID();
1305 mesh::Vertex &vertex = outMesh->createVertex(Vector2d(0, 0));
1306 outMesh->allocateDataValues();
1307 addGlobalIndex(outMesh);
1308
1309 // Setup mapping with mapping coordinates and geometry used
1310 mapping.setMeshes(inMesh, outMesh);
1311 BOOST_TEST(mapping.hasComputedMapping() == false);
1312
1313 vertex.setCoords(Vector2d(0.0, 0.0));
1314 mapping.computeMapping();
1315 Eigen::VectorXd guess;
1316 mapping.map(inDataID, outDataID, guess);
1317 double value1 = outData->values()(0);
1318 double value2 = outData->values()(1);
1319 BOOST_TEST(mapping.hasComputedMapping() == true);
1320 BOOST_TEST(value1 == 1.0);
1321 BOOST_TEST(value2 == 4.0);
1322 BOOST_TEST(guess.size() > 0);
1323 guess.resize(0);
1324
1325 vertex.setCoords(Vector2d(0.0, 0.5));
1326 mapping.computeMapping();
1327 mapping.map(inDataID, outDataID, guess);
1328 value1 = outData->values()(0);
1329 value2 = outData->values()(1);
1330 BOOST_TEST(mapping.hasComputedMapping() == true);
1331 BOOST_TEST(value1 == 1.0);
1332 BOOST_TEST(value2 == 4.0);
1333 BOOST_TEST(guess.size() > 0);
1334 guess.resize(0);
1335
1336 vertex.setCoords(Vector2d(0.0, 1.0));
1337 mapping.computeMapping();
1338 mapping.map(inDataID, outDataID, guess);
1339 value1 = outData->values()(0);
1340 value2 = outData->values()(1);
1341 BOOST_TEST(mapping.hasComputedMapping() == true);
1342 BOOST_TEST(value1 == 1.0);
1343 BOOST_TEST(value2 == 4.0);
1344 BOOST_TEST(guess.size() > 0);
1345 guess.resize(0);
1346
1347 vertex.setCoords(Vector2d(1.0, 0.0));
1348 mapping.computeMapping();
1349 mapping.map(inDataID, outDataID, guess);
1350 value1 = outData->values()(0);
1351 value2 = outData->values()(1);
1352 BOOST_TEST(mapping.hasComputedMapping() == true);
1353 BOOST_TEST(value1 == 2.0);
1354 BOOST_TEST(value2 == 5.0);
1355 BOOST_TEST(guess.size() > 0);
1356 guess.resize(0);
1357
1358 vertex.setCoords(Vector2d(1.0, 0.5));
1359 mapping.computeMapping();
1360 mapping.map(inDataID, outDataID, guess);
1361 value1 = outData->values()(0);
1362 value2 = outData->values()(1);
1363 BOOST_TEST(mapping.hasComputedMapping() == true);
1364 BOOST_TEST(value1 == 2.0);
1365 BOOST_TEST(value2 == 5.0);
1366 BOOST_TEST(guess.size() > 0);
1367 guess.resize(0);
1368
1369 vertex.setCoords(Vector2d(1.0, 1.0));
1370 mapping.computeMapping();
1371 mapping.map(inDataID, outDataID, guess);
1372 value1 = outData->values()(0);
1373 value2 = outData->values()(1);
1374 BOOST_TEST(mapping.hasComputedMapping() == true);
1375 BOOST_TEST(value1 == 2.0);
1376 BOOST_TEST(value2 == 5.0);
1377 BOOST_TEST(guess.size() > 0);
1378 guess.resize(0);
1379
1380 vertex.setCoords(Vector2d(0.5, 0.0));
1381 mapping.computeMapping();
1382 mapping.map(inDataID, outDataID, guess);
1383 value1 = outData->values()(0);
1384 value2 = outData->values()(1);
1385 BOOST_TEST(mapping.hasComputedMapping() == true);
1386 BOOST_TEST(value1 == 1.5);
1387 BOOST_TEST(value2 == 4.5);
1388 BOOST_TEST(guess.size() > 0);
1389 guess.resize(0);
1390
1391 vertex.setCoords(Vector2d(0.5, 0.5));
1392 mapping.computeMapping();
1393 mapping.map(inDataID, outDataID, guess);
1394 value1 = outData->values()(0);
1395 value2 = outData->values()(1);
1396 BOOST_TEST(mapping.hasComputedMapping() == true);
1397 BOOST_TEST(value1 == 1.5);
1398 BOOST_TEST(value2 == 4.5);
1399 BOOST_TEST(guess.size() > 0);
1400 guess.resize(0);
1401
1402 vertex.setCoords(Vector2d(0.5, 1.0));
1403 mapping.computeMapping();
1404 mapping.map(inDataID, outDataID, guess);
1405 value1 = outData->values()(0);
1406 value2 = outData->values()(1);
1407 BOOST_TEST(mapping.hasComputedMapping() == true);
1408 BOOST_TEST(value1 == 1.5);
1409 BOOST_TEST(value2 == 4.5);
1410 BOOST_TEST(guess.size() > 0);
1411}
1412
1414{
1415 int dimensions = 3;
1416
1417 // Create mesh to map from
1418 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1419 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1420 int inDataID = inData->getID();
1421 inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1422 inMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
1423 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
1424 inMesh->createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
1425 inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 1.0));
1426 inMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 1.0));
1427 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 1.0));
1428 inMesh->createVertex(Eigen::Vector3d(1.0, 1.0, 1.0));
1429 inMesh->allocateDataValues();
1430 addGlobalIndex(inMesh);
1431
1432 auto &values = inData->values();
1433 values << 1.0, 1.0, 1.0, 1.0, 2.0, 2.0, 2.0, 2.0;
1434
1435 // Create mesh to map to
1436 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1437 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1438 int outDataID = outData->getID();
1439 mesh::Vertex &vertex = outMesh->createVertex(Eigen::Vector3d::Zero());
1440 outMesh->allocateDataValues();
1441 addGlobalIndex(outMesh);
1442
1443 // Setup mapping with mapping coordinates and geometry used
1444 mapping.setMeshes(inMesh, outMesh);
1445 BOOST_TEST(mapping.hasComputedMapping() == false);
1446
1447 vertex.setCoords(Eigen::Vector3d(0.0, 0.0, 0.0));
1448 mapping.computeMapping();
1449 Eigen::VectorXd guess;
1450 mapping.map(inDataID, outDataID, guess);
1451 double value = outData->values()(0);
1452 BOOST_TEST(mapping.hasComputedMapping() == true);
1453 BOOST_TEST(value == 1.0);
1454 BOOST_TEST(guess.size() > 0);
1455 guess.resize(0);
1456
1457 vertex.setCoords(Eigen::Vector3d(0.0, 0.5, 0.0));
1458 mapping.computeMapping();
1459 mapping.map(inDataID, outDataID, guess);
1460 value = outData->values()(0);
1461 BOOST_TEST(mapping.hasComputedMapping() == true);
1462 BOOST_TEST(value == 1.0);
1463 BOOST_TEST(guess.size() > 0);
1464 guess.resize(0);
1465
1466 vertex.setCoords(Eigen::Vector3d(0.5, 0.5, 0.0));
1467 mapping.computeMapping();
1468 mapping.map(inDataID, outDataID, guess);
1469 value = outData->values()(0);
1470 BOOST_TEST(mapping.hasComputedMapping() == true);
1471 BOOST_TEST(value == 1.0);
1472 BOOST_TEST(guess.size() > 0);
1473 guess.resize(0);
1474
1475 vertex.setCoords(Eigen::Vector3d(1.0, 0.0, 0.0));
1476 mapping.computeMapping();
1477 mapping.map(inDataID, outDataID, guess);
1478 value = outData->values()(0);
1479 BOOST_TEST(mapping.hasComputedMapping() == true);
1480 BOOST_TEST(value == 1.0);
1481 BOOST_TEST(guess.size() > 0);
1482 guess.resize(0);
1483
1484 vertex.setCoords(Eigen::Vector3d(1.0, 1.0, 0.0));
1485 mapping.computeMapping();
1486 mapping.map(inDataID, outDataID, guess);
1487 value = outData->values()(0);
1488 BOOST_TEST(mapping.hasComputedMapping() == true);
1489 BOOST_TEST(value == 1.0);
1490 BOOST_TEST(guess.size() > 0);
1491 guess.resize(0);
1492
1493 vertex.setCoords(Eigen::Vector3d(0.0, 0.0, 1.0));
1494 mapping.computeMapping();
1495 mapping.map(inDataID, outDataID, guess);
1496 value = outData->values()(0);
1497 BOOST_TEST(mapping.hasComputedMapping() == true);
1498 BOOST_TEST(value == 2.0);
1499 BOOST_TEST(guess.size() > 0);
1500 guess.resize(0);
1501
1502 vertex.setCoords(Eigen::Vector3d(1.0, 0.0, 1.0));
1503 mapping.computeMapping();
1504 mapping.map(inDataID, outDataID, guess);
1505 value = outData->values()(0);
1506 BOOST_TEST(mapping.hasComputedMapping() == true);
1507 BOOST_TEST(value == 2.0);
1508 BOOST_TEST(guess.size() > 0);
1509 guess.resize(0);
1510
1511 vertex.setCoords(Eigen::Vector3d(1.0, 1.0, 1.0));
1512 mapping.computeMapping();
1513 mapping.map(inDataID, outDataID, guess);
1514 value = outData->values()(0);
1515 BOOST_TEST(mapping.hasComputedMapping() == true);
1516 BOOST_TEST(value == 2.0);
1517 BOOST_TEST(guess.size() > 0);
1518 guess.resize(0);
1519
1520 vertex.setCoords(Eigen::Vector3d(0.5, 0.5, 1.0));
1521 mapping.computeMapping();
1522 mapping.map(inDataID, outDataID, guess);
1523 value = outData->values()(0);
1524 BOOST_TEST(mapping.hasComputedMapping() == true);
1525 BOOST_TEST(value == 2.0);
1526 BOOST_TEST(guess.size() > 0);
1527 guess.resize(0);
1528
1529 vertex.setCoords(Eigen::Vector3d(0.0, 0.0, 0.5));
1530 mapping.computeMapping();
1531 mapping.map(inDataID, outDataID, guess);
1532 value = outData->values()(0);
1533 BOOST_TEST(mapping.hasComputedMapping() == true);
1534 BOOST_TEST(value, 1.5);
1535 BOOST_TEST(guess.size() > 0);
1536 guess.resize(0);
1537
1538 vertex.setCoords(Eigen::Vector3d(1.0, 0.0, 0.5));
1539 mapping.computeMapping();
1540 mapping.map(inDataID, outDataID, guess);
1541 value = outData->values()(0);
1542 BOOST_TEST(mapping.hasComputedMapping() == true);
1543 BOOST_TEST(value == 1.5);
1544 BOOST_TEST(guess.size() > 0);
1545 guess.resize(0);
1546
1547 vertex.setCoords(Eigen::Vector3d(0.0, 1.0, 0.5));
1548 mapping.computeMapping();
1549 mapping.map(inDataID, outDataID, guess);
1550 value = outData->values()(0);
1551 BOOST_TEST(mapping.hasComputedMapping() == true);
1552 BOOST_TEST(value == 1.5);
1553 BOOST_TEST(guess.size() > 0);
1554 guess.resize(0);
1555
1556 vertex.setCoords(Eigen::Vector3d(1.0, 1.0, 0.5));
1557 mapping.computeMapping();
1558 mapping.map(inDataID, outDataID, guess);
1559 value = outData->values()(0);
1560 BOOST_TEST(mapping.hasComputedMapping() == true);
1561 BOOST_TEST(value == 1.5);
1562 BOOST_TEST(guess.size() > 0);
1563 guess.resize(0);
1564
1565 vertex.setCoords(Eigen::Vector3d(0.5, 0.5, 0.5));
1566 mapping.computeMapping();
1567 mapping.map(inDataID, outDataID, guess);
1568 value = outData->values()(0);
1569 BOOST_TEST(mapping.hasComputedMapping() == true);
1570 BOOST_TEST(value == 1.5);
1571 BOOST_TEST(guess.size() > 0);
1572}
1573
1575{
1576 int dimensions = 2;
1577 using Eigen::Vector2d;
1578
1579 // Create mesh to map from
1580 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1581 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1582 int inDataID = inData->getID();
1583 auto &inV1 = inMesh->createVertex(Vector2d(0.0, 0.0));
1584 auto &inV2 = inMesh->createVertex(Vector2d(1.0, 0.0));
1585 auto &inV3 = inMesh->createVertex(Vector2d(1.0, 1.0));
1586 auto &inV4 = inMesh->createVertex(Vector2d(0.0, 1.0));
1587
1588 inMesh->createEdge(inV1, inV2);
1589 inMesh->createEdge(inV2, inV3);
1590 inMesh->createEdge(inV3, inV4);
1591 inMesh->createEdge(inV1, inV4);
1592
1593 inMesh->allocateDataValues();
1594 addGlobalIndex(inMesh);
1595
1596 auto &inValues = inData->values();
1597 inValues << 1.0, 2.0, 2.0, 1.0;
1598
1599 // Create mesh to map to
1600 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1601 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1602 int outDataID = outData->getID();
1603 auto &outV1 = outMesh->createVertex(Vector2d(0.0, 0.0));
1604 auto &outV2 = outMesh->createVertex(Vector2d(0.0, 1.0));
1605 auto &outV3 = outMesh->createVertex(Vector2d(1.1, 1.1));
1606 auto &outV4 = outMesh->createVertex(Vector2d(0.1, 1.1));
1607 outMesh->createEdge(outV1, outV2);
1608 outMesh->createEdge(outV2, outV3);
1609 outMesh->createEdge(outV3, outV4);
1610 outMesh->createEdge(outV1, outV4);
1611 outMesh->allocateDataValues();
1612 addGlobalIndex(outMesh);
1613
1614 // Setup mapping with mapping coordinates and geometry used
1615 mapping.setMeshes(inMesh, outMesh);
1616 BOOST_TEST(mapping.hasComputedMapping() == false);
1617
1618 mapping.computeMapping();
1619 Eigen::VectorXd guess;
1620 mapping.map(inDataID, outDataID, guess);
1621 BOOST_TEST(guess.size() > 0);
1622
1623 testSerialScaledConsistent(inMesh, outMesh, inData, outData);
1624}
1625
1627{
1628 int dimensions = 3;
1629
1630 // Create mesh to map from
1631 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1632 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1633 int inDataID = inData->getID();
1634 auto &inV1 = inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1635 auto &inV2 = inMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
1636 auto &inV3 = inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.5));
1637 auto &inV4 = inMesh->createVertex(Eigen::Vector3d(2.0, 0.0, 0.0));
1638 auto &inV5 = inMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 0.0));
1639 auto &inV6 = inMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 1.0));
1640 auto &inE1 = inMesh->createEdge(inV1, inV2);
1641 auto &inE2 = inMesh->createEdge(inV2, inV3);
1642 auto &inE3 = inMesh->createEdge(inV1, inV3);
1643 auto &inE4 = inMesh->createEdge(inV4, inV5);
1644 auto &inE5 = inMesh->createEdge(inV5, inV6);
1645 auto &inE6 = inMesh->createEdge(inV4, inV6);
1646 inMesh->createTriangle(inE1, inE2, inE3);
1647 inMesh->createTriangle(inE4, inE5, inE6);
1648
1649 inMesh->allocateDataValues();
1650 addGlobalIndex(inMesh);
1651
1652 auto &inValues = inData->values();
1653 inValues << 1.0, 2.0, 4.0, 6.0, 8.0, 9.0;
1654
1655 // Create mesh to map to
1656 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1657 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1658 int outDataID = outData->getID();
1659 auto &outV1 = outMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1660 auto &outV2 = outMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
1661 auto &outV3 = outMesh->createVertex(Eigen::Vector3d(0.0, 1.1, 0.6));
1662 auto &outE1 = outMesh->createEdge(outV1, outV2);
1663 auto &outE2 = outMesh->createEdge(outV2, outV3);
1664 auto &outE3 = outMesh->createEdge(outV1, outV3);
1665 outMesh->createTriangle(outE1, outE2, outE3);
1666
1667 outMesh->allocateDataValues();
1668 addGlobalIndex(outMesh);
1669
1670 // Setup mapping with mapping coordinates and geometry used
1671 mapping.setMeshes(inMesh, outMesh);
1672 BOOST_TEST(mapping.hasComputedMapping() == false);
1673 mapping.computeMapping();
1674 BOOST_TEST(mapping.hasComputedMapping() == true);
1675 Eigen::VectorXd guess;
1676 mapping.map(inDataID, outDataID, guess);
1677 BOOST_TEST(guess.size() > 0);
1678
1679 testSerialScaledConsistent(inMesh, outMesh, inData, outData);
1680}
1681
1683{
1684 const int dimensions = 2;
1685 const double tolerance = 1e-6;
1686 using Eigen::Vector2d;
1687
1688 // Create mesh to map from
1689 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1690 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1691 int inDataID = inData->getID();
1692 mesh::Vertex &vertex0 = inMesh->createVertex(Vector2d(0, 0));
1693 mesh::Vertex &vertex1 = inMesh->createVertex(Vector2d(0, 0));
1694 inMesh->allocateDataValues();
1695 inData->values() << 1.0, 2.0;
1696 addGlobalIndex(inMesh);
1697
1698 // Create mesh to map to
1699 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1700 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1701 int outDataID = outData->getID();
1702 outMesh->createVertex(Vector2d(0.0, 0.0));
1703 outMesh->createVertex(Vector2d(1.0, 0.0));
1704 outMesh->createVertex(Vector2d(1.0, 1.0));
1705 outMesh->createVertex(Vector2d(0.0, 1.0));
1706 outMesh->allocateDataValues();
1707 addGlobalIndex(outMesh);
1708
1709 auto &values = outData->values();
1710
1711 mapping.setMeshes(inMesh, outMesh);
1712 BOOST_TEST(mapping.hasComputedMapping() == false);
1713
1714 vertex0.setCoords(Vector2d(0.5, 0.0));
1715 vertex1.setCoords(Vector2d(0.5, 1.0));
1716 mapping.computeMapping();
1717 Eigen::VectorXd guess;
1718 mapping.map(inDataID, outDataID, guess);
1719 BOOST_TEST(mapping.hasComputedMapping() == true);
1720 BOOST_TEST(testing::equals(values, Eigen::Vector4d(0.5, 0.5, 1.0, 1.0), tolerance));
1721 BOOST_TEST(guess.size() > 0);
1722 guess.resize(0);
1723
1724 vertex0.setCoords(Vector2d(0.0, 0.5));
1725 vertex1.setCoords(Vector2d(1.0, 0.5));
1726 mapping.computeMapping();
1727 mapping.map(inDataID, outDataID, guess);
1728 BOOST_TEST(mapping.hasComputedMapping() == true);
1729 BOOST_TEST(testing::equals(values, Eigen::Vector4d(0.5, 1.0, 1.0, 0.5), tolerance));
1730 BOOST_TEST(guess.size() > 0);
1731 guess.resize(0);
1732
1733 vertex0.setCoords(Vector2d(0.0, 1.0));
1734 vertex1.setCoords(Vector2d(1.0, 0.0));
1735 mapping.computeMapping();
1736 mapping.map(inDataID, outDataID, guess);
1737 BOOST_TEST(mapping.hasComputedMapping() == true);
1738 BOOST_TEST(testing::equals(values, Eigen::Vector4d(0.0, 2.0, 0.0, 1.0), tolerance));
1739 BOOST_TEST(guess.size() > 0);
1740 guess.resize(0);
1741
1742 vertex0.setCoords(Vector2d(0.0, 0.0));
1743 vertex1.setCoords(Vector2d(1.0, 1.0));
1744 mapping.computeMapping();
1745 mapping.map(inDataID, outDataID, guess);
1746 BOOST_TEST(mapping.hasComputedMapping() == true);
1747 BOOST_TEST(testing::equals(values, Eigen::Vector4d(1.0, 0.0, 2.0, 0.0), tolerance));
1748 BOOST_TEST(guess.size() > 0);
1749 guess.resize(0);
1750
1751 vertex0.setCoords(Vector2d(0.4, 0.5));
1752 vertex1.setCoords(Vector2d(0.6, 0.5));
1753 mapping.computeMapping();
1754 mapping.map(inDataID, outDataID, guess);
1755 BOOST_TEST(mapping.hasComputedMapping() == true);
1756 BOOST_TEST(values.sum() == 3.0);
1757 BOOST_TEST(guess.size() > 0);
1758}
1759
1761{
1762 const int dimensions = 2;
1763 const double tolerance = 1e-6;
1764 using Eigen::Vector2d;
1765
1766 // Create mesh to map from
1767 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1768 mesh::PtrData inData = inMesh->createData("InData", 2, 0_dataID);
1769 int inDataID = inData->getID();
1770 mesh::Vertex &vertex0 = inMesh->createVertex(Vector2d(0, 0));
1771 mesh::Vertex &vertex1 = inMesh->createVertex(Vector2d(0, 0));
1772 inMesh->allocateDataValues();
1773 inData->values() << 1.0, 4.0, 2.0, 5.0;
1774 addGlobalIndex(inMesh);
1775
1776 // Create mesh to map to
1777 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1778 mesh::PtrData outData = outMesh->createData("OutData", 2, 1_dataID);
1779 int outDataID = outData->getID();
1780 outMesh->createVertex(Vector2d(0.0, 0.0));
1781 outMesh->createVertex(Vector2d(1.0, 0.0));
1782 outMesh->createVertex(Vector2d(1.0, 1.0));
1783 outMesh->createVertex(Vector2d(0.0, 1.0));
1784 outMesh->allocateDataValues();
1785 addGlobalIndex(outMesh);
1786
1787 auto &values = outData->values();
1788
1789 mapping.setMeshes(inMesh, outMesh);
1790 BOOST_TEST(mapping.hasComputedMapping() == false);
1791
1792 vertex0.setCoords(Vector2d(0.5, 0.0));
1793 vertex1.setCoords(Vector2d(0.5, 1.0));
1794 mapping.computeMapping();
1795 Eigen::VectorXd guess;
1796 mapping.map(inDataID, outDataID, guess);
1797 BOOST_TEST(mapping.hasComputedMapping() == true);
1798 Eigen::VectorXd refValues(8);
1799 refValues << 0.5, 2, 0.5, 2, 1.0, 2.5, 1.0, 2.5;
1800 BOOST_TEST(testing::equals(values, refValues, tolerance));
1801 BOOST_TEST(guess.size() > 0);
1802 guess.resize(0);
1803
1804 vertex0.setCoords(Vector2d(0.0, 0.5));
1805 vertex1.setCoords(Vector2d(1.0, 0.5));
1806 mapping.computeMapping();
1807 mapping.map(inDataID, outDataID, guess);
1808 BOOST_TEST(mapping.hasComputedMapping() == true);
1809 refValues << 0.5, 2, 1.0, 2.5, 1.0, 2.5, 0.5, 2;
1810 BOOST_TEST(testing::equals(values, refValues, tolerance));
1811 BOOST_TEST(guess.size() > 0);
1812 guess.resize(0);
1813
1814 vertex0.setCoords(Vector2d(0.0, 1.0));
1815 vertex1.setCoords(Vector2d(1.0, 0.0));
1816 mapping.computeMapping();
1817 mapping.map(inDataID, outDataID, guess);
1818 BOOST_TEST(mapping.hasComputedMapping() == true);
1819 refValues << 0.0, 0.0, 2.0, 5.0, 0.0, 0.0, 1.0, 4.0;
1820 BOOST_TEST(testing::equals(values, refValues, tolerance));
1821 BOOST_TEST(guess.size() > 0);
1822 guess.resize(0);
1823
1824 vertex0.setCoords(Vector2d(0.0, 0.0));
1825 vertex1.setCoords(Vector2d(1.0, 1.0));
1826 mapping.computeMapping();
1827 mapping.map(inDataID, outDataID, guess);
1828 BOOST_TEST(mapping.hasComputedMapping() == true);
1829 refValues << 1.0, 4.0, 0.0, 0.0, 2.0, 5.0, 0.0, 0.0;
1830 BOOST_TEST(testing::equals(values, refValues, tolerance));
1831 BOOST_TEST(guess.size() > 0);
1832 guess.resize(0);
1833
1834 vertex0.setCoords(Vector2d(0.4, 0.5));
1835 vertex1.setCoords(Vector2d(0.6, 0.5));
1836 mapping.computeMapping();
1837 mapping.map(inDataID, outDataID, guess);
1838 BOOST_TEST(mapping.hasComputedMapping() == true);
1839 BOOST_TEST(values.sum() == 12.0);
1840 BOOST_TEST(guess.size() > 0);
1841}
1842
1844{
1845 using Eigen::Vector3d;
1846 int dimensions = 3;
1847
1848 // Create mesh to map from
1849 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
1850 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
1851 int inDataID = inData->getID();
1852 mesh::Vertex &vertex0 = inMesh->createVertex(Vector3d(0, 0, 0));
1853 mesh::Vertex &vertex1 = inMesh->createVertex(Vector3d(0, 0, 0));
1854 inMesh->allocateDataValues();
1855 inData->values() << 1.0, 2.0;
1856 addGlobalIndex(inMesh);
1857
1858 // Create mesh to map to
1859 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
1860 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
1861 int outDataID = outData->getID();
1862 outMesh->createVertex(Vector3d(0.0, 0.0, 0.0));
1863 outMesh->createVertex(Vector3d(1.0, 0.0, 0.0));
1864 outMesh->createVertex(Vector3d(1.0, 1.0, 0.0));
1865 outMesh->createVertex(Vector3d(0.0, 1.0, 0.0));
1866 outMesh->createVertex(Vector3d(0.0, 0.0, 1.0));
1867 outMesh->createVertex(Vector3d(1.0, 0.0, 1.0));
1868 outMesh->createVertex(Vector3d(1.0, 1.0, 1.0));
1869 outMesh->createVertex(Vector3d(0.0, 1.0, 1.0));
1870 outMesh->allocateDataValues();
1871 addGlobalIndex(outMesh);
1872
1873 auto &values = outData->values();
1874 double expectedSum = inData->values().sum();
1875
1876 mapping.setMeshes(inMesh, outMesh);
1877 BOOST_TEST(mapping.hasComputedMapping() == false);
1878
1879 vertex0.setCoords(Vector3d(0.5, 0.0, 0.0));
1880 vertex1.setCoords(Vector3d(0.5, 1.0, 0.0));
1881 mapping.computeMapping();
1882 Eigen::VectorXd guess;
1883 mapping.map(inDataID, outDataID, guess);
1884 BOOST_TEST(mapping.hasComputedMapping());
1885 BOOST_TEST(values.sum() == expectedSum);
1886 BOOST_TEST(guess.size() > 0);
1887}
1888
1889PRECICE_TEST_SETUP(1_rank, Require::PETSc)
1890BOOST_AUTO_TEST_CASE(MapThinPlateSplines)
1891{
1892 PRECICE_TEST();
1893 bool xDead = false;
1894 bool yDead = false;
1895 bool zDead = false;
1896 ThinPlateSplines fct;
1897 PetRadialBasisFctMapping<ThinPlateSplines> consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1898 perform2DTestConsistentMapping(consistentMap2D);
1899 PetRadialBasisFctMapping<ThinPlateSplines> consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
1900 perform3DTestConsistentMapping(consistentMap3D);
1901 PetRadialBasisFctMapping<ThinPlateSplines> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
1902 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
1903 PetRadialBasisFctMapping<ThinPlateSplines> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
1904 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
1905 PetRadialBasisFctMapping<ThinPlateSplines> conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1906 perform2DTestConservativeMapping(conservativeMap2D);
1907 PetRadialBasisFctMapping<ThinPlateSplines> conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
1908 perform3DTestConservativeMapping(conservativeMap3D);
1909}
1910
1911PRECICE_TEST_SETUP(1_rank, Require::PETSc)
1912BOOST_AUTO_TEST_CASE(MapMultiquadrics)
1913{
1914 PRECICE_TEST();
1915 bool xDead = false;
1916 bool yDead = false;
1917 bool zDead = false;
1918 Multiquadrics fct(1e-3);
1919 PetRadialBasisFctMapping<Multiquadrics> consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1920 perform2DTestConsistentMapping(consistentMap2D);
1921 PetRadialBasisFctMapping<Multiquadrics> consistentMap2DVector(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1922 perform2DTestConsistentMappingVector(consistentMap2DVector);
1923 PetRadialBasisFctMapping<Multiquadrics> consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
1924 perform3DTestConsistentMapping(consistentMap3D);
1925 PetRadialBasisFctMapping<Multiquadrics> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
1926 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
1927 PetRadialBasisFctMapping<Multiquadrics> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
1928 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
1929 PetRadialBasisFctMapping<Multiquadrics> conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1930 perform2DTestConservativeMapping(conservativeMap2D);
1931 PetRadialBasisFctMapping<Multiquadrics> conservativeMap2DVector(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1932 perform2DTestConservativeMappingVector(conservativeMap2DVector);
1933 PetRadialBasisFctMapping<Multiquadrics> conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
1934 perform3DTestConservativeMapping(conservativeMap3D);
1935}
1936
1937PRECICE_TEST_SETUP(1_rank, Require::PETSc)
1938BOOST_AUTO_TEST_CASE(MapInverseMultiquadrics)
1939{
1940 PRECICE_TEST();
1941 bool xDead = false;
1942 bool yDead = false;
1943 bool zDead = false;
1944 InverseMultiquadrics fct(1e-3);
1945 PetRadialBasisFctMapping<InverseMultiquadrics> consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1946 perform2DTestConsistentMapping(consistentMap2D);
1947 PetRadialBasisFctMapping<InverseMultiquadrics> consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
1948 perform3DTestConsistentMapping(consistentMap3D);
1949 PetRadialBasisFctMapping<InverseMultiquadrics> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
1950 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
1951 PetRadialBasisFctMapping<InverseMultiquadrics> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
1952 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
1953 PetRadialBasisFctMapping<InverseMultiquadrics> conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1954 perform2DTestConservativeMapping(conservativeMap2D);
1955 PetRadialBasisFctMapping<InverseMultiquadrics> conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
1956 perform3DTestConservativeMapping(conservativeMap3D);
1957}
1958
1959PRECICE_TEST_SETUP(1_rank, Require::PETSc)
1960BOOST_AUTO_TEST_CASE(MapVolumeSplines)
1961{
1962 PRECICE_TEST();
1963 bool xDead = false;
1964 bool yDead = false;
1965 bool zDead = false;
1966 VolumeSplines fct;
1967 PetRadialBasisFctMapping<VolumeSplines> consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1968 perform2DTestConsistentMapping(consistentMap2D);
1969 PetRadialBasisFctMapping<VolumeSplines> consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
1970 perform3DTestConsistentMapping(consistentMap3D);
1971 PetRadialBasisFctMapping<VolumeSplines> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
1972 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
1973 PetRadialBasisFctMapping<VolumeSplines> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
1974 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
1975 PetRadialBasisFctMapping<VolumeSplines> conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1976 perform2DTestConservativeMapping(conservativeMap2D);
1977 PetRadialBasisFctMapping<VolumeSplines> conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
1978 perform3DTestConservativeMapping(conservativeMap3D);
1979}
1980
1981PRECICE_TEST_SETUP(1_rank, Require::PETSc)
1983{
1984 PRECICE_TEST();
1985 bool xDead = false;
1986 bool yDead = false;
1987 bool zDead = false;
1988 Gaussian fct(1.0);
1989 PetRadialBasisFctMapping<Gaussian> consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
1990 perform2DTestConsistentMapping(consistentMap2D);
1991 PetRadialBasisFctMapping<Gaussian> consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
1992 perform3DTestConsistentMapping(consistentMap3D);
1993 PetRadialBasisFctMapping<Gaussian> scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
1994 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
1995 PetRadialBasisFctMapping<Gaussian> scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
1996 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
1997 PetRadialBasisFctMapping<Gaussian> conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
1998 perform2DTestConservativeMapping(conservativeMap2D);
1999 PetRadialBasisFctMapping<Gaussian> conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
2000 perform3DTestConservativeMapping(conservativeMap3D);
2001}
2002
2003PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2004BOOST_AUTO_TEST_CASE(MapCompactThinPlateSplinesC2)
2005{
2006 PRECICE_TEST();
2007 double supportRadius = 1.2;
2008 bool xDead = false;
2009 bool yDead = false;
2010 bool zDead = false;
2011 CompactThinPlateSplinesC2 fct(supportRadius);
2013 Mapping consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
2014 perform2DTestConsistentMapping(consistentMap2D);
2015 Mapping consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
2016 perform3DTestConsistentMapping(consistentMap3D);
2017 Mapping scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
2018 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
2019 Mapping scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
2020 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
2021 Mapping conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
2022 perform2DTestConservativeMapping(conservativeMap2D);
2023 Mapping conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
2024 perform3DTestConservativeMapping(conservativeMap3D);
2025}
2026
2027PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2028BOOST_AUTO_TEST_CASE(MapPetCompactPolynomialC0)
2029{
2030 PRECICE_TEST();
2031 double supportRadius = 1.2;
2032 bool xDead = false;
2033 bool yDead = false;
2034 bool zDead = false;
2035 CompactPolynomialC0 fct(supportRadius);
2037 Mapping consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
2038 perform2DTestConsistentMapping(consistentMap2D);
2039 Mapping consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
2040 perform3DTestConsistentMapping(consistentMap3D);
2041 Mapping scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
2042 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
2043 Mapping scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
2044 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
2045 Mapping conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
2046 perform2DTestConservativeMapping(conservativeMap2D);
2047 Mapping conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
2048 perform3DTestConservativeMapping(conservativeMap3D);
2049}
2050
2051PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2052BOOST_AUTO_TEST_CASE(MapPetCompactPolynomialC6)
2053{
2054 PRECICE_TEST();
2055 double supportRadius = 1.2;
2056 bool xDead = false;
2057 bool yDead = false;
2058 bool zDead = false;
2059 CompactPolynomialC6 fct(supportRadius);
2061 Mapping consistentMap2D(Mapping::CONSISTENT, 2, fct, {{xDead, yDead, zDead}});
2062 perform2DTestConsistentMapping(consistentMap2D);
2063 Mapping consistentMap3D(Mapping::CONSISTENT, 3, fct, {{xDead, yDead, zDead}});
2064 perform3DTestConsistentMapping(consistentMap3D);
2065 Mapping scaledConsistentMap2D(Mapping::SCALED_CONSISTENT_SURFACE, 2, fct, {{xDead, yDead, zDead}});
2066 perform2DTestScaledConsistentMapping(scaledConsistentMap2D);
2067 Mapping scaledConsistentMap3D(Mapping::SCALED_CONSISTENT_SURFACE, 3, fct, {{xDead, yDead, zDead}});
2068 perform3DTestScaledConsistentMapping(scaledConsistentMap3D);
2069 Mapping conservativeMap2D(Mapping::CONSERVATIVE, 2, fct, {{xDead, yDead, zDead}});
2070 perform2DTestConservativeMapping(conservativeMap2D);
2071 Mapping conservativeMap3D(Mapping::CONSERVATIVE, 3, fct, {{xDead, yDead, zDead}});
2072 perform3DTestConservativeMapping(conservativeMap3D);
2073}
2074
2075PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2077{
2078 PRECICE_TEST();
2079 using Eigen::Vector2d;
2080 int dimensions = 2;
2081
2082 bool xDead = false;
2083 bool yDead = true;
2084 bool zDead = false;
2085
2086 ThinPlateSplines fct;
2088 {{xDead, yDead, zDead}});
2089
2090 // Create mesh to map from
2091 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
2092 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2093 int inDataID = inData->getID();
2094 inMesh->createVertex(Vector2d(0.0, 1.0));
2095 inMesh->createVertex(Vector2d(1.0, 1.0));
2096 inMesh->createVertex(Vector2d(2.0, 1.0));
2097 inMesh->createVertex(Vector2d(3.0, 1.0));
2098 inMesh->allocateDataValues();
2099 addGlobalIndex(inMesh);
2100
2101 auto &values = inData->values();
2102 values << 1.0, 2.0, 2.0, 1.0;
2103
2104 // Create mesh to map to
2105 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
2106 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2107 int outDataID = outData->getID();
2108 mesh::Vertex &vertex = outMesh->createVertex(Vector2d(0, 0));
2109 outMesh->allocateDataValues();
2110 addGlobalIndex(outMesh);
2111
2112 // Setup mapping with mapping coordinates and geometry used
2113 mapping.setMeshes(inMesh, outMesh);
2114 BOOST_TEST(mapping.hasComputedMapping() == false);
2115
2116 vertex.setCoords(Vector2d(0.0, 3.0));
2117 mapping.computeMapping();
2118 Eigen::VectorXd guess;
2119 mapping.map(inDataID, outDataID, guess);
2120 double value = outData->values()(0);
2121 BOOST_TEST(mapping.hasComputedMapping() == true);
2122 BOOST_TEST(value == 1.0);
2123 BOOST_TEST(guess.size() > 0);
2124}
2125
2126PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2128{
2129 PRECICE_TEST();
2130 using Eigen::Vector3d;
2131 int dimensions = 3;
2132
2133 double supportRadius = 1.2;
2134 CompactPolynomialC6 fct(supportRadius);
2135 bool xDead = false;
2136 bool yDead = true;
2137 bool zDead = false;
2139 Mapping mapping(Mapping::CONSISTENT, dimensions, fct, {{xDead, yDead, zDead}});
2140
2141 // Create mesh to map from
2142 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
2143 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2144 int inDataID = inData->getID();
2145 inMesh->createVertex(Vector3d(0.0, 3.0, 0.0));
2146 inMesh->createVertex(Vector3d(1.0, 3.0, 0.0));
2147 inMesh->createVertex(Vector3d(0.0, 3.0, 1.0));
2148 inMesh->createVertex(Vector3d(1.0, 3.0, 1.0));
2149 inMesh->allocateDataValues();
2150 addGlobalIndex(inMesh);
2151
2152 auto &values = inData->values();
2153 values << 1.0, 2.0, 3.0, 4.0;
2154
2155 // Create mesh to map to
2156 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
2157 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2158 int outDataID = outData->getID();
2159 outMesh->createVertex(Vector3d(0.0, 2.9, 0.0));
2160 outMesh->createVertex(Vector3d(0.8, 2.9, 0.1));
2161 outMesh->createVertex(Vector3d(0.1, 2.9, 0.9));
2162 outMesh->createVertex(Vector3d(1.1, 2.9, 1.1));
2163 outMesh->allocateDataValues();
2164 addGlobalIndex(outMesh);
2165
2166 // Setup mapping with mapping coordinates and geometry used
2167 mapping.setMeshes(inMesh, outMesh);
2168 BOOST_TEST(mapping.hasComputedMapping() == false);
2169
2170 mapping.computeMapping();
2171 Eigen::VectorXd guess;
2172 mapping.map(inDataID, outDataID, guess);
2173 BOOST_TEST(mapping.hasComputedMapping() == true);
2174 BOOST_TEST(guess.size() > 0);
2175 guess.resize(0);
2176
2177 BOOST_TEST(outData->values()(0) == 1.0);
2178 BOOST_TEST(outData->values()(1) == 2.0);
2179 BOOST_TEST(outData->values()(2) == 2.9);
2180 BOOST_TEST(outData->values()(3) == 4.3);
2181}
2182
2183PRECICE_TEST_SETUP(1_rank, Require::PETSc);
2184BOOST_AUTO_TEST_CASE(ConsistentPolynomialSwitch,
2185 *boost::unit_test::tolerance(1e-6))
2186{
2187 PRECICE_TEST();
2188 using Eigen::Vector2d;
2189 int dimensions = 2;
2190
2191 bool xDead = false, yDead = false, zDead = false;
2192
2193 Gaussian fct(1); // supportRadius = 4.55
2194
2195 // Create mesh to map from
2196 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
2197 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2198 int inDataID = inData->getID();
2199 inMesh->createVertex(Vector2d(1, 1));
2200 inMesh->createVertex(Vector2d(1, 0));
2201 inMesh->createVertex(Vector2d(0, 0));
2202 inMesh->createVertex(Vector2d(0, 1));
2203 inMesh->allocateDataValues();
2204 addGlobalIndex(inMesh);
2205 inData->values() << 1, 1, 1, 1;
2206
2207 // Create mesh to map to
2208 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
2209 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2210 int outDataID = outData->getID();
2211 outMesh->createVertex(Vector2d(3, 3)); // Point is outside the inMesh
2212
2213 outMesh->allocateDataValues();
2214 addGlobalIndex(outMesh);
2215
2216 // Test deactivated polynomial
2217 PetRadialBasisFctMapping<Gaussian> mappingOff(Mapping::CONSISTENT, dimensions, fct,
2218 {{xDead, yDead, zDead}},
2219 1e-9, Polynomial::OFF);
2220 mappingOff.setMeshes(inMesh, outMesh);
2221 mappingOff.computeMapping();
2222 Eigen::VectorXd guess;
2223 mappingOff.map(inDataID, outDataID, guess);
2224
2225 BOOST_TEST(outData->values()(0) <= 0.01); // Mapping to almost 0 since almost no basis function at (3,3) and no polynomial
2226 BOOST_TEST(guess.size() > 0);
2227 guess.resize(0);
2228
2229 // Test integrated polynomial
2231 {{xDead, yDead, zDead}},
2232 1e-9, Polynomial::ON);
2233
2234 mappingOn.setMeshes(inMesh, outMesh);
2235 mappingOn.computeMapping();
2236 mappingOn.map(inDataID, outDataID, guess);
2237
2238 BOOST_TEST(outData->values()(0) == 1.0); // Mapping to 1 since there is the polynomial
2239 BOOST_TEST(guess.size() > 0);
2240 guess.resize(0);
2241
2242 // Test separated polynomial
2243 PetRadialBasisFctMapping<Gaussian> mappingSep(Mapping::CONSISTENT, dimensions, fct,
2244 {{xDead, yDead, zDead}},
2245 1e-9, Polynomial::SEPARATE);
2246
2247 mappingSep.setMeshes(inMesh, outMesh);
2248 mappingSep.computeMapping();
2249 mappingSep.map(inDataID, outDataID, guess);
2250
2251 BOOST_TEST(outData->values()(0) == 1.0); // Mapping to 1 since there is the polynomial
2252 BOOST_TEST(guess.size() > 0);
2253}
2254
2255PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2256BOOST_AUTO_TEST_CASE(ConservativePolynomialSwitch,
2257 *boost::unit_test::tolerance(1e-6))
2258{
2259 PRECICE_TEST();
2260 using Eigen::Vector2d;
2261 int dimensions = 2;
2262
2263 bool xDead = false, yDead = false, zDead = false;
2264
2265 Gaussian fct(1); // supportRadius = 4.55
2266
2267 // Create mesh to map from
2268 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
2269 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2270 int inDataID = inData->getID();
2271 inMesh->createVertex(Vector2d(0, 0));
2272 inMesh->createVertex(Vector2d(1, 0));
2273 inMesh->createVertex(Vector2d(1, 1));
2274 inMesh->createVertex(Vector2d(0, 1));
2275 inMesh->allocateDataValues();
2276 addGlobalIndex(inMesh);
2277 inData->values() << 1, 1, 1, 1;
2278
2279 // Create mesh to map to
2280 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
2281 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2282 int outDataID = outData->getID();
2283 outMesh->createVertex(Vector2d(0.4, 0));
2284 outMesh->createVertex(Vector2d(6, 6));
2285 outMesh->createVertex(Vector2d(7, 7));
2286
2287 outMesh->allocateDataValues();
2288 addGlobalIndex(outMesh);
2289
2290 // Test deactivated polynomial
2292 {{xDead, yDead, zDead}},
2293 1e-9, Polynomial::OFF);
2294 mappingOff.setMeshes(inMesh, outMesh);
2295 mappingOff.computeMapping();
2296 Eigen::VectorXd guess;
2297 mappingOff.map(inDataID, outDataID, guess);
2298 BOOST_TEST(guess.size() > 0);
2299 guess.resize(0);
2300
2301 BOOST_TEST(outData->values()(0) == 2.119967); // Conservativness is not retained, because no polynomial
2302 BOOST_TEST(outData->values()(1) == 0.0); // Mapping to 0 since no basis function at (5,5) and no polynomial
2303 BOOST_TEST(outData->values()(2) == 0.0); // Mapping to 0 since no basis function at (5,5) and no polynomial
2304
2305 // Test integrated polynomial
2307 {{xDead, yDead, zDead}},
2308 1e-9, Polynomial::ON);
2309
2310 mappingOn.setMeshes(inMesh, outMesh);
2311 mappingOn.computeMapping();
2312 mappingOn.map(inDataID, outDataID, guess);
2313 BOOST_TEST(guess.size() > 0);
2314 guess.resize(0);
2315
2316 BOOST_TEST(outData->values()(0) == 0);
2317 BOOST_TEST(outData->values()(1) == 26.0);
2318 BOOST_TEST(outData->values()(2) == -22.0);
2319
2320 // Test separated polynomial
2322 {{xDead, yDead, zDead}},
2323 1e-9, Polynomial::SEPARATE);
2324
2325 mappingSep.setMeshes(inMesh, outMesh);
2326 mappingSep.computeMapping();
2327 mappingSep.map(inDataID, outDataID, guess);
2328 BOOST_TEST(guess.size() > 0);
2329 guess.resize(0);
2330
2331 BOOST_TEST(outData->values()(0) == 0);
2332 BOOST_TEST(outData->values()(1) == 26.0);
2333 BOOST_TEST(outData->values()(2) == -22.0);
2334}
2335
2336PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2338{
2339 PRECICE_TEST();
2340 /*
2341 * RATIONALE: Correctly destroying PETSc objects in OOP context can be a bit
2342 * tricky. We test if an RBF object can be destroyed right after creation
2343 * and if only computeMapping (not map) is called.
2344 */
2345
2346 ThinPlateSplines fct;
2347
2348 // Call neither computeMapping nor map
2349 {
2351 {{false, false, false}});
2352 }
2353
2354 {
2355 // Call only computeMapping
2356 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", 2, testing::nextMeshID()));
2357 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2358 inMesh->createVertex(Eigen::Vector2d(0, 0));
2359 inMesh->allocateDataValues();
2360 addGlobalIndex(inMesh);
2361
2362 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", 2, testing::nextMeshID()));
2363 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2364 outMesh->createVertex(Eigen::Vector2d(0, 0));
2365 outMesh->allocateDataValues();
2366 addGlobalIndex(outMesh);
2367
2369 {{false, false, false}});
2370
2371 mapping2.setMeshes(inMesh, outMesh);
2372 mapping2.computeMapping();
2373 }
2374}
2375
2376PRECICE_TEST_SETUP(1_rank, Require::PETSc)
2377BOOST_AUTO_TEST_CASE(TestNonHomongenousGlobalIndex)
2378{
2379 PRECICE_TEST();
2380 using Eigen::Vector2d;
2381 int dimensions = 2;
2382
2383 bool xDead = false, yDead = false, zDead = false;
2384
2385 Gaussian fct(1); // supportRadius = 4.55
2386
2387 // Create mesh to map from
2388 mesh::PtrMesh inMesh(new mesh::Mesh("InMesh", dimensions, testing::nextMeshID()));
2389 mesh::PtrData inData = inMesh->createData("InData", 1, 0_dataID);
2390 int inDataID = inData->getID();
2391 inMesh->createVertex(Vector2d(1, 1)).setGlobalIndex(2);
2392 inMesh->createVertex(Vector2d(1, 0)).setGlobalIndex(3);
2393 inMesh->createVertex(Vector2d(0, 0)).setGlobalIndex(6);
2394 inMesh->createVertex(Vector2d(0, 1)).setGlobalIndex(5);
2395 inMesh->allocateDataValues();
2396
2397 inData->values() << 1, 1, 1, 1;
2398
2399 // Create mesh to map to
2400 mesh::PtrMesh outMesh(new mesh::Mesh("OutMesh", dimensions, testing::nextMeshID()));
2401 mesh::PtrData outData = outMesh->createData("OutData", 1, 1_dataID);
2402 int outDataID = outData->getID();
2403 outMesh->createVertex(Vector2d(0.5, 0.5));
2404
2405 outMesh->allocateDataValues();
2406 addGlobalIndex(outMesh);
2407
2409 {{xDead, yDead, zDead}});
2410 mapping1.setMeshes(inMesh, outMesh);
2411 mapping1.computeMapping();
2412 Eigen::VectorXd guess;
2413 mapping1.map(inDataID, outDataID, guess);
2414 BOOST_TEST(guess.size() > 0);
2415 guess.resize(0);
2416
2417 BOOST_TEST(outData->values()(0) == 1);
2418
2420 {{xDead, yDead, zDead}});
2421 inData->values() << 0, 0, 0, 0; // reset
2422 outData->values() << 4; // used as inData here
2423 mapping2.setMeshes(outMesh, inMesh);
2424 mapping2.computeMapping();
2425 mapping2.map(outDataID, inDataID, guess);
2426
2427 BOOST_TEST(inData->values()(0) == 1.0);
2428 BOOST_TEST(inData->values()(1) == 1.0);
2429 BOOST_TEST(inData->values()(2) == 1.0);
2430 BOOST_TEST(inData->values()(3) == 1.0);
2431}
2432
2433BOOST_AUTO_TEST_SUITE_END() // Serial
2434
2435BOOST_AUTO_TEST_SUITE_END() // PetRadialBasisFunctionMapping
2437
2438#endif
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
void testTagging(const TestContext &context, MeshSpecification inMeshSpec, MeshSpecification outMeshSpec, MeshSpecification shouldTagFirstRound, MeshSpecification shouldTagSecondRound, bool consistent)
void perform3DTestScaledConsistentMapping(Mapping &mapping)
void addGlobalIndex(mesh::PtrMesh &mesh, int offset=0)
void perform2DTestConservativeMappingVector(Mapping &mapping)
void perform2DTestScaledConsistentMapping(Mapping &mapping)
void perform2DTestConsistentMappingVector(Mapping &mapping)
void perform2DTestConsistentMapping(Mapping &mapping)
void perform3DTestConsistentMapping(Mapping &mapping)
void perform3DTestConservativeMapping(Mapping &mapping)
void perform2DTestConservativeMapping(Mapping &mapping)
void testDistributed(const TestContext &context, Mapping &mapping, MeshSpecification inMeshSpec, MeshSpecification outMeshSpec, ReferenceSpecification referenceSpec, int inGlobalIndexOffset=0)
void getDistributedMesh(const TestContext &context, MeshSpecification const &vertices, mesh::PtrMesh &mesh, mesh::PtrData &data, int globalIndexOffset=0)
void testSerialScaledConsistent(mesh::PtrMesh inMesh, mesh::PtrMesh outMesh, PtrData inData, PtrData outData)
BOOST_AUTO_TEST_CASE(MapThinPlateSplines)
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 back(T... args)
T begin(T... args)
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 using the Petsc library to solve the resulting system.
Radial basis function with global support.
Radial basis function with global support.
DataID getID() const
Returns the ID of the data set (supposed to be unique).
Definition Data.cpp:110
Eigen::VectorXd & values()
Returns a reference to the data values.
Definition Data.cpp:28
Container and creator for meshes.
Definition Mesh.hpp:38
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:121
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
PtrData & createData(const std::string &name, int dimension, DataID id, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Create only data for vertex.
Definition Mesh.cpp:153
void computeBoundingBox()
Computes the boundingBox for the vertices.
Definition Mesh.cpp:244
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:113
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:235
Vertex of a mesh.
Definition Vertex.hpp:16
void setCoords(const VECTOR_T &coordinates)
Sets the coordinates of the vertex.
Definition Vertex.hpp:101
void setGlobalIndex(int globalIndex)
Definition Vertex.cpp:17
void setOwner(bool owner)
Definition Vertex.cpp:27
Rank rank
the rank of the current participant
T count(T... args)
T data(T... args)
T emplace(T... args)
T end(T... args)
contains data mapping from points to meshes.
provides Mesh, Data and primitives.
Eigen::VectorXd integrateSurface(const PtrMesh &mesh, const Eigen::VectorXd &input)
Given the data and the mesh, this function returns the surface integral. Assumes no overlap exists fo...
Definition Utils.cpp:11
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
Main namespace of the precice library.
T size(T... args)
Holds rank, owner, position and value of a single vertex.
const Eigen::Map< const Eigen::VectorXd > asEigen() const
static constexpr bool hasCompactSupport()