preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
NearestProjectionMappingTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <memory>
4#include <ostream>
5#include "mapping/Mapping.hpp"
7#include "math/constants.hpp"
8#include "mesh/Data.hpp"
9#include "mesh/Mesh.hpp"
11#include "mesh/Utils.hpp"
13#include "testing/Testing.hpp"
14#include "utils/assertion.hpp"
15
16namespace precice::mesh {
17class Edge;
18class Vertex;
19} // namespace precice::mesh
20
21using namespace precice;
22
23BOOST_AUTO_TEST_SUITE(MappingTests)
24BOOST_AUTO_TEST_SUITE(NearestProjectionMapping)
25
27BOOST_AUTO_TEST_CASE(testConservativeNonIncremental)
28{
30 using namespace mesh;
31 int dimensions = 2;
32
33 // Setup geometry to map to
34 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
35 Vertex &v1 = outMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
36 Vertex &v2 = outMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
37 outMesh->createEdge(v1, v2);
38 outMesh->allocateDataValues();
39
40 // Base-value for tests
41 double value = 1.0;
42
43 {
44 // Setup mapping with mapping coordinates and geometry used
46 PtrMesh inMesh(new Mesh("InMesh0", dimensions, testing::nextMeshID()));
47 // Map value 1.0 from middle of edge to geometry. Expect half of the
48 // value to be added to vertex1 and half of it to vertex2.
49 inMesh->createVertex(Eigen::Vector2d(0.5, 0.5));
50 // Map value 1.0 from below edge to geometry. Expect vertex1 to get the
51 // full data value, i.e. 1.0 and in addition the value from before. In total
52 // v1 should have 1.5 * dataValue then.
53 inMesh->createVertex(Eigen::Vector2d(-0.5, -0.5));
54 // Do the same thing from above, expect vertex2 to get the full value now.
55 inMesh->createVertex(Eigen::Vector2d(1.5, 1.5));
56
57 Eigen::VectorXd inValues(3);
58 inValues = Eigen::VectorXd::Constant(inValues.size(), value);
59 Eigen::VectorXd values(2);
60 values = Eigen::VectorXd::Constant(values.size(), 0.0);
61
62 mapping.setMeshes(inMesh, outMesh);
63 mapping.computeMapping();
64 time::Sample inSample(1, inValues);
65 mapping.map(inSample, values);
66 BOOST_TEST_CONTEXT(*inMesh)
67 {
68 BOOST_TEST(values(0) == value * 1.5);
69 BOOST_TEST(values(1) == value * 1.5);
70 }
71 }
72 {
73 // Setup mapping with mapping coordinates and geometry used
75 PtrMesh inMesh(new Mesh("InMesh1", dimensions, testing::nextMeshID()));
76
77 inMesh->createVertex(Eigen::Vector2d(-1.0, -1.0));
78 inMesh->createVertex(Eigen::Vector2d(-1.0, -1.0));
79 inMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
80
81 Eigen::VectorXd inValues(3);
82 inValues = Eigen::VectorXd::Constant(inValues.size(), value);
83 Eigen::VectorXd values(2);
84 values = Eigen::VectorXd::Constant(values.size(), 0.0);
85
86 mapping.setMeshes(inMesh, outMesh);
87 mapping.computeMapping();
88 time::Sample inSample(1, inValues);
89 mapping.map(inSample, values);
90 BOOST_TEST_CONTEXT(*inMesh)
91 {
92 BOOST_TEST(values(0) == value * 2.0);
93 BOOST_TEST(values(1) == value * 1.0);
94 }
95
96 // reset output value and remap
97 // assign(values) = 0.0;
98 values = Eigen::VectorXd::Constant(values.size(), 0.0);
99
100 mapping.map(inSample, values);
101 BOOST_TEST_CONTEXT(*inMesh)
102 {
103 BOOST_TEST(values(0) == value * 2.0);
104 BOOST_TEST(values(1) == value * 1.0);
105 }
106 }
107}
108
109PRECICE_TEST_SETUP(1_rank)
110BOOST_AUTO_TEST_CASE(ConsistentNonIncremental2D)
111{
112 PRECICE_TEST();
113 using namespace mesh;
114 int dimensions = 2;
115
116 // Create mesh to map from
117 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
118 Vertex &v1 = inMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
119 Vertex &v2 = inMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
120 inMesh->createEdge(v1, v2);
121 double valueVertex1 = 1.0;
122 double valueVertex2 = 2.0;
123 Eigen::VectorXd inValues(2);
124 inValues(0) = valueVertex1;
125 inValues(1) = valueVertex2;
126
127 {
128 // Create mesh to map to
129 PtrMesh outMesh(new Mesh("OutMesh0", dimensions, testing::nextMeshID()));
130
131 // Setup mapping with mapping coordinates and geometry used
133 mapping.setMeshes(inMesh, outMesh);
134 BOOST_TEST(mapping.hasComputedMapping() == false);
135
136 outMesh->createVertex(Eigen::Vector2d(0.5, 0.5));
137 outMesh->createVertex(Eigen::Vector2d(-0.5, -0.5));
138 outMesh->createVertex(Eigen::Vector2d(1.5, 1.5));
139
140 // Compute and perform mapping
141 mapping.computeMapping();
142 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
143 time::Sample inSample(1, inValues);
144 mapping.map(inSample, outValues);
145
146 // Validate results
147 BOOST_TEST(mapping.hasComputedMapping() == true);
148 BOOST_TEST(outValues(0) == (valueVertex1 + valueVertex2) * 0.5);
149 BOOST_TEST(outValues(1) == valueVertex1);
150 BOOST_TEST(outValues(2) == valueVertex2);
151
152 // Redo mapping, results should be
153 // assign(outData->values()) = 0.0;
154 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
155
156 mapping.map(inSample, outValues);
157 BOOST_TEST(outValues(0) == (valueVertex1 + valueVertex2) * 0.5);
158 BOOST_TEST(outValues(1) == valueVertex1);
159 BOOST_TEST(outValues(2) == valueVertex2);
160 }
161
162 {
163 // Create mesh to map to
164 PtrMesh outMesh(new Mesh("OutMesh1", dimensions, testing::nextMeshID()));
165
166 // Setup mapping with mapping coordinates and geometry used
168 mapping.setMeshes(inMesh, outMesh);
169 BOOST_TEST(mapping.hasComputedMapping() == false);
170
171 outMesh->createVertex(Eigen::Vector2d(-0.5, -0.5));
172 outMesh->createVertex(Eigen::Vector2d(1.5, 1.5));
173 outMesh->createVertex(Eigen::Vector2d(0.5, 0.5));
174
175 // assign(outData->values()) = 0.0;
176 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
177 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
178
179 mapping.computeMapping();
180 time::Sample inSample(1, inValues);
181 mapping.map(inSample, outValues);
182 BOOST_TEST(outValues(0) == valueVertex1);
183 BOOST_TEST(outValues(1) == valueVertex2);
184 BOOST_TEST(outValues(2) == (valueVertex1 + valueVertex2) * 0.5);
185
186 // Reset output data to zero and redo the mapping
187 // assign(outData->values()) = 0.0;
188 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
189
190 mapping.map(inSample, outValues);
191 BOOST_TEST(outValues(0) == valueVertex1);
192 BOOST_TEST(outValues(1) == valueVertex2);
193 BOOST_TEST(outValues(2) == (valueVertex1 + valueVertex2) * 0.5);
194 }
195}
196
197PRECICE_TEST_SETUP(1_rank)
198BOOST_AUTO_TEST_CASE(ScaleConsistentNonIncremental2DCase1)
199{
200 PRECICE_TEST();
201 using namespace mesh;
202 int dimensions = 2;
203
204 // Create mesh to map from
205 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
206 Vertex &v1 = inMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
207 Vertex &v2 = inMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
208 inMesh->createEdge(v1, v2);
209 double valueVertex1 = 1.0;
210 double valueVertex2 = 2.0;
211 Eigen::VectorXd inValues(2);
212 inValues(0) = valueVertex1;
213 inValues(1) = valueVertex2;
214
215 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
216 // Create mesh to map to
217 PtrMesh outMesh(new Mesh("OutMesh0", dimensions, testing::nextMeshID()));
218 Eigen::VectorXd outValues(3);
219 // Setup mapping with mapping coordinates and geometry used
221 mapping.setMeshes(inMesh, outMesh);
222 BOOST_TEST(mapping.hasComputedMapping() == false);
223
224 Vertex &outV1 = outMesh->createVertex(Eigen::Vector2d(0.5, 0.5));
225 Vertex &outV2 = outMesh->createVertex(Eigen::Vector2d(-0.5, -0.5));
226 Vertex &outV3 = outMesh->createVertex(Eigen::Vector2d(1.5, 1.5));
227
228 outMesh->createEdge(outV1, outV2);
229 outMesh->createEdge(outV1, outV3);
230
231 outMesh->allocateDataValues();
232 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
233
234 // Compute and perform mapping
235 mapping.computeMapping();
236
237 time::Sample inSample(1, inValues);
238 mapping.map(inSample, outValues);
239
240 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
241 double scaleFactor = outValues(1) / inValues(0);
242 BOOST_TEST(scaleFactor != 1.0);
243
244 // Validate results
245 BOOST_TEST(mapping.hasComputedMapping() == true);
246 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
247 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
248 }
249 BOOST_TEST(outValues(0) == (inValues(0) + inValues(1)) * 0.5 * scaleFactor);
250 BOOST_TEST(outValues(1) == inValues(0) * scaleFactor);
251 BOOST_TEST(outValues(2) == inValues(1) * scaleFactor);
252}
253
254PRECICE_TEST_SETUP(1_rank)
255BOOST_AUTO_TEST_CASE(ScaleConsistentNonIncremental2DCase2)
256{
257 PRECICE_TEST();
258 using namespace mesh;
259 int dimensions = 2;
260
261 // Create mesh to map from
262 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
263 Vertex &v1 = inMesh->createVertex(Eigen::Vector2d(0.0, 0.0));
264 Vertex &v2 = inMesh->createVertex(Eigen::Vector2d(1.0, 1.0));
265 inMesh->createEdge(v1, v2);
266 double valueVertex1 = 1.0;
267 double valueVertex2 = 2.0;
268 Eigen::VectorXd inValues(2);
269 inValues(0) = valueVertex1;
270 inValues(1) = valueVertex2;
271
272 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
273
274 // Create mesh to map to
275 PtrMesh outMesh(new Mesh("OutMesh1", dimensions, testing::nextMeshID()));
276 Eigen::VectorXd outValues(3);
277 // Setup mapping with mapping coordinates and geometry used
279 mapping.setMeshes(inMesh, outMesh);
280 BOOST_TEST(mapping.hasComputedMapping() == false);
281
282 Vertex &outV1 = outMesh->createVertex(Eigen::Vector2d(-0.5, -0.5));
283 Vertex &outV2 = outMesh->createVertex(Eigen::Vector2d(1.5, 1.5));
284 Vertex &outV3 = outMesh->createVertex(Eigen::Vector2d(0.5, 0.5));
285
286 outMesh->createEdge(outV3, outV1);
287 outMesh->createEdge(outV3, outV2);
288
289 outValues = Eigen::VectorXd::Constant(outValues.size(), 0.0);
290
291 mapping.computeMapping();
292 time::Sample inSample(1, inValues);
293 mapping.map(inSample, outValues);
294
295 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
296 double scaleFactor = outValues(0) / inValues(0);
297 BOOST_TEST(scaleFactor != 1.0);
298
299 // Validate results
300 BOOST_TEST(mapping.hasComputedMapping() == true);
301 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
302 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
303 }
304 BOOST_TEST(outValues(0) == inValues(0) * scaleFactor);
305 BOOST_TEST(outValues(1) == inValues(1) * scaleFactor);
306 BOOST_TEST(outValues(2) == (inValues(0) + inValues(1)) * 0.5 * scaleFactor);
307}
308
309PRECICE_TEST_SETUP(1_rank)
310BOOST_AUTO_TEST_CASE(Consistent3DFalbackOnEdges)
311{
312 PRECICE_TEST();
313 using namespace mesh;
314 int dimensions = 3;
315
316 // Create mesh to map from
317 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
318 Vertex &v1 = inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
319 Vertex &v2 = inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
320 Vertex &v3 = inMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
321 inMesh->createEdge(v1, v2);
322 inMesh->createEdge(v2, v3);
323 inMesh->createEdge(v3, v1);
324
325 double valueVertex1 = 1.0;
326 double valueVertex2 = 2.0;
327 double valueVertex3 = 3.0;
328 Eigen::VectorXd values(3);
329 values(0) = valueVertex1;
330 values(1) = valueVertex2;
331 values(2) = valueVertex3;
332
333 // Create mesh to map to
334 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
335
336 // Setup mapping with mapping coordinates and geometry used
338 mapping.setMeshes(inMesh, outMesh);
339 BOOST_TEST(mapping.hasComputedMapping() == false);
340
341 outMesh->createVertex(Eigen::Vector3d(0.0, 0.5, 0.0));
342 outMesh->createVertex(Eigen::Vector3d(0.5, 0.0, 0.0));
343 outMesh->createVertex(Eigen::Vector3d(0.5, 0.5, 0.0));
344
345 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
346 // Compute and perform mapping
347 mapping.computeMapping();
348 time::Sample inSample(1, values);
349 mapping.map(inSample, outValues);
350
351 // Validate results
352 BOOST_TEST(mapping.hasComputedMapping() == true);
353 BOOST_TEST_CONTEXT(*inMesh)
354 {
355 BOOST_TEST(outValues(0) == (valueVertex1 + valueVertex2) * 0.5);
356 BOOST_TEST(outValues(1) == (valueVertex1 + valueVertex3) * 0.5);
357 BOOST_TEST(outValues(2) == (valueVertex2 + valueVertex3) * 0.5);
358 }
359}
360
361PRECICE_TEST_SETUP(1_rank)
362BOOST_AUTO_TEST_CASE(Consistent3DFalbackOnVertices)
363{
364 PRECICE_TEST();
365 using namespace mesh;
366 int dimensions = 3;
367
368 // Create mesh to map from
369 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
370 inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
371 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
372 inMesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
373
374 double valueVertex1 = 1.0;
375 double valueVertex2 = 2.0;
376 double valueVertex3 = 3.0;
377 Eigen::VectorXd values(3);
378 values(0) = valueVertex1;
379 values(1) = valueVertex2;
380 values(2) = valueVertex3;
381
382 // Create mesh to map to
383 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
384
385 // Setup mapping with mapping coordinates and geometry used
387 mapping.setMeshes(inMesh, outMesh);
388 BOOST_TEST(mapping.hasComputedMapping() == false);
389
390 outMesh->createVertex(Eigen::Vector3d(0.1, 0.1, 0.0));
391 outMesh->createVertex(Eigen::Vector3d(0.1, 1.1, 0.0));
392 outMesh->createVertex(Eigen::Vector3d(0.9, 0.1, 0.0));
393
394 // Compute and perform mapping
395 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(3);
396 mapping.computeMapping();
397 time::Sample inSample(1, values);
398 mapping.map(inSample, outValues);
399
400 // Validate results
401 BOOST_TEST(mapping.hasComputedMapping() == true);
402 BOOST_TEST_CONTEXT(*inMesh)
403 {
404 BOOST_TEST(outValues(0) == valueVertex1);
405 BOOST_TEST(outValues(1) == valueVertex2);
406 BOOST_TEST(outValues(2) == valueVertex3);
407 }
408}
409
410PRECICE_TEST_SETUP(1_rank)
411BOOST_AUTO_TEST_CASE(AxisAlignedTriangles)
412{
413 PRECICE_TEST();
414 using namespace precice::mesh;
415 constexpr int dimensions = 3;
416
417 // Create mesh to map from with Triangles A-B-D and B-D-C
418 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
419 Vertex &inVA = inMesh->createVertex(Eigen::Vector3d{0, 0, 0});
420 Vertex &inVB = inMesh->createVertex(Eigen::Vector3d{0, 1, 0});
421 Vertex &inVC = inMesh->createVertex(Eigen::Vector3d{1, 1, 0});
422 Vertex &inVD = inMesh->createVertex(Eigen::Vector3d{1, 0, 0});
423
424 Edge &inEDA = inMesh->createEdge(inVD, inVA);
425 Edge &inEAB = inMesh->createEdge(inVA, inVB);
426 Edge &inEBD = inMesh->createEdge(inVB, inVD);
427 Edge &inEDC = inMesh->createEdge(inVD, inVC);
428 Edge &inECB = inMesh->createEdge(inVC, inVB);
429
430 inMesh->createTriangle(inEAB, inEBD, inEDA);
431 inMesh->createTriangle(inEBD, inEDC, inECB);
432 Eigen::VectorXd inValues = Eigen::VectorXd::Zero(4);
433 inValues << 1.0, 1.0, 1.0, 1.0;
434
435 // Create mesh to map to with one vertex per defined triangle
436 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
437 outMesh->createVertex(Eigen::Vector3d{0.33, 0.33, 0});
438 outMesh->createVertex(Eigen::Vector3d{0.66, 0.66, 0});
439 Eigen::VectorXd outValues = Eigen::VectorXd::Zero(2);
440 outValues << 0.0, 0.0;
441
442 // Setup mapping with mapping coordinates and geometry used
444 mapping.setMeshes(inMesh, outMesh);
445 BOOST_TEST(mapping.hasComputedMapping() == false);
446
447 mapping.computeMapping();
448 BOOST_TEST(mapping.hasComputedMapping() == true);
449 BOOST_TEST_INFO("In Data:" << inValues);
450 BOOST_TEST_INFO("Out Data before Mapping:" << outValues);
451 time::Sample inSample(1, inValues);
452 mapping.map(inSample, outValues);
453 BOOST_TEST_INFO("Out Data after Mapping:" << outValues);
454 BOOST_TEST(outValues == outValues.cwiseAbs());
455}
456
457PRECICE_TEST_SETUP(1_rank)
458BOOST_AUTO_TEST_CASE(Query_3D_FullMesh)
459{
460 PRECICE_TEST();
461 using namespace precice::mesh;
462 constexpr int dimensions = 3;
463
464 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
465 // PtrData inData = inMesh->createData("InData", 1, 0_dataID);
466 const double z1 = 0.1;
467 const double z2 = -0.1;
468 auto &v00 = inMesh->createVertex(Eigen::Vector3d(0, 0, 0));
469 auto &v01 = inMesh->createVertex(Eigen::Vector3d(0, 1, 0));
470 auto &v10 = inMesh->createVertex(Eigen::Vector3d(1, 0, z1));
471 auto &v11 = inMesh->createVertex(Eigen::Vector3d(1, 1, z1));
472 auto &v20 = inMesh->createVertex(Eigen::Vector3d(2, 0, z2));
473 auto &v21 = inMesh->createVertex(Eigen::Vector3d(2, 1, z2));
474 auto &ell = inMesh->createEdge(v00, v01);
475 auto &elt = inMesh->createEdge(v01, v11);
476 auto &elr = inMesh->createEdge(v11, v10);
477 auto &elb = inMesh->createEdge(v10, v00);
478 auto &eld = inMesh->createEdge(v00, v11);
479 auto &erl = elr;
480 auto &ert = inMesh->createEdge(v11, v21);
481 auto &err = inMesh->createEdge(v21, v20);
482 auto &erb = inMesh->createEdge(v20, v10);
483 auto &erd = inMesh->createEdge(v10, v21);
484 inMesh->createTriangle(ell, elt, eld);
485 inMesh->createTriangle(eld, elb, elr);
486 inMesh->createTriangle(erl, ert, erd);
487 inMesh->createTriangle(erd, erb, err);
488
489 Eigen::VectorXd inValues = Eigen::VectorXd::Constant(6, 1.0);
490
491 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
492 outMesh->createVertex(Eigen::Vector3d{0.7, 0.5, 0.0});
493
494 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(1, 0.0);
495
496 // Setup mapping with mapping coordinates and geometry used
498 mapping.setMeshes(inMesh, outMesh);
499 BOOST_TEST(mapping.hasComputedMapping() == false);
500
501 mapping.computeMapping();
502 BOOST_TEST(mapping.hasComputedMapping() == true);
503
504 BOOST_TEST_INFO("In Data:" << inValues);
505 BOOST_TEST_INFO("Out Data before Mapping:" << outValues);
506 time::Sample inSample(1, inValues);
507 mapping.map(inSample, outValues);
508 BOOST_TEST_INFO("Out Data after Mapping:" << outValues);
509 BOOST_TEST(outValues(0) == 1.0);
510}
511
512PRECICE_TEST_SETUP(1_rank)
513BOOST_AUTO_TEST_CASE(ScaledConsistentQuery3DFullMesh)
514{
515 PRECICE_TEST();
516 using namespace precice::mesh;
517 constexpr int dimensions = 3;
518
519 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
520 const double z1 = 0.1;
521 const double z2 = -0.1;
522 auto &v00 = inMesh->createVertex(Eigen::Vector3d(0, 0, 0));
523 auto &v01 = inMesh->createVertex(Eigen::Vector3d(0, 1, 0));
524 auto &v10 = inMesh->createVertex(Eigen::Vector3d(1, 0, z1));
525 auto &v11 = inMesh->createVertex(Eigen::Vector3d(1, 1, z1));
526 auto &v20 = inMesh->createVertex(Eigen::Vector3d(2, 0, z2));
527 auto &v21 = inMesh->createVertex(Eigen::Vector3d(2, 1, z2));
528 auto &ell = inMesh->createEdge(v00, v01);
529 auto &elt = inMesh->createEdge(v01, v11);
530 auto &elr = inMesh->createEdge(v11, v10);
531 auto &elb = inMesh->createEdge(v10, v00);
532 auto &eld = inMesh->createEdge(v00, v11);
533 auto &erl = elr;
534 auto &ert = inMesh->createEdge(v11, v21);
535 auto &err = inMesh->createEdge(v21, v20);
536 auto &erb = inMesh->createEdge(v20, v10);
537 auto &erd = inMesh->createEdge(v10, v21);
538 inMesh->createTriangle(ell, elt, eld);
539 inMesh->createTriangle(eld, elb, elr);
540 inMesh->createTriangle(erl, ert, erd);
541 inMesh->createTriangle(erd, erb, err);
542
543 Eigen::VectorXd inValues = Eigen::VectorXd::Constant(6, 1.0);
544
545 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
546 auto &outV1 = outMesh->createVertex(Eigen::Vector3d{0.7, 0.5, 0.0});
547 auto &outV2 = outMesh->createVertex(Eigen::Vector3d{0.5, 0.0, 0.05});
548 auto &outV3 = outMesh->createVertex(Eigen::Vector3d{0.5, 0.0, 0.0});
549 auto &outE1 = outMesh->createEdge(outV1, outV2);
550 auto &outE2 = outMesh->createEdge(outV2, outV3);
551 auto &outE3 = outMesh->createEdge(outV1, outV3);
552 outMesh->createTriangle(outE1, outE2, outE3);
553 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(3, 0.0);
554
555 // Setup mapping with mapping coordinates and geometry used
557 mapping.setMeshes(inMesh, outMesh);
558 BOOST_TEST(mapping.hasComputedMapping() == false);
559
560 mapping.computeMapping();
561 BOOST_TEST(mapping.hasComputedMapping() == true);
562
563 time::Sample inSample(1, inValues);
564 mapping.map(inSample, outValues);
565
566 auto inputIntegral = mesh::integrateSurface(inMesh, inValues);
567 auto outputIntegral = mesh::integrateSurface(outMesh, outValues);
568
569 for (int dim = 0; dim < inputIntegral.size(); ++dim) {
570 BOOST_TEST(inputIntegral(dim) == outputIntegral(dim));
571 }
572}
573
574namespace {
575using namespace precice::mesh;
576const Eigen::VectorXd &runNPMapping(mapping::Mapping::Constraint constraint, PtrMesh &inMesh, Eigen::VectorXd *inData, PtrMesh &outMesh, Eigen::VectorXd *outData)
577{
578 BOOST_REQUIRE(inMesh->getDimensions() == outMesh->getDimensions());
579 precice::mapping::NearestProjectionMapping mapping(constraint, inMesh->getDimensions());
580 mapping.setMeshes(inMesh, outMesh);
581 BOOST_REQUIRE(mapping.hasComputedMapping() == false);
582 mapping.computeMapping();
583 BOOST_REQUIRE(mapping.hasComputedMapping() == true);
584 time::Sample inSample(1, *inData);
585 mapping.map(inSample, *outData);
586 return *outData;
587}
588
589void makeTriangle(PtrMesh &inMesh, Vertex &a, Vertex &b, Vertex &c)
590{
591 auto &ab = inMesh->createEdge(a, b);
592 auto &bc = inMesh->createEdge(b, c);
593 auto &ca = inMesh->createEdge(c, a);
594 inMesh->createTriangle(ab, bc, ca);
595}
596} // namespace
597
598PRECICE_TEST_SETUP(1_rank)
599BOOST_AUTO_TEST_CASE(AvoidClosestTriangle)
600{
601 PRECICE_TEST();
602 using namespace precice::mesh;
603 constexpr int dimensions = 3;
604
605 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
606 // Close triangle - extrapolating
607 auto &vc0 = inMesh->createVertex(Eigen::Vector3d(3, 0, 0));
608 auto &vc1 = inMesh->createVertex(Eigen::Vector3d(3, 2, 0));
609 auto &vc2 = inMesh->createVertex(Eigen::Vector3d(4, 1, 0));
610 makeTriangle(inMesh, vc0, vc1, vc2);
611
612 // Far triangle - interpolating
613 auto &vf0 = inMesh->createVertex(Eigen::Vector3d(0, 0, -1));
614 auto &vf1 = inMesh->createVertex(Eigen::Vector3d(0, 2, -1));
615 auto &vf2 = inMesh->createVertex(Eigen::Vector3d(0, 1, 1));
616 makeTriangle(inMesh, vf0, vf1, vf2);
617
618 Eigen::VectorXd inValues(5);
619 inValues << 0, 0, 0, 1, 1;
620
621 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
622 outMesh->createVertex(Eigen::Vector3d{2, 1, 0});
623 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(1, 0.0);
624
625 const auto &values = runNPMapping(mapping::Mapping::CONSISTENT, inMesh, &inValues, outMesh, &outValues);
626
627 // Interpolating triangle is further than NN => fall back on NN
628 BOOST_TEST(values(0) == 0.0);
629}
630
631PRECICE_TEST_SETUP(1_rank)
632BOOST_AUTO_TEST_CASE(PickClosestTriangle)
633{
634 PRECICE_TEST();
635 using namespace precice::mesh;
636
637 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
638 // Far triangle - interpolating
639 auto &vf0 = inMesh->createVertex(Eigen::Vector3d(0, 0, -1));
640 auto &vf1 = inMesh->createVertex(Eigen::Vector3d(0, 1, 1));
641 auto &vf2 = inMesh->createVertex(Eigen::Vector3d(0, 2, -1));
642 makeTriangle(inMesh, vf0, vf1, vf2);
643
644 // Close triangle - extrapolating
645 auto &vc0 = inMesh->createVertex(Eigen::Vector3d(3, 0, 0));
646 auto &vc1 = inMesh->createVertex(Eigen::Vector3d(3, 2, 0));
647 auto &vc2 = inMesh->createVertex(Eigen::Vector3d(4, 1, 0));
648 makeTriangle(inMesh, vc0, vc1, vc2);
649
650 Eigen::VectorXd inValues(6);
651 inValues << 1, 1, 1, 0, 0, 0;
652
653 PtrMesh outMesh(new Mesh("OutMesh", 3, testing::nextMeshID()));
654 outMesh->createVertex(Eigen::Vector3d{1, 1, 0});
655 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(1, 0.0);
656
657 const auto &values = runNPMapping(mapping::Mapping::CONSISTENT, inMesh, &inValues, outMesh, &outValues);
658
659 BOOST_TEST(values(0) == 1.0);
660}
661
662PRECICE_TEST_SETUP(1_rank)
663BOOST_AUTO_TEST_CASE(PreferTriangleOverEdge)
664{
665 PRECICE_TEST();
666 using namespace precice::mesh;
667 constexpr int dimensions = 3;
668
669 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
670 // Close edge ->
671 auto &vc0 = inMesh->createVertex(Eigen::Vector3d(0, 0, 0));
672 auto &vc1 = inMesh->createVertex(Eigen::Vector3d(0, 2, 0));
673 inMesh->createEdge(vc0, vc1);
674
675 // Far triangle - interpolating
676 auto &vf0 = inMesh->createVertex(Eigen::Vector3d(3, 0, 0));
677 auto &vf1 = inMesh->createVertex(Eigen::Vector3d(3, 2, 2));
678 auto &vf2 = inMesh->createVertex(Eigen::Vector3d(3, 1, 0));
679 makeTriangle(inMesh, vf0, vf1, vf2);
680
681 Eigen::VectorXd inValues(5);
682 inValues << 0, 1, 2, 2, 2;
683
684 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
685 outMesh->createVertex(Eigen::Vector3d{1, 1, 1});
686 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(1, 0.0);
687
688 const auto &values = runNPMapping(mapping::Mapping::CONSISTENT, inMesh, &inValues, outMesh, &outValues);
689
690 // Distance to triangle > distance to NN > distance to edge => Interpolation on edge
691 // Projection is on the middle of the edge => (0+1)/2 = 0.5
692 BOOST_TEST(values(0) == 0.5);
693}
694
695PRECICE_TEST_SETUP(1_rank)
696BOOST_AUTO_TEST_CASE(TriangleDistances)
697{
698 PRECICE_TEST();
699 using namespace precice::mesh;
700 constexpr int dimensions = 3;
701
702 PtrMesh inMesh(new mesh::Mesh("InMesh", 3, testing::nextMeshID()));
703
704 // Close triangle
705 auto &vc0 = inMesh->createVertex(Eigen::Vector3d(0, 0, 0));
706 auto &vc1 = inMesh->createVertex(Eigen::Vector3d(0, 2, 2));
707 auto &vc2 = inMesh->createVertex(Eigen::Vector3d(0, 1, 0));
708 makeTriangle(inMesh, vc0, vc1, vc2);
709
710 // Far triangle
711 auto &vf0 = inMesh->createVertex(Eigen::Vector3d(3, 0, 0));
712 auto &vf1 = inMesh->createVertex(Eigen::Vector3d(3, 2, 2));
713 auto &vf2 = inMesh->createVertex(Eigen::Vector3d(3, 1, 0));
714 makeTriangle(inMesh, vf0, vf1, vf2);
715
716 Eigen::VectorXd inValues(6);
717 inValues << 1, 1, 1, 0, 0, 0;
718
719 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
720 outMesh->createVertex(Eigen::Vector3d{1, 1, 1});
721 Eigen::VectorXd outValues = Eigen::VectorXd::Constant(1, 0.0);
722
723 const auto &values = runNPMapping(mapping::Mapping::CONSISTENT, inMesh, &inValues, outMesh, &outValues);
724
725 BOOST_TEST(values(0) == 1.0);
726}
727
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(testConservativeNonIncremental)
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
void setMeshes(const mesh::PtrMesh &input, const mesh::PtrMesh &output)
Sets input and output meshes carrying data to be mapped.
Definition Mapping.cpp:28
bool hasComputedMapping() const
Returns true, if the mapping has been computed.
Definition Mapping.cpp:253
void map(int inputDataID, int outputDataID)
Definition Mapping.cpp:127
Mapping using orthogonal projection to nearest triangle/edge/vertex and linear interpolation from pro...
void computeMapping() final override
Computes the projections and interpolation relations.
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:15
Container and creator for meshes.
Definition Mesh.hpp:38
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:121
int getDimensions() const
Definition Mesh.cpp:100
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:113
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:235
Vertex of a mesh.
Definition Vertex.hpp:16
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
Main namespace of the precice library.