preCICE v3.1.1
Loading...
Searching...
No Matches
GeometryTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include "../geometry.hpp"
4#include "logging/Logger.hpp"
7#include "testing/Testing.hpp"
8#include "utils/algorithm.hpp"
9
10using namespace precice;
11using namespace precice::math;
12
13BOOST_AUTO_TEST_SUITE(MathTests)
15
17{
18 PRECICE_TEST(1_rank);
19 // 2D test setup
20 Eigen::Vector2d a2D(0, 0);
21 Eigen::Vector2d b2D(1, 1);
22 Eigen::Vector2d collinearPoint2D(0.5, 0.5);
23 Eigen::Vector2d notCollinearPoint2D(0.5, 0.6);
24
25 // 3D test setup
26 Eigen::Vector3d a3D(0, 0, 0);
27 Eigen::Vector3d b3D(1, 1, 1);
28 Eigen::Vector3d collinearPoint3D(0.5, 0.5, 0.5);
29 Eigen::Vector3d notCollinearPoint3D(0.5, 0.6, 0.5);
30
31 // 2D test validations
32 BOOST_CHECK(geometry::collinear(a2D, b2D, collinearPoint2D));
33 BOOST_CHECK(!geometry::collinear(a2D, b2D, notCollinearPoint2D));
34
35 // 3D test validations
36 BOOST_CHECK(geometry::collinear(a3D, b3D, collinearPoint3D));
37 BOOST_CHECK(!geometry::collinear(a3D, b3D, notCollinearPoint3D));
38}
39
41 *boost::unit_test::tolerance(1e-3))
42{
43 PRECICE_TEST(1_rank);
44 Eigen::Vector3d a(1, 2, 3);
45 Eigen::Vector3d b(3, 2, 1);
46 Eigen::Vector3d c(4, 5, 6);
47 Eigen::Vector3d d(6, 5, 4);
48 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 0);
49
50 a << 5, 0, 0;
51 b << 0, -3, -3;
52 c << 0, 3, 4;
53 d << -1, -2, 6;
54 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 38.6666);
55
56 d << -1.47, -4.1, 8.3;
57 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 62.1816);
58}
59
61{
62 PRECICE_TEST(1_rank);
63 for (int dim = 2; dim <= 3; dim++) {
64 Eigen::VectorXd a(dim);
65 Eigen::VectorXd b(dim);
66 Eigen::VectorXd betweenPoint(dim);
67 Eigen::VectorXd betweenLimitPoint(dim);
68 Eigen::VectorXd collinearOutsidePoint(dim);
69 Eigen::VectorXd outsidePoint(dim);
70 if (dim == 2) {
71 a << 0.0, 0.0;
72 b << 1.0, 1.0;
73 betweenPoint << 0.5, 0.5;
74 betweenLimitPoint << 1.0, 1.0;
75 collinearOutsidePoint << 2.0, 2.0;
76 outsidePoint << 0.5, 0.4;
77 } else {
78 a << 0.0, 0.0, 0.0;
79 b << 1.0, 1.0, 1.0;
80 betweenPoint << 0.5, 0.5, 0.5;
81 betweenLimitPoint << 1.0, 1.0, 1.0;
82 collinearOutsidePoint << 2.0, 2.0, 2.0;
83 outsidePoint << 0.5, 0.4, 0.5;
84 }
85 BOOST_CHECK(geometry::between(a, b, betweenPoint));
86 BOOST_CHECK(geometry::between(a, b, betweenLimitPoint));
87 BOOST_CHECK(!geometry::between(a, b, collinearOutsidePoint));
88 BOOST_CHECK(!geometry::between(a, b, outsidePoint));
89 }
90}
91
93{
94 PRECICE_TEST(1_rank);
95 { // 2D
96 Eigen::Vector2d a, b, c;
97 double area;
98 a << 0.0, 0.0;
99 b << 1.0, 0.0;
100 c << 0.0, 1.0;
101 area = geometry::triangleArea(a, b, c);
102 BOOST_TEST(area == 0.5);
103
104 b << 0.0, 1.0;
105 c << 1.0, 0.0;
106 area = geometry::triangleArea(a, b, c);
107 // Areas are not signed.
108 BOOST_TEST(area == 0.5);
109 }
110 { // 3D
111 Eigen::Vector3d a, b, c;
112 double area;
113 a << 0.0, 0.0, 0.0;
114 b << 1.0, 0.0, 1.0;
115 c << 1.0, 1.0, 1.0;
116 area = geometry::triangleArea(a, b, c);
117 BOOST_CHECK(area == 0.5 * sqrt(2.0));
118 }
119}
120
121BOOST_AUTO_TEST_CASE(SegmentPlaneIntersection)
122{
123 PRECICE_TEST(1_rank);
124 using Eigen::Vector3d;
125 Vector3d planeNormal = Vector3d::Constant(1.0);
126 Vector3d pointOnPlane = Vector3d::Constant(0.0);
127 Vector3d firstSegmentPoint = Vector3d::Constant(1.0);
128 Vector3d secondSegmentPoint = Vector3d::Constant(-1.0);
129 Vector3d intersectionPoint = Vector3d::Constant(1.0);
130 Vector3d expected = Vector3d::Constant(0.0);
131
132 // True intersection
134 pointOnPlane, planeNormal, firstSegmentPoint,
135 secondSegmentPoint, intersectionPoint);
136 BOOST_TEST(result == geometry::INTERSECTION);
137 BOOST_CHECK(equals(intersectionPoint, expected));
138
139 // Touching second segment vertex
140 secondSegmentPoint = Vector3d::Constant(0.0);
142 pointOnPlane, planeNormal, firstSegmentPoint,
143 secondSegmentPoint, intersectionPoint);
144 BOOST_TEST(result == geometry::TOUCHING);
145 BOOST_CHECK(equals(intersectionPoint, expected));
146
147 // Touching first segment vertex
148 firstSegmentPoint = Vector3d::Constant(0.0);
149 secondSegmentPoint = Vector3d::Constant(-1.0);
151 pointOnPlane, planeNormal, firstSegmentPoint,
152 secondSegmentPoint, intersectionPoint);
153 BOOST_TEST(result == geometry::TOUCHING);
154 BOOST_CHECK(equals(intersectionPoint, expected));
155
156 // Parallel segment with distance to plain
157 firstSegmentPoint << 0.0, 0.0, -3.0;
158 intersectionPoint << 1.0, 2.0, 3.0; // should not be modified
159 expected = intersectionPoint;
161 pointOnPlane, planeNormal, firstSegmentPoint,
162 secondSegmentPoint, intersectionPoint);
163 BOOST_TEST(result == geometry::NO_INTERSECTION);
164 BOOST_CHECK(equals(intersectionPoint, expected));
165
166 // Parallel segment contained in plane
167 firstSegmentPoint << 0.0, 0.0, 0.0;
168 secondSegmentPoint << 1.0, 1.0, -2.0;
170 pointOnPlane, planeNormal, firstSegmentPoint,
171 secondSegmentPoint, intersectionPoint);
172 BOOST_TEST(result == geometry::CONTAINED);
173 BOOST_CHECK(equals(intersectionPoint, expected));
174
175 // Segment ending before intersection
176 firstSegmentPoint << -2.0, -2.0, -2.0;
177 secondSegmentPoint << -1.0, -1.0, -1.0;
179 pointOnPlane, planeNormal, firstSegmentPoint,
180 secondSegmentPoint, intersectionPoint);
181 BOOST_TEST(result == geometry::NO_INTERSECTION);
182 BOOST_CHECK(equals(intersectionPoint, expected));
183
184 // Segment ending before intersection (inversed segment points)
185 firstSegmentPoint << -1.0, -1.0, -1.0;
186 secondSegmentPoint << -2.0, -2.0, -2.0;
188 pointOnPlane, planeNormal, firstSegmentPoint,
189 secondSegmentPoint, intersectionPoint);
190 BOOST_TEST(result == geometry::NO_INTERSECTION);
191 BOOST_CHECK(equals(intersectionPoint, expected));
192}
193
195{
196 PRECICE_TEST(1_rank);
197 Eigen::Vector3d vector3D(1.0, 2.0, 3.0);
198 Eigen::Vector2d vector2D;
199 Eigen::Vector2d vectorExpected(1.0, 2.0);
200
201 vector2D = geometry::projectVector(vector3D, 2);
202 BOOST_CHECK(equals(vector2D, vectorExpected));
203
204 vector2D = geometry::projectVector(vector3D, 1);
205 vectorExpected << 1.0, 3.0;
206 BOOST_CHECK(equals(vector2D, vectorExpected));
207
208 vector2D = geometry::projectVector(vector3D, 0);
209 vectorExpected << 2.0, 3.0;
210 BOOST_CHECK(equals(vector2D, vectorExpected));
211}
212
213BOOST_AUTO_TEST_CASE(ContainedInHyperrectangle)
214{
215 PRECICE_TEST(1_rank);
216 // 2D
217 Eigen::Vector2d center2D(0, 0);
218 Eigen::Vector2d sidelengths2D(1, 1);
219
220 // Not contained 2D
221 Eigen::Vector2d testPoint2D(2, 2);
223 sidelengths2D, center2D, testPoint2D);
224 BOOST_TEST(result == geometry::NOT_CONTAINED);
225
226 testPoint2D << -2.0, -2.0;
228 sidelengths2D, center2D, testPoint2D);
229 BOOST_TEST(result == geometry::NOT_CONTAINED);
230
231 testPoint2D << 2.0, -2.0;
233 sidelengths2D, center2D, testPoint2D);
234 BOOST_TEST(result == geometry::NOT_CONTAINED);
235
236 testPoint2D << -2.0, 2.0;
238 sidelengths2D, center2D, testPoint2D);
239 BOOST_TEST(result == geometry::NOT_CONTAINED);
240
241 testPoint2D << 2.0, 0.0;
243 sidelengths2D, center2D, testPoint2D);
244 BOOST_TEST(result == geometry::NOT_CONTAINED);
245
246 testPoint2D << -2.0, 0.0;
248 sidelengths2D, center2D, testPoint2D);
249 BOOST_TEST(result == geometry::NOT_CONTAINED);
250
251 testPoint2D << 0.0, 2.0;
253 sidelengths2D, center2D, testPoint2D);
254 BOOST_TEST(result == geometry::NOT_CONTAINED);
255
256 testPoint2D << 0.0, -2.0;
258 sidelengths2D, center2D, testPoint2D);
259 BOOST_TEST(result == geometry::NOT_CONTAINED);
260
261 // Contained 2D
262 testPoint2D << 0.0, 0.0;
264 sidelengths2D, center2D, testPoint2D);
265 BOOST_TEST(result == geometry::CONTAINED);
266
267 testPoint2D << 0.25, 0.25;
269 sidelengths2D, center2D, testPoint2D);
270 BOOST_TEST(result == geometry::CONTAINED);
271
272 testPoint2D << -0.25, 0.25;
274 sidelengths2D, center2D, testPoint2D);
275 BOOST_TEST(result == geometry::CONTAINED);
276
277 testPoint2D << 0.25, -0.25;
279 sidelengths2D, center2D, testPoint2D);
280 BOOST_TEST(result == geometry::CONTAINED);
281
282 testPoint2D << -0.25, -0.25;
284 sidelengths2D, center2D, testPoint2D);
285 BOOST_TEST(result == geometry::CONTAINED);
286
287 testPoint2D << 0.49999999999, 0.49999999999;
289 sidelengths2D, center2D, testPoint2D);
290 BOOST_TEST(result == geometry::CONTAINED);
291
292 testPoint2D << -0.49999999999, -0.49999999999;
294 sidelengths2D, center2D, testPoint2D);
295 BOOST_TEST(result == geometry::CONTAINED);
296
297 // Touching 2D
298 testPoint2D << 0.5, 0.5;
300 sidelengths2D, center2D, testPoint2D);
301 BOOST_TEST(result == geometry::TOUCHING);
302
303 testPoint2D << -0.5, 0.5;
305 sidelengths2D, center2D, testPoint2D);
306 BOOST_TEST(result == geometry::TOUCHING);
307
308 testPoint2D << 0.5, -0.5;
310 sidelengths2D, center2D, testPoint2D);
311 BOOST_TEST(result == geometry::TOUCHING);
312
313 testPoint2D << -0.5, -0.5;
315 sidelengths2D, center2D, testPoint2D);
316 BOOST_TEST(result == geometry::TOUCHING);
317
318 testPoint2D << 0.4999999999999999, 0.4999999999999999;
320 sidelengths2D, center2D, testPoint2D);
321 BOOST_TEST(result == geometry::TOUCHING);
322
323 testPoint2D << 0.500000000000001, 0.500000000000001;
325 sidelengths2D, center2D, testPoint2D);
326 BOOST_TEST(result == geometry::TOUCHING);
327
328 testPoint2D << -0.500000000000001, -0.500000000000001;
330 sidelengths2D, center2D, testPoint2D);
331 BOOST_TEST(result == geometry::TOUCHING);
332
333 testPoint2D << 0.5, 0.0;
335 sidelengths2D, center2D, testPoint2D);
336 BOOST_TEST(result == geometry::TOUCHING);
337
338 testPoint2D << -0.5, 0.0;
340 sidelengths2D, center2D, testPoint2D);
341 BOOST_TEST(result == geometry::TOUCHING);
342
343 testPoint2D << 0.0, 0.5;
345 sidelengths2D, center2D, testPoint2D);
346 BOOST_TEST(result == geometry::TOUCHING);
347
348 testPoint2D << 0.0, -0.5;
350 sidelengths2D, center2D, testPoint2D);
351 BOOST_TEST(result == geometry::TOUCHING);
352
353 testPoint2D << 0.5, 0.0;
355 sidelengths2D, center2D, testPoint2D);
356 BOOST_TEST(result == geometry::TOUCHING);
357
358 testPoint2D << -0.5, 0.0;
360 sidelengths2D, center2D, testPoint2D);
361 BOOST_TEST(result == geometry::TOUCHING);
362
363 testPoint2D << 0.0, 0.499999999999999;
365 sidelengths2D, center2D, testPoint2D);
366 BOOST_TEST(result == geometry::TOUCHING);
367
368 testPoint2D << 0.0, -0.499999999999999;
370 sidelengths2D, center2D, testPoint2D);
371 BOOST_TEST(result == geometry::TOUCHING);
372
373 testPoint2D << 0.0, 0.500000000000001;
375 sidelengths2D, center2D, testPoint2D);
376 BOOST_TEST(result == geometry::TOUCHING);
377
378 testPoint2D << 0.0, -0.500000000000001;
380 sidelengths2D, center2D, testPoint2D);
381 BOOST_TEST(result == geometry::TOUCHING);
382
383 testPoint2D << 0.499999999999999, 0.0;
385 sidelengths2D, center2D, testPoint2D);
386 BOOST_TEST(result == geometry::TOUCHING);
387
388 testPoint2D << -0.499999999999999, 0.0;
390 sidelengths2D, center2D, testPoint2D);
391 BOOST_TEST(result == geometry::TOUCHING);
392
393 testPoint2D << 0.500000000000001, 0.0;
395 sidelengths2D, center2D, testPoint2D);
396 BOOST_TEST(result == geometry::TOUCHING);
397
398 testPoint2D << -0.500000000000001, 0.0;
400 sidelengths2D, center2D, testPoint2D);
401 BOOST_TEST(result == geometry::TOUCHING);
402
403 // 3D
404 Eigen::Vector3d center3D = Eigen::Vector3d::Zero();
405 Eigen::Vector3d sidelengths3D = Eigen::Vector3d::Constant(1.0);
406
407 // Not contained 3D
408 Eigen::Vector3d testPoint3D(2, 2, 2);
410 sidelengths3D, center3D, testPoint3D);
411 BOOST_TEST(result == geometry::NOT_CONTAINED);
412
413 testPoint3D << -2.0, -2.0, -2.0;
415 sidelengths3D, center3D, testPoint3D);
416 BOOST_TEST(result == geometry::NOT_CONTAINED);
417
418 testPoint3D << 2.0, -2.0, 2.0;
420 sidelengths3D, center3D, testPoint3D);
421 BOOST_TEST(result == geometry::NOT_CONTAINED);
422
423 testPoint3D << -2.0, 2.0, 2.0;
425 sidelengths3D, center3D, testPoint3D);
426 BOOST_TEST(result == geometry::NOT_CONTAINED);
427
428 testPoint3D << 2.0, 0.0, 0.0;
430 sidelengths3D, center3D, testPoint3D);
431 BOOST_TEST(result == geometry::NOT_CONTAINED);
432
433 testPoint3D << -2.0, 0.0, 0.0;
435 sidelengths3D, center3D, testPoint3D);
436 BOOST_TEST(result == geometry::NOT_CONTAINED);
437
438 testPoint3D << 0.0, 2.0, 0.0;
440 sidelengths3D, center3D, testPoint3D);
441 BOOST_TEST(result == geometry::NOT_CONTAINED);
442
443 testPoint3D << 0.0, 0.0, 2.0;
445 sidelengths3D, center3D, testPoint3D);
446 BOOST_TEST(result == geometry::NOT_CONTAINED);
447
448 testPoint3D << 0.0, 0.0, -2.0;
450 sidelengths3D, center3D, testPoint3D);
451 BOOST_TEST(result == geometry::NOT_CONTAINED);
452
453 // Contained 3D
454 testPoint3D << 0.0, 0.0, 0.0;
456 sidelengths3D, center3D, testPoint3D);
457 BOOST_TEST(result == geometry::CONTAINED);
458
459 testPoint3D << 0.25, 0.25, 0.25;
461 sidelengths3D, center3D, testPoint3D);
462 BOOST_TEST(result == geometry::CONTAINED);
463
464 testPoint3D << -0.25, 0.25, 0.25;
466 sidelengths3D, center3D, testPoint3D);
467 BOOST_TEST(result == geometry::CONTAINED);
468
469 testPoint3D << 0.25, -0.25, 0.25;
471 sidelengths3D, center3D, testPoint3D);
472 BOOST_TEST(result == geometry::CONTAINED);
473
474 testPoint3D << -0.25, -0.25, 0.25;
476 sidelengths3D, center3D, testPoint3D);
477 BOOST_TEST(result == geometry::CONTAINED);
478
479 testPoint3D << 0.25, 0.25, -0.25;
481 sidelengths3D, center3D, testPoint3D);
482 BOOST_TEST(result == geometry::CONTAINED);
483
484 testPoint3D << 0.49999999999, 0.49999999999, 0.49999999999;
486 sidelengths3D, center3D, testPoint3D);
487 BOOST_TEST(result == geometry::CONTAINED);
488
489 testPoint3D << -0.49999999999, -0.49999999999, -0.49999999999;
491 sidelengths3D, center3D, testPoint3D);
492 BOOST_TEST(result == geometry::CONTAINED);
493
494 // Touching 3D
495 testPoint3D << 0.5, 0.5, 0.5;
497 sidelengths3D, center3D, testPoint3D);
498 BOOST_TEST(result == geometry::TOUCHING);
499
500 testPoint3D << -0.5, 0.5, 0.5;
502 sidelengths3D, center3D, testPoint3D);
503 BOOST_TEST(result == geometry::TOUCHING);
504
505 testPoint3D << 0.5, -0.5, 0.5;
507 sidelengths3D, center3D, testPoint3D);
508 BOOST_TEST(result == geometry::TOUCHING);
509
510 testPoint3D << -0.5, -0.5, -0.5;
512 sidelengths3D, center3D, testPoint3D);
513 BOOST_TEST(result == geometry::TOUCHING);
514
515 testPoint3D << 0.4999999999999999, 0.4999999999999999, 0.4999999999999999;
517 sidelengths3D, center3D, testPoint3D);
518 BOOST_TEST(result == geometry::TOUCHING);
519
520 testPoint3D << 0.500000000000001, 0.500000000000001, 0.500000000000001;
522 sidelengths3D, center3D, testPoint3D);
523 BOOST_TEST(result == geometry::TOUCHING);
524
525 testPoint3D << -0.500000000000001, -0.500000000000001, -0.500000000000001;
527 sidelengths3D, center3D, testPoint3D);
528 BOOST_TEST(result == geometry::TOUCHING);
529
530 testPoint3D << 0.5, 0.0, 0.0;
532 sidelengths3D, center3D, testPoint3D);
533 BOOST_TEST(result == geometry::TOUCHING);
534
535 testPoint3D << -0.5, 0.0, 0.0;
537 sidelengths3D, center3D, testPoint3D);
538 BOOST_TEST(result == geometry::TOUCHING);
539
540 testPoint3D << 0.0, 0.5, 0.0;
542 sidelengths3D, center3D, testPoint3D);
543 BOOST_TEST(result == geometry::TOUCHING);
544
545 testPoint3D << 0.0, -0.5, 0.0;
547 sidelengths3D, center3D, testPoint3D);
548 BOOST_TEST(result == geometry::TOUCHING);
549
550 testPoint3D << 0.5, 0.0, -0.5;
552 sidelengths3D, center3D, testPoint3D);
553 BOOST_TEST(result == geometry::TOUCHING);
554
555 testPoint3D << -0.5, 0.0, 0.5;
557 sidelengths3D, center3D, testPoint3D);
558 BOOST_TEST(result == geometry::TOUCHING);
559
560 testPoint3D << 0.0, 0.499999999999999, 0.0;
562 sidelengths3D, center3D, testPoint3D);
563 BOOST_TEST(result == geometry::TOUCHING);
564
565 testPoint3D << 0.0, -0.499999999999999, 0.0;
567 sidelengths3D, center3D, testPoint3D);
568 BOOST_TEST(result == geometry::TOUCHING);
569
570 testPoint3D << 0.0, 0.500000000000001, 0.0;
572 sidelengths3D, center3D, testPoint3D);
573 BOOST_TEST(result == geometry::TOUCHING);
574
575 testPoint3D << 0.0, -0.500000000000001, 0.0;
577 sidelengths3D, center3D, testPoint3D);
578 BOOST_TEST(result == geometry::TOUCHING);
579
580 testPoint3D << 0.499999999999999, 0.0, 0.0;
582 sidelengths3D, center3D, testPoint3D);
583 BOOST_TEST(result == geometry::TOUCHING);
584
585 testPoint3D << -0.499999999999999, 0.0, 0.0;
587 sidelengths3D, center3D, testPoint3D);
588 BOOST_TEST(result == geometry::TOUCHING);
589
590 testPoint3D << 0.500000000000001, 0.0, 0.0;
592 sidelengths3D, center3D, testPoint3D);
593 BOOST_TEST(result == geometry::TOUCHING);
594
595 testPoint3D << -0.500000000000001, 0.0, 0.0;
597 sidelengths3D, center3D, testPoint3D);
598 BOOST_TEST(result == geometry::TOUCHING);
599
600 testPoint3D << 0.0, 0.0, 0.499999999999999;
602 sidelengths3D, center3D, testPoint3D);
603 BOOST_TEST(result == geometry::TOUCHING);
604
605 testPoint3D << 0.0, 0.0, -0.499999999999999;
607 sidelengths3D, center3D, testPoint3D);
608 BOOST_TEST(result == geometry::TOUCHING);
609
610 testPoint3D << 0.0, 0.0, 0.500000000000001;
612 sidelengths3D, center3D, testPoint3D);
613 BOOST_TEST(result == geometry::TOUCHING);
614
615 testPoint3D << 0.0, 0.0, -0.500000000000001;
617 sidelengths3D, center3D, testPoint3D);
618 BOOST_TEST(result == geometry::TOUCHING);
619}
620
621BOOST_AUTO_TEST_SUITE(Convexity)
622
623BOOST_AUTO_TEST_CASE(ComputeUnitQuadConvexity)
624{
625 PRECICE_TEST(1_rank);
626 int dim = 3;
627 Eigen::VectorXd coords0(dim);
628 Eigen::VectorXd coords1(dim);
629 Eigen::VectorXd coords2(dim);
630 Eigen::VectorXd coords3(dim);
631 coords0 << 0.0, 0.0, 0.0;
632 coords1 << 1.0, 0.0, 0.0;
633 coords2 << 1.0, 1.0, 0.0;
634 coords3 << 0.0, 1.0, 0.0;
635
636 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
637 auto result = geometry::isConvexQuad(vertexList);
638
639 BOOST_TEST(result.convex);
640 BOOST_TEST(utils::unique_elements(result.vertexOrder));
641 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
642 BOOST_TEST(result.vertexOrder.at(0) == 0);
643 BOOST_TEST(result.vertexOrder.at(1) == 3);
644 BOOST_TEST(result.vertexOrder.at(2) == 2);
645 BOOST_TEST(result.vertexOrder.at(3) == 1);
646}
647
648BOOST_AUTO_TEST_CASE(ComputeReversedUnitQuadConvexity)
649{
650 PRECICE_TEST(1_rank);
651 int dim = 3;
652 Eigen::VectorXd coords0(dim);
653 Eigen::VectorXd coords1(dim);
654 Eigen::VectorXd coords2(dim);
655 Eigen::VectorXd coords3(dim);
656 coords0 << 0.0, 1.0, 0.0;
657 coords1 << 1.0, 1.0, 0.0;
658 coords2 << 1.0, 0.0, 0.0;
659 coords3 << 0.0, 0.0, 0.0;
660
661 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
662 auto result = geometry::isConvexQuad(vertexList);
663
664 BOOST_TEST(result.convex);
665 BOOST_TEST(utils::unique_elements(result.vertexOrder));
666 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
667 BOOST_TEST(result.vertexOrder.at(0) == 0);
668 BOOST_TEST(result.vertexOrder.at(1) == 3);
669 BOOST_TEST(result.vertexOrder.at(2) == 2);
670 BOOST_TEST(result.vertexOrder.at(3) == 1);
671}
672
673BOOST_AUTO_TEST_CASE(ComputeValidQuadConvexity)
674{
675 PRECICE_TEST(1_rank);
676 int dim = 3;
677 Eigen::VectorXd coords0(dim);
678 Eigen::VectorXd coords1(dim);
679 Eigen::VectorXd coords2(dim);
680 Eigen::VectorXd coords3(dim);
681 coords0 << 0.5, 0.34, 0.0;
682 coords1 << 0.62, 0.32, 0.0;
683 coords2 << 0.6, 0.24, 0.0;
684 coords3 << 0.3, 0.22, 0.0;
685
686 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
687 auto result = geometry::isConvexQuad(vertexList);
688
689 BOOST_TEST(result.convex);
690 BOOST_TEST(utils::unique_elements(result.vertexOrder));
691 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
692 BOOST_TEST(result.vertexOrder.at(0) == 3);
693 BOOST_TEST(result.vertexOrder.at(1) == 2);
694 BOOST_TEST(result.vertexOrder.at(2) == 1);
695 BOOST_TEST(result.vertexOrder.at(3) == 0);
696}
697
698BOOST_AUTO_TEST_CASE(ComputeValidQuadConvexityWithOffPlane)
699{
700 PRECICE_TEST(1_rank);
701 int dim = 3;
702 Eigen::VectorXd coords0(dim);
703 Eigen::VectorXd coords1(dim);
704 Eigen::VectorXd coords2(dim);
705 Eigen::VectorXd coords3(dim);
706 coords0 << 0.5, 0.34, 0.0;
707 coords1 << 0.62, 0.32, 0.0;
708 coords2 << 0.6, 0.24, 0.0;
709 coords3 << 0.3, 0.22, 0.5;
710
711 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
712 // This test should result in an error
713 // auto result = geometry::isConvexQuad(vertexList);
714 //
715 // BOOST_TEST(result.convex);
716 // BOOST_TEST(utils::unique_elements(result.vertexOrder));
717 // BOOST_TEST_MESSAGE("Vertex Order" << result.vertexOrder);
718 // BOOST_TEST(result.vertexOrder.at(0) == 3);
719 // BOOST_TEST(result.vertexOrder.at(1) == 2);
720 // BOOST_TEST(result.vertexOrder.at(2) == 1);
721 // BOOST_TEST(result.vertexOrder.at(3) == 0);
722}
723
724BOOST_AUTO_TEST_CASE(ComputeInvalidUnitQuadConvexity)
725{
726 PRECICE_TEST(1_rank);
727 int dim = 3;
728 Eigen::VectorXd coords0(dim);
729 Eigen::VectorXd coords1(dim);
730 Eigen::VectorXd coords2(dim);
731 Eigen::VectorXd coords3(dim);
732 coords2 << 1.0, 0.0, 0.0;
733 coords1 << 1.0, 1.0, 0.0;
734 coords0 << 0.5, 0.9, 0.0;
735 coords3 << 0.0, 1.0, 0.0;
736
737 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
738 auto result = geometry::isConvexQuad(vertexList);
739
740 BOOST_TEST(!result.convex);
741}
742
743BOOST_AUTO_TEST_CASE(ComputeInvalidQuadConvexity)
744{
745 PRECICE_TEST(1_rank);
746 int dim = 3;
747 Eigen::VectorXd coords0(dim);
748 Eigen::VectorXd coords1(dim);
749 Eigen::VectorXd coords2(dim);
750 Eigen::VectorXd coords3(dim);
751 coords0 << 0.5, 0.34, 0.0;
752 coords1 << 0.62, 0.32, 0.0;
753 coords2 << 0.52, 0.31, 0.0;
754 coords3 << 0.51, 0.22, 0.0;
755
756 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
757 auto result = geometry::isConvexQuad(vertexList);
758
759 BOOST_TEST(!result.convex);
760}
761
762BOOST_AUTO_TEST_SUITE_END() // convexity
763
764BOOST_AUTO_TEST_SUITE_END() // geometry
765
BOOST_AUTO_TEST_CASE(Collinear)
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
#define PRECICE_TEST(...)
Definition Testing.hpp:27
bool between(const VECTORA_T &a, const VECTORB_T &b, const VECTORC_T &toCheck)
Determines, if a point lies on the line segment defined by a, b.
Definition geometry.hpp:151
double triangleArea(const Eigen::VectorXd &a, const Eigen::VectorXd &b, const Eigen::VectorXd &c)
Computes the signed area of a triangle in 2D.
Definition geometry.cpp:93
double tetraVolume(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const Eigen::Vector3d &d)
Computes the (unsigned) area of a triangle in 3D.
Definition geometry.cpp:116
Eigen::Vector2d projectVector(const Eigen::Vector3d &vector, int indexDimensionToRemove)
Projects a 3D vector to a 2D one by removing one dimension.
Definition geometry.cpp:125
int containedInHyperrectangle(const Eigen::MatrixBase< Derived > &sidelengths, const Eigen::MatrixBase< Derived > &center, const Eigen::MatrixBase< Derived > &testPoint)
Tests, if a vertex is contained in a hyperrectangle.
Definition geometry.hpp:198
ResultConstants segmentPlaneIntersection(const Eigen::Vector3d &pointOnPlane, const Eigen::Vector3d &planeNormal, const Eigen::Vector3d &firstPointSegment, const Eigen::Vector3d &secondPointSegment, Eigen::Vector3d &intersectionPoint)
Determines the intersection point of a segment with a plane in 3D.
Definition geometry.cpp:42
bool collinear(const Eigen::MatrixBase< Derived > &a, const Eigen::MatrixBase< Derived > &b, const Eigen::MatrixBase< Derived > &c)
Determines, if three points are collinear (on one line)
Definition geometry.hpp:170
ConvexityResult isConvexQuad(std::array< Eigen::VectorXd, 4 > coords)
Definition geometry.cpp:143
provides general mathematical constants and functions.
Definition barycenter.cpp:9
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.
auto make_array(Elements &&... elements) -> std::array< typename std::common_type< Elements... >::type, sizeof...(Elements)>
Function that generates an array from given elements.
Definition algorithm.hpp:51
bool unique_elements(const Container &c, BinaryPredicate p={})
Definition algorithm.hpp:64
Main namespace of the precice library.