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