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