preCICE v3.1.2
Loading...
Searching...
No Matches
NearestNeighborMappingTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <memory>
4#include <string_view>
6#include "mapping/Mapping.hpp"
8#include "math/constants.hpp"
9#include "mesh/Data.hpp"
10#include "mesh/Mesh.hpp"
12#include "mesh/Utils.hpp"
13#include "mesh/Vertex.hpp"
15#include "testing/Testing.hpp"
16#include "time/Sample.hpp"
17
18using namespace precice;
19using namespace precice::mesh;
20
21BOOST_AUTO_TEST_SUITE(MappingTests)
22BOOST_AUTO_TEST_SUITE(NearestNeighborMapping)
23
24BOOST_AUTO_TEST_CASE(ConsistentNonIncremental)
25{
26 PRECICE_TEST(1_rank);
27 int dimensions = 2;
28 using testing::equals;
29
30 // Create mesh to map from
31 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
32 Vertex &inVertex0 = inMesh->createVertex(Eigen::Vector2d::Constant(0.0));
33 Vertex &inVertex1 = inMesh->createVertex(Eigen::Vector2d::Constant(1.0));
34
35 Eigen::VectorXd inValuesScalar = Eigen::VectorXd::Zero(2);
36 Eigen::VectorXd inValuesVector = Eigen::VectorXd::Zero(4);
37 inValuesScalar << 1.0, 2.0;
38 inValuesVector << 1.0, 2.0, 3.0, 4.0;
39
40 // Create mesh to map to
41 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
42 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector2d::Constant(0.0));
43 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector2d::Constant(1.0));
44
45 // Setup mapping with mapping coordinates and geometry used
47 mapping.setMeshes(inMesh, outMesh);
48 BOOST_TEST(mapping.hasComputedMapping() == false);
49
50 // Map data with coinciding vertices, has to result in equal values.
51 mapping.computeMapping();
52 Eigen::VectorXd outValuesScalar = Eigen::VectorXd::Zero(2);
53 time::Sample inSample(1, inValuesScalar);
54 mapping.map(inSample, outValuesScalar);
55 BOOST_TEST(mapping.hasComputedMapping() == true);
56 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0));
57 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1));
58 Eigen::VectorXd outValuesVector = Eigen::VectorXd::Zero(4);
59 inSample = time::Sample(2, inValuesVector);
60 mapping.map(inSample, outValuesVector);
61 BOOST_CHECK(equals(inValuesVector, outValuesVector));
62
63 // Map data with almost coinciding vertices, has to result in equal values.
64 inVertex0.setCoords(outVertex0.getCoords() + Eigen::Vector2d::Constant(0.1));
65 inVertex1.setCoords(outVertex1.getCoords() + Eigen::Vector2d::Constant(0.1));
66 mapping.computeMapping();
67 inSample = time::Sample(1, inValuesScalar);
68 mapping.map(inSample, outValuesScalar);
69 BOOST_TEST(mapping.hasComputedMapping() == true);
70 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0));
71 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1));
72 inSample = time::Sample(2, inValuesVector);
73 mapping.map(inSample, outValuesVector);
74 BOOST_CHECK(equals(inValuesVector, outValuesVector));
75
76 // Map data with exchanged vertices, has to result in exchanged values.
77 inVertex0.setCoords(outVertex1.getCoords());
78 inVertex1.setCoords(outVertex0.getCoords());
79 mapping.computeMapping();
80 inSample = time::Sample(1, inValuesScalar);
81 mapping.map(inSample, outValuesScalar);
82 BOOST_TEST(mapping.hasComputedMapping() == true);
83 BOOST_TEST(outValuesScalar(1) == inValuesScalar(0));
84 BOOST_TEST(outValuesScalar(0) == inValuesScalar(1));
85 inSample = time::Sample(2, inValuesVector);
86 mapping.map(inSample, outValuesVector);
87 Eigen::Vector4d expected(3.0, 4.0, 1.0, 2.0);
88 BOOST_CHECK(equals(expected, outValuesVector));
89
90 // Map data with coinciding output vertices, has to result in same values.
91 outVertex1.setCoords(outVertex0.getCoords());
92 mapping.computeMapping();
93 inSample = time::Sample(1, inValuesScalar);
94 mapping.map(inSample, outValuesScalar);
95 BOOST_TEST(mapping.hasComputedMapping() == true);
96 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1));
97 BOOST_TEST(outValuesScalar(0) == inValuesScalar(1));
98 inSample = time::Sample(2, inValuesVector);
99 mapping.map(inSample, outValuesVector);
100 expected << 3.0, 4.0, 3.0, 4.0;
101 BOOST_CHECK(equals(expected, outValuesVector));
102}
103
104BOOST_AUTO_TEST_CASE(ConservativeNonIncremental)
105{
106 PRECICE_TEST(1_rank);
107 int dimensions = 2;
108
109 // Create mesh to map from
110 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
111 Vertex & inVertex0 = inMesh->createVertex(Eigen::Vector2d::Constant(0.0));
112 Vertex & inVertex1 = inMesh->createVertex(Eigen::Vector2d::Constant(1.0));
113 Eigen::VectorXd inValues = Eigen::VectorXd::Zero(2);
114 inValues(0) = 1.0;
115 inValues(1) = 2.0;
116
117 // Create mesh to map to
118 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
119 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector2d::Constant(0.0));
120 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector2d::Constant(1.0));
121
122 // Setup mapping with mapping coordinates and geometry used
124 mapping.setMeshes(inMesh, outMesh);
125 BOOST_TEST(mapping.hasComputedMapping() == false);
126
127 // Map data with coinciding vertices, has to result in equal values.
128 mapping.computeMapping();
129 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(2);
130 time::Sample inSample(1, inValues);
131 mapping.map(inSample, outValues);
132 BOOST_TEST(mapping.hasComputedMapping() == true);
133 BOOST_TEST(outValues(0) == inValues(0));
134 BOOST_TEST(outValues(1) == inValues(1));
135 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
136
137 // Map data with almost coinciding vertices, has to result in equal values.
138 inVertex0.setCoords(outVertex0.getCoords() + Eigen::Vector2d::Constant(0.1));
139 inVertex1.setCoords(outVertex1.getCoords() + Eigen::Vector2d::Constant(0.1));
140 mapping.computeMapping();
141 inSample = time::Sample(1, inValues);
142 mapping.map(inSample, outValues);
143 BOOST_TEST(mapping.hasComputedMapping() == true);
144 BOOST_TEST(outValues(0) == inValues(0));
145 BOOST_TEST(outValues(1) == inValues(1));
146 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
147
148 // Map data with exchanged vertices, has to result in exchanged values.
149 inVertex0.setCoords(outVertex1.getCoords());
150 inVertex1.setCoords(outVertex0.getCoords());
151 mapping.computeMapping();
152 inSample = time::Sample(1, inValues);
153 mapping.map(inSample, outValues);
154 BOOST_TEST(mapping.hasComputedMapping() == true);
155 BOOST_TEST(outValues(1) == inValues(0));
156 BOOST_TEST(outValues(0) == inValues(1));
157 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
158
159 // Map data with coinciding output vertices, has to result in double values.
160 outVertex1.setCoords(Eigen::Vector2d::Constant(-1.0));
161 mapping.computeMapping();
162 inSample = time::Sample(1, inValues);
163 mapping.map(inSample, outValues);
164 BOOST_TEST(mapping.hasComputedMapping() == true);
165 BOOST_TEST(outValues(0) == inValues(0) + inValues(1));
166 BOOST_TEST(outValues(1) == 0.0);
167}
168
169BOOST_AUTO_TEST_CASE(ScaledConsistentNonIncremental)
170{
171 PRECICE_TEST(1_rank);
172 int dimensions = 2;
173
174 // Create mesh to map from
175 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
176 Vertex &inVertex0 = inMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
177 Vertex &inVertex1 = inMesh->createVertex(Eigen::Vector2d{1.0, 0.0});
178 Vertex &inVertex2 = inMesh->createVertex(Eigen::Vector2d{3.0, 0.0});
179 Vertex &inVertex3 = inMesh->createVertex(Eigen::Vector2d{6.0, 0.0});
180
181 inMesh->createEdge(inVertex0, inVertex1);
182 inMesh->createEdge(inVertex1, inVertex2);
183 inMesh->createEdge(inVertex2, inVertex3);
184
185 Eigen::VectorXd inValues = Eigen::VectorXd::Zero(4);
186 inValues(0) = 1.0;
187 inValues(1) = 2.0;
188 inValues(2) = 3.0;
189 inValues(3) = 4.0;
190
191 // Create mesh to map to
192 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
193 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
194 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector2d(0.8, 0.0));
195 Vertex &outVertex2 = outMesh->createVertex(Eigen::Vector2d(3.0, 0.0));
196 Vertex &outVertex3 = outMesh->createVertex(Eigen::Vector2d(6.2, 0.0));
197
198 outMesh->createEdge(outVertex0, outVertex1);
199 outMesh->createEdge(outVertex1, outVertex2);
200 outMesh->createEdge(outVertex2, outVertex3);
201
202 // Setup mapping with mapping coordinates and geometry used
204
205 mapping.setMeshes(inMesh, outMesh);
206 BOOST_TEST(mapping.hasComputedMapping() == false);
207
208 mapping.computeMapping();
209
210 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(4);
211 BOOST_TEST(mapping.hasComputedMapping() == true);
212 time::Sample inSample(1, inValues);
213 mapping.map(inSample, outValues);
214
215 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
216 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
217
218 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
219 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
220 }
221
222 double scaleFactor = outValues(0) / inValues(0);
223 BOOST_TEST(scaleFactor != 1.0);
224
225 BOOST_TEST(inValues(0) * scaleFactor == outValues(0));
226 BOOST_TEST(inValues(1) * scaleFactor == outValues(1));
227 BOOST_TEST(inValues(2) * scaleFactor == outValues(2));
228 BOOST_TEST(inValues(3) * scaleFactor == outValues(3));
229}
230
231namespace {
232template <int N>
233PtrMesh create2DLinSpaceMesh(std::string_view name, int dims, double x0, double x1)
234{
235 static_assert(N > 1, "More than 1 vertex required");
236 PtrMesh mesh(new Mesh("InMesh", dims, testing::nextMeshID()));
237 std::vector<Vertex *> vertices;
238 double dx = (x1 - x0) / (double) (N - 1);
239 for (int i = 0; i < N; ++i) {
240 double x = x0 + dx * i;
241 vertices.push_back(&mesh->createVertex(Eigen::Vector2d{x, 0.0}));
242 }
243 auto first = vertices.begin();
244 auto second = ++vertices.begin();
245 while (second != vertices.end()) {
246 mesh->createEdge(**first, **second);
247 ++first;
248 ++second;
249 }
250 return mesh;
251}
252} // namespace
253
254BOOST_AUTO_TEST_CASE(ScaledConsistentZeroData)
255{
256 PRECICE_TEST(1_rank);
257 int dimensions = 2;
258
259 // Create mesh to map from
260 PtrMesh inMesh = create2DLinSpaceMesh<4>("InMesh", dimensions, 0, 1);
261
262 // Create mesh to map to
263 PtrMesh outMesh = create2DLinSpaceMesh<3>("OutMesh", dimensions, 0, 1);
264
265 // Setup mapping with mapping coordinates and geometry used
267
268 mapping.setMeshes(inMesh, outMesh);
269 BOOST_TEST(mapping.hasComputedMapping() == false);
270 mapping.computeMapping();
271 BOOST_TEST(mapping.hasComputedMapping() == true);
272
273 Eigen::VectorXd inValues = Eigen::VectorXd::Zero(4);
274 time::Sample inSample(1, inValues);
275 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
276 mapping.map(inSample, outValues);
277
278 BOOST_TEST(!outValues.hasNaN());
279 BOOST_TEST(outValues.isZero());
280
281 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
282 BOOST_TEST(!inputIntegral.hasNaN());
283 BOOST_TEST(inputIntegral.isZero());
284
285 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
286 BOOST_TEST(!outputIntegral.hasNaN());
287 BOOST_TEST(outputIntegral.isZero());
288}
289
290BOOST_AUTO_TEST_CASE(ScaledConsistentZeroIntegral)
291{
292 PRECICE_TEST(1_rank);
293 int dimensions = 2;
294
295 // Create mesh to map from
296 PtrMesh inMesh = create2DLinSpaceMesh<4>("InMesh", dimensions, 0, 1);
297
298 // Create mesh to map to
299 PtrMesh outMesh = create2DLinSpaceMesh<3>("OutMesh", dimensions, 0, 1);
300
301 // Setup mapping with mapping coordinates and geometry used
303
304 mapping.setMeshes(inMesh, outMesh);
305 BOOST_TEST(mapping.hasComputedMapping() == false);
306 mapping.computeMapping();
307 BOOST_TEST(mapping.hasComputedMapping() == true);
308
309 Eigen::Vector4d inValues{1, 1, -1, -1};
310 time::Sample inSample(1, inValues);
311 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
312 mapping.map(inSample, outValues);
313
314 BOOST_TEST(!outValues.hasNaN());
315 BOOST_TEST((outValues.array() == 0.0).count() == 0);
316
317 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
318 BOOST_TEST(!inputIntegral.hasNaN());
319 BOOST_TEST(inputIntegral(0) == 0.0);
320
321 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
322 BOOST_TEST(!outputIntegral.hasNaN());
323 BOOST_TEST(outputIntegral(0) == 0.0);
324}
325
326BOOST_AUTO_TEST_CASE(ScaledConsistentZeroDataComponent)
327{
328 PRECICE_TEST(1_rank);
329 int dimensions = 2;
330
331 // Create mesh to map from
332 PtrMesh inMesh = create2DLinSpaceMesh<4>("InMesh", dimensions, 0, 1);
333
334 // Create mesh to map to
335 PtrMesh outMesh = create2DLinSpaceMesh<3>("OutMesh", dimensions, 0, 1);
336
337 // Setup mapping with mapping coordinates and geometry used
339
340 mapping.setMeshes(inMesh, outMesh);
341 BOOST_TEST(mapping.hasComputedMapping() == false);
342 mapping.computeMapping();
343 BOOST_TEST(mapping.hasComputedMapping() == true);
344
345 Eigen::VectorXd inValues(8);
346 inValues << 1, 0, 0.5, 0, 1.5, 0, 1, 0;
347 time::Sample inSample(2, inValues);
348 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(6);
349 mapping.map(inSample, outValues);
350
351 BOOST_TEST(!outValues.hasNaN());
352 BOOST_TEST((outValues.array() == 0.0).count() == 3);
353
354 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
355 BOOST_TEST(!inputIntegral.hasNaN());
356 BOOST_TEST(inputIntegral(0) == 1.0);
357 BOOST_TEST(inputIntegral(1) == 0.0);
358
359 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
360 BOOST_TEST(!outputIntegral.hasNaN());
361 BOOST_TEST(outputIntegral(0) == 1.0);
362 BOOST_TEST(outputIntegral(1) == 0.0);
363}
364
365BOOST_AUTO_TEST_CASE(ScaledConsistentZeroIntegralComponent)
366{
367 PRECICE_TEST(1_rank);
368 int dimensions = 2;
369
370 // Create mesh to map from
371 PtrMesh inMesh = create2DLinSpaceMesh<4>("InMesh", dimensions, 0, 1);
372
373 // Create mesh to map to
374 PtrMesh outMesh = create2DLinSpaceMesh<3>("OutMesh", dimensions, 0, 1);
375
376 // Setup mapping with mapping coordinates and geometry used
378
379 mapping.setMeshes(inMesh, outMesh);
380 BOOST_TEST(mapping.hasComputedMapping() == false);
381 mapping.computeMapping();
382 BOOST_TEST(mapping.hasComputedMapping() == true);
383
384 Eigen::VectorXd inValues(8);
385 inValues << 2, 1, 3, 1, 1, -1, 2, -1;
386 time::Sample inSample(2, inValues);
387 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(6);
388 mapping.map(inSample, outValues);
389
390 BOOST_TEST(!outValues.hasNaN());
391 BOOST_TEST((outValues.array() == 0.0).count() == 0);
392
393 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
394 BOOST_TEST(!inputIntegral.hasNaN());
395 BOOST_TEST(inputIntegral(0) == 2.0);
396 BOOST_TEST(inputIntegral(1) == 0.0);
397
398 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
399 BOOST_TEST(!outputIntegral.hasNaN());
400 BOOST_TEST(outputIntegral(0) == 2.0);
401 BOOST_TEST(outputIntegral(1) == 0.0);
402}
403
404BOOST_AUTO_TEST_CASE(ScaledConsistentVolume2D)
405{
406 PRECICE_TEST(1_rank);
407 int dimensions = 2;
408
409 // Create mesh to map from
410 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
411
412 // One square with 3 triangles
413 Vertex &inVertex0 = inMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
414 Vertex &inVertex1 = inMesh->createVertex(Eigen::Vector2d{1.0, 0.0});
415 Vertex &inVertex2 = inMesh->createVertex(Eigen::Vector2d{1.0, 1.0});
416 Vertex &inVertex3 = inMesh->createVertex(Eigen::Vector2d{0.0, 1.0});
417 Vertex &inVertex4 = inMesh->createVertex(Eigen::Vector2d{0.5, 1.0});
418
419 inMesh->createTriangle(inVertex0, inVertex1, inVertex4);
420 inMesh->createTriangle(inVertex0, inVertex3, inVertex4);
421 inMesh->createTriangle(inVertex1, inVertex2, inVertex4);
422
423 Eigen::VectorXd inValues(5);
424 inValues(0) = 1.0;
425 inValues(1) = 2.0;
426 inValues(2) = 3.0;
427 inValues(3) = 4.0;
428 inValues(4) = 5.0;
429
430 // Create mesh to map to
431 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
432
433 // Unit square as 2 triangles
434 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
435 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector2d(1.0, 0.0));
436 Vertex &outVertex2 = outMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
437 Vertex &outVertex3 = outMesh->createVertex(Eigen::Vector2d(0.0, 1.0));
438
439 outMesh->createTriangle(outVertex0, outVertex1, outVertex2);
440 outMesh->createTriangle(outVertex0, outVertex2, outVertex3);
441
442 // Setup mapping with mapping coordinates and geometry used
444
445 mapping.setMeshes(inMesh, outMesh);
446 BOOST_TEST(mapping.hasComputedMapping() == false);
447
448 Eigen::VectorXd outValues(5);
449 mapping.computeMapping();
450 time::Sample inSample(1, inValues);
451 mapping.map(inSample, outValues);
452
453 BOOST_TEST(mapping.hasComputedMapping() == true);
454
455 auto inputIntegral = mesh::integrateVolume(inMesh, inValues);
456 auto outputIntegral = mesh::integrateVolume(outMesh, outValues);
457
458 Eigen::VectorXd expectedIntegral(1);
459 expectedIntegral << 3.0;
460
461 BOOST_TEST(inputIntegral(0) == expectedIntegral(0));
462
463 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
464 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
465 }
466
467 double scaleFactor = outValues(0) / inValues(0);
468 BOOST_TEST(scaleFactor != 1.0);
469
470 BOOST_TEST(math::equals(inValues(0) * scaleFactor, outValues(0)));
471 BOOST_TEST(inValues(1) * scaleFactor == outValues(1));
472 BOOST_TEST(inValues(2) * scaleFactor == outValues(2));
473 BOOST_TEST(inValues(3) * scaleFactor == outValues(3));
474}
475
476BOOST_AUTO_TEST_CASE(ScaledConsistentVolume3D)
477{
478 PRECICE_TEST(1_rank);
479 int dimensions = 3;
480
481 // Create mesh to map from
482 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
483
484 // One tetra on "out". The "in" has the same but split into two
485 Vertex &inVertex0 = inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
486 Vertex &inVertex1 = inMesh->createVertex(Eigen::Vector3d{0.5, 1.0, 0.0});
487 Vertex &inVertex2 = inMesh->createVertex(Eigen::Vector3d{1.0, 0.0, 0.0});
488 Vertex &inVertex3 = inMesh->createVertex(Eigen::Vector3d{0.5, 0.0, 1.0});
489 Vertex &inVertex4 = inMesh->createVertex(Eigen::Vector3d{0.5, 0.0, 0.0});
490
491 inMesh->createTetrahedron(inVertex0, inVertex1, inVertex3, inVertex4);
492 inMesh->createTetrahedron(inVertex1, inVertex2, inVertex3, inVertex4);
493
494 Eigen::VectorXd inValues(5);
495 inValues(0) = 1.0;
496 inValues(1) = 2.0;
497 inValues(2) = 3.0;
498 inValues(3) = 4.0;
499 inValues(4) = 5.0;
500
501 // Create mesh to map to
502 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
503
504 // One big tetra
505 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
506 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector3d{0.5, 1.0, 0.0});
507 Vertex &outVertex2 = outMesh->createVertex(Eigen::Vector3d{1.0, 0.0, 0.0});
508 Vertex &outVertex3 = outMesh->createVertex(Eigen::Vector3d{0.5, 0.0, 1.0});
509
510 outMesh->createTetrahedron(outVertex0, outVertex1, outVertex2, outVertex3);
511
512 // Setup mapping with mapping coordinates and geometry used
514
515 mapping.setMeshes(inMesh, outMesh);
516 BOOST_TEST(mapping.hasComputedMapping() == false);
517
518 Eigen::VectorXd outValues(5);
519 mapping.computeMapping();
520 time::Sample inSample(1, inValues);
521 mapping.map(inSample, outValues);
522
523 BOOST_TEST(mapping.hasComputedMapping() == true);
524
525 auto inputIntegral = mesh::integrateVolume(inMesh, inValues);
526 auto outputIntegral = mesh::integrateVolume(outMesh, outValues);
527
528 Eigen::VectorXd expectedIntegral(1);
529 expectedIntegral << 6.5 * 1. / 12;
530
531 BOOST_TEST(inputIntegral(0) == expectedIntegral(0));
532
533 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
534 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
535 }
536
537 double scaleFactor = outValues(0) / inValues(0);
538 BOOST_TEST(scaleFactor != 1.0);
539
540 BOOST_TEST(math::equals(inValues(0) * scaleFactor, outValues(0)));
541 BOOST_TEST(inValues(1) * scaleFactor == outValues(1));
542 BOOST_TEST(inValues(2) * scaleFactor == outValues(2));
543 BOOST_TEST(inValues(3) * scaleFactor == outValues(3));
544}
545
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(ConsistentNonIncremental)
std::string name
#define PRECICE_TEST(...)
Definition Testing.hpp:27
T begin(T... args)
void setMeshes(const mesh::PtrMesh &input, const mesh::PtrMesh &output)
Sets input and output meshes carrying data to be mapped.
Definition Mapping.cpp:27
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
void computeMapping() final override
Computes the mapping coefficients from the in- and output mesh.
Mapping using nearest neighboring vertices.
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
Tetrahedron & createTetrahedron(Vertex &vertexOne, Vertex &vertexTwo, Vertex &vertexThree, Vertex &vertexFour)
Creates and initializes a Tetrahedron object.
Definition Mesh.cpp:141
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
Vertex of a mesh.
Definition Vertex.hpp:16
void setCoords(const VECTOR_T &coordinates)
Sets the coordinates of the vertex.
Definition Vertex.hpp:102
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:116
T end(T... args)
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
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
Eigen::VectorXd integrateVolume(const PtrMesh &mesh, const Eigen::VectorXd &input)
Given the data and the mesh, this function returns the volume integral. Assumes no overlap exists for...
Definition Utils.cpp:41
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 push_back(T... args)