preCICE v3.2.0
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
18{
20 // 2D test setup
21 Eigen::Vector2d a2D(0, 0);
22 Eigen::Vector2d b2D(1, 1);
23 Eigen::Vector2d collinearPoint2D(0.5, 0.5);
24 Eigen::Vector2d notCollinearPoint2D(0.5, 0.6);
25
26 // 3D test setup
27 Eigen::Vector3d a3D(0, 0, 0);
28 Eigen::Vector3d b3D(1, 1, 1);
29 Eigen::Vector3d collinearPoint3D(0.5, 0.5, 0.5);
30 Eigen::Vector3d notCollinearPoint3D(0.5, 0.6, 0.5);
31
32 // 2D test validations
33 BOOST_CHECK(geometry::collinear(a2D, b2D, collinearPoint2D));
34 BOOST_CHECK(!geometry::collinear(a2D, b2D, notCollinearPoint2D));
35
36 // 3D test validations
37 BOOST_CHECK(geometry::collinear(a3D, b3D, collinearPoint3D));
38 BOOST_CHECK(!geometry::collinear(a3D, b3D, notCollinearPoint3D));
39}
40
43 *boost::unit_test::tolerance(1e-3))
44{
46 Eigen::Vector3d a(1, 2, 3);
47 Eigen::Vector3d b(3, 2, 1);
48 Eigen::Vector3d c(4, 5, 6);
49 Eigen::Vector3d d(6, 5, 4);
50 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 0);
51
52 a << 5, 0, 0;
53 b << 0, -3, -3;
54 c << 0, 3, 4;
55 d << -1, -2, 6;
56 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 38.6666);
57
58 d << -1.47, -4.1, 8.3;
59 BOOST_TEST(geometry::tetraVolume(a, b, c, d) == 62.1816);
60}
61
64{
66 for (int dim = 2; dim <= 3; dim++) {
67 Eigen::VectorXd a(dim);
68 Eigen::VectorXd b(dim);
69 Eigen::VectorXd betweenPoint(dim);
70 Eigen::VectorXd betweenLimitPoint(dim);
71 Eigen::VectorXd collinearOutsidePoint(dim);
72 Eigen::VectorXd outsidePoint(dim);
73 if (dim == 2) {
74 a << 0.0, 0.0;
75 b << 1.0, 1.0;
76 betweenPoint << 0.5, 0.5;
77 betweenLimitPoint << 1.0, 1.0;
78 collinearOutsidePoint << 2.0, 2.0;
79 outsidePoint << 0.5, 0.4;
80 } else {
81 a << 0.0, 0.0, 0.0;
82 b << 1.0, 1.0, 1.0;
83 betweenPoint << 0.5, 0.5, 0.5;
84 betweenLimitPoint << 1.0, 1.0, 1.0;
85 collinearOutsidePoint << 2.0, 2.0, 2.0;
86 outsidePoint << 0.5, 0.4, 0.5;
87 }
88 BOOST_CHECK(geometry::between(a, b, betweenPoint));
89 BOOST_CHECK(geometry::between(a, b, betweenLimitPoint));
90 BOOST_CHECK(!geometry::between(a, b, collinearOutsidePoint));
91 BOOST_CHECK(!geometry::between(a, b, outsidePoint));
92 }
93}
94
97{
99 { // 2D
100 Eigen::Vector2d a, b, c;
101 double area;
102 a << 0.0, 0.0;
103 b << 1.0, 0.0;
104 c << 0.0, 1.0;
105 area = geometry::triangleArea(a, b, c);
106 BOOST_TEST(area == 0.5);
107
108 b << 0.0, 1.0;
109 c << 1.0, 0.0;
110 area = geometry::triangleArea(a, b, c);
111 // Areas are not signed.
112 BOOST_TEST(area == 0.5);
113 }
114 { // 3D
115 Eigen::Vector3d a, b, c;
116 double area;
117 a << 0.0, 0.0, 0.0;
118 b << 1.0, 0.0, 1.0;
119 c << 1.0, 1.0, 1.0;
120 area = geometry::triangleArea(a, b, c);
121 BOOST_CHECK(area == 0.5 * sqrt(2.0));
122 }
123}
124
125PRECICE_TEST_SETUP(1_rank)
126BOOST_AUTO_TEST_CASE(SegmentPlaneIntersection)
127{
128 PRECICE_TEST();
129 using Eigen::Vector3d;
130 Vector3d planeNormal = Vector3d::Constant(1.0);
131 Vector3d pointOnPlane = Vector3d::Constant(0.0);
132 Vector3d firstSegmentPoint = Vector3d::Constant(1.0);
133 Vector3d secondSegmentPoint = Vector3d::Constant(-1.0);
134 Vector3d intersectionPoint = Vector3d::Constant(1.0);
135 Vector3d expected = Vector3d::Constant(0.0);
136
137 // True intersection
139 pointOnPlane, planeNormal, firstSegmentPoint,
140 secondSegmentPoint, intersectionPoint);
141 BOOST_TEST(result == geometry::INTERSECTION);
142 BOOST_CHECK(equals(intersectionPoint, expected));
143
144 // Touching second segment vertex
145 secondSegmentPoint = Vector3d::Constant(0.0);
147 pointOnPlane, planeNormal, firstSegmentPoint,
148 secondSegmentPoint, intersectionPoint);
149 BOOST_TEST(result == geometry::TOUCHING);
150 BOOST_CHECK(equals(intersectionPoint, expected));
151
152 // Touching first segment vertex
153 firstSegmentPoint = Vector3d::Constant(0.0);
154 secondSegmentPoint = Vector3d::Constant(-1.0);
156 pointOnPlane, planeNormal, firstSegmentPoint,
157 secondSegmentPoint, intersectionPoint);
158 BOOST_TEST(result == geometry::TOUCHING);
159 BOOST_CHECK(equals(intersectionPoint, expected));
160
161 // Parallel segment with distance to plain
162 firstSegmentPoint << 0.0, 0.0, -3.0;
163 intersectionPoint << 1.0, 2.0, 3.0; // should not be modified
164 expected = intersectionPoint;
166 pointOnPlane, planeNormal, firstSegmentPoint,
167 secondSegmentPoint, intersectionPoint);
168 BOOST_TEST(result == geometry::NO_INTERSECTION);
169 BOOST_CHECK(equals(intersectionPoint, expected));
170
171 // Parallel segment contained in plane
172 firstSegmentPoint << 0.0, 0.0, 0.0;
173 secondSegmentPoint << 1.0, 1.0, -2.0;
175 pointOnPlane, planeNormal, firstSegmentPoint,
176 secondSegmentPoint, intersectionPoint);
177 BOOST_TEST(result == geometry::CONTAINED);
178 BOOST_CHECK(equals(intersectionPoint, expected));
179
180 // Segment ending before intersection
181 firstSegmentPoint << -2.0, -2.0, -2.0;
182 secondSegmentPoint << -1.0, -1.0, -1.0;
184 pointOnPlane, planeNormal, firstSegmentPoint,
185 secondSegmentPoint, intersectionPoint);
186 BOOST_TEST(result == geometry::NO_INTERSECTION);
187 BOOST_CHECK(equals(intersectionPoint, expected));
188
189 // Segment ending before intersection (inversed segment points)
190 firstSegmentPoint << -1.0, -1.0, -1.0;
191 secondSegmentPoint << -2.0, -2.0, -2.0;
193 pointOnPlane, planeNormal, firstSegmentPoint,
194 secondSegmentPoint, intersectionPoint);
195 BOOST_TEST(result == geometry::NO_INTERSECTION);
196 BOOST_CHECK(equals(intersectionPoint, expected));
197}
198
199PRECICE_TEST_SETUP(1_rank)
201{
202 PRECICE_TEST();
203 Eigen::Vector3d vector3D(1.0, 2.0, 3.0);
204 Eigen::Vector2d vector2D;
205 Eigen::Vector2d vectorExpected(1.0, 2.0);
206
207 vector2D = geometry::projectVector(vector3D, 2);
208 BOOST_CHECK(equals(vector2D, vectorExpected));
209
210 vector2D = geometry::projectVector(vector3D, 1);
211 vectorExpected << 1.0, 3.0;
212 BOOST_CHECK(equals(vector2D, vectorExpected));
213
214 vector2D = geometry::projectVector(vector3D, 0);
215 vectorExpected << 2.0, 3.0;
216 BOOST_CHECK(equals(vector2D, vectorExpected));
217}
218
219PRECICE_TEST_SETUP(1_rank)
220BOOST_AUTO_TEST_CASE(ContainedInHyperrectangle)
221{
222 PRECICE_TEST();
223 // 2D
224 Eigen::Vector2d center2D(0, 0);
225 Eigen::Vector2d sidelengths2D(1, 1);
226
227 // Not contained 2D
228 Eigen::Vector2d testPoint2D(2, 2);
230 sidelengths2D, center2D, testPoint2D);
231 BOOST_TEST(result == geometry::NOT_CONTAINED);
232
233 testPoint2D << -2.0, -2.0;
235 sidelengths2D, center2D, testPoint2D);
236 BOOST_TEST(result == geometry::NOT_CONTAINED);
237
238 testPoint2D << 2.0, -2.0;
240 sidelengths2D, center2D, testPoint2D);
241 BOOST_TEST(result == geometry::NOT_CONTAINED);
242
243 testPoint2D << -2.0, 2.0;
245 sidelengths2D, center2D, testPoint2D);
246 BOOST_TEST(result == geometry::NOT_CONTAINED);
247
248 testPoint2D << 2.0, 0.0;
250 sidelengths2D, center2D, testPoint2D);
251 BOOST_TEST(result == geometry::NOT_CONTAINED);
252
253 testPoint2D << -2.0, 0.0;
255 sidelengths2D, center2D, testPoint2D);
256 BOOST_TEST(result == geometry::NOT_CONTAINED);
257
258 testPoint2D << 0.0, 2.0;
260 sidelengths2D, center2D, testPoint2D);
261 BOOST_TEST(result == geometry::NOT_CONTAINED);
262
263 testPoint2D << 0.0, -2.0;
265 sidelengths2D, center2D, testPoint2D);
266 BOOST_TEST(result == geometry::NOT_CONTAINED);
267
268 // Contained 2D
269 testPoint2D << 0.0, 0.0;
271 sidelengths2D, center2D, testPoint2D);
272 BOOST_TEST(result == geometry::CONTAINED);
273
274 testPoint2D << 0.25, 0.25;
276 sidelengths2D, center2D, testPoint2D);
277 BOOST_TEST(result == geometry::CONTAINED);
278
279 testPoint2D << -0.25, 0.25;
281 sidelengths2D, center2D, testPoint2D);
282 BOOST_TEST(result == geometry::CONTAINED);
283
284 testPoint2D << 0.25, -0.25;
286 sidelengths2D, center2D, testPoint2D);
287 BOOST_TEST(result == geometry::CONTAINED);
288
289 testPoint2D << -0.25, -0.25;
291 sidelengths2D, center2D, testPoint2D);
292 BOOST_TEST(result == geometry::CONTAINED);
293
294 testPoint2D << 0.49999999999, 0.49999999999;
296 sidelengths2D, center2D, testPoint2D);
297 BOOST_TEST(result == geometry::CONTAINED);
298
299 testPoint2D << -0.49999999999, -0.49999999999;
301 sidelengths2D, center2D, testPoint2D);
302 BOOST_TEST(result == geometry::CONTAINED);
303
304 // Touching 2D
305 testPoint2D << 0.5, 0.5;
307 sidelengths2D, center2D, testPoint2D);
308 BOOST_TEST(result == geometry::TOUCHING);
309
310 testPoint2D << -0.5, 0.5;
312 sidelengths2D, center2D, testPoint2D);
313 BOOST_TEST(result == geometry::TOUCHING);
314
315 testPoint2D << 0.5, -0.5;
317 sidelengths2D, center2D, testPoint2D);
318 BOOST_TEST(result == geometry::TOUCHING);
319
320 testPoint2D << -0.5, -0.5;
322 sidelengths2D, center2D, testPoint2D);
323 BOOST_TEST(result == geometry::TOUCHING);
324
325 testPoint2D << 0.4999999999999999, 0.4999999999999999;
327 sidelengths2D, center2D, testPoint2D);
328 BOOST_TEST(result == geometry::TOUCHING);
329
330 testPoint2D << 0.500000000000001, 0.500000000000001;
332 sidelengths2D, center2D, testPoint2D);
333 BOOST_TEST(result == geometry::TOUCHING);
334
335 testPoint2D << -0.500000000000001, -0.500000000000001;
337 sidelengths2D, center2D, testPoint2D);
338 BOOST_TEST(result == geometry::TOUCHING);
339
340 testPoint2D << 0.5, 0.0;
342 sidelengths2D, center2D, testPoint2D);
343 BOOST_TEST(result == geometry::TOUCHING);
344
345 testPoint2D << -0.5, 0.0;
347 sidelengths2D, center2D, testPoint2D);
348 BOOST_TEST(result == geometry::TOUCHING);
349
350 testPoint2D << 0.0, 0.5;
352 sidelengths2D, center2D, testPoint2D);
353 BOOST_TEST(result == geometry::TOUCHING);
354
355 testPoint2D << 0.0, -0.5;
357 sidelengths2D, center2D, testPoint2D);
358 BOOST_TEST(result == geometry::TOUCHING);
359
360 testPoint2D << 0.5, 0.0;
362 sidelengths2D, center2D, testPoint2D);
363 BOOST_TEST(result == geometry::TOUCHING);
364
365 testPoint2D << -0.5, 0.0;
367 sidelengths2D, center2D, testPoint2D);
368 BOOST_TEST(result == geometry::TOUCHING);
369
370 testPoint2D << 0.0, 0.499999999999999;
372 sidelengths2D, center2D, testPoint2D);
373 BOOST_TEST(result == geometry::TOUCHING);
374
375 testPoint2D << 0.0, -0.499999999999999;
377 sidelengths2D, center2D, testPoint2D);
378 BOOST_TEST(result == geometry::TOUCHING);
379
380 testPoint2D << 0.0, 0.500000000000001;
382 sidelengths2D, center2D, testPoint2D);
383 BOOST_TEST(result == geometry::TOUCHING);
384
385 testPoint2D << 0.0, -0.500000000000001;
387 sidelengths2D, center2D, testPoint2D);
388 BOOST_TEST(result == geometry::TOUCHING);
389
390 testPoint2D << 0.499999999999999, 0.0;
392 sidelengths2D, center2D, testPoint2D);
393 BOOST_TEST(result == geometry::TOUCHING);
394
395 testPoint2D << -0.499999999999999, 0.0;
397 sidelengths2D, center2D, testPoint2D);
398 BOOST_TEST(result == geometry::TOUCHING);
399
400 testPoint2D << 0.500000000000001, 0.0;
402 sidelengths2D, center2D, testPoint2D);
403 BOOST_TEST(result == geometry::TOUCHING);
404
405 testPoint2D << -0.500000000000001, 0.0;
407 sidelengths2D, center2D, testPoint2D);
408 BOOST_TEST(result == geometry::TOUCHING);
409
410 // 3D
411 Eigen::Vector3d center3D = Eigen::Vector3d::Zero();
412 Eigen::Vector3d sidelengths3D = Eigen::Vector3d::Constant(1.0);
413
414 // Not contained 3D
415 Eigen::Vector3d testPoint3D(2, 2, 2);
417 sidelengths3D, center3D, testPoint3D);
418 BOOST_TEST(result == geometry::NOT_CONTAINED);
419
420 testPoint3D << -2.0, -2.0, -2.0;
422 sidelengths3D, center3D, testPoint3D);
423 BOOST_TEST(result == geometry::NOT_CONTAINED);
424
425 testPoint3D << 2.0, -2.0, 2.0;
427 sidelengths3D, center3D, testPoint3D);
428 BOOST_TEST(result == geometry::NOT_CONTAINED);
429
430 testPoint3D << -2.0, 2.0, 2.0;
432 sidelengths3D, center3D, testPoint3D);
433 BOOST_TEST(result == geometry::NOT_CONTAINED);
434
435 testPoint3D << 2.0, 0.0, 0.0;
437 sidelengths3D, center3D, testPoint3D);
438 BOOST_TEST(result == geometry::NOT_CONTAINED);
439
440 testPoint3D << -2.0, 0.0, 0.0;
442 sidelengths3D, center3D, testPoint3D);
443 BOOST_TEST(result == geometry::NOT_CONTAINED);
444
445 testPoint3D << 0.0, 2.0, 0.0;
447 sidelengths3D, center3D, testPoint3D);
448 BOOST_TEST(result == geometry::NOT_CONTAINED);
449
450 testPoint3D << 0.0, 0.0, 2.0;
452 sidelengths3D, center3D, testPoint3D);
453 BOOST_TEST(result == geometry::NOT_CONTAINED);
454
455 testPoint3D << 0.0, 0.0, -2.0;
457 sidelengths3D, center3D, testPoint3D);
458 BOOST_TEST(result == geometry::NOT_CONTAINED);
459
460 // Contained 3D
461 testPoint3D << 0.0, 0.0, 0.0;
463 sidelengths3D, center3D, testPoint3D);
464 BOOST_TEST(result == geometry::CONTAINED);
465
466 testPoint3D << 0.25, 0.25, 0.25;
468 sidelengths3D, center3D, testPoint3D);
469 BOOST_TEST(result == geometry::CONTAINED);
470
471 testPoint3D << -0.25, 0.25, 0.25;
473 sidelengths3D, center3D, testPoint3D);
474 BOOST_TEST(result == geometry::CONTAINED);
475
476 testPoint3D << 0.25, -0.25, 0.25;
478 sidelengths3D, center3D, testPoint3D);
479 BOOST_TEST(result == geometry::CONTAINED);
480
481 testPoint3D << -0.25, -0.25, 0.25;
483 sidelengths3D, center3D, testPoint3D);
484 BOOST_TEST(result == geometry::CONTAINED);
485
486 testPoint3D << 0.25, 0.25, -0.25;
488 sidelengths3D, center3D, testPoint3D);
489 BOOST_TEST(result == geometry::CONTAINED);
490
491 testPoint3D << 0.49999999999, 0.49999999999, 0.49999999999;
493 sidelengths3D, center3D, testPoint3D);
494 BOOST_TEST(result == geometry::CONTAINED);
495
496 testPoint3D << -0.49999999999, -0.49999999999, -0.49999999999;
498 sidelengths3D, center3D, testPoint3D);
499 BOOST_TEST(result == geometry::CONTAINED);
500
501 // Touching 3D
502 testPoint3D << 0.5, 0.5, 0.5;
504 sidelengths3D, center3D, testPoint3D);
505 BOOST_TEST(result == geometry::TOUCHING);
506
507 testPoint3D << -0.5, 0.5, 0.5;
509 sidelengths3D, center3D, testPoint3D);
510 BOOST_TEST(result == geometry::TOUCHING);
511
512 testPoint3D << 0.5, -0.5, 0.5;
514 sidelengths3D, center3D, testPoint3D);
515 BOOST_TEST(result == geometry::TOUCHING);
516
517 testPoint3D << -0.5, -0.5, -0.5;
519 sidelengths3D, center3D, testPoint3D);
520 BOOST_TEST(result == geometry::TOUCHING);
521
522 testPoint3D << 0.4999999999999999, 0.4999999999999999, 0.4999999999999999;
524 sidelengths3D, center3D, testPoint3D);
525 BOOST_TEST(result == geometry::TOUCHING);
526
527 testPoint3D << 0.500000000000001, 0.500000000000001, 0.500000000000001;
529 sidelengths3D, center3D, testPoint3D);
530 BOOST_TEST(result == geometry::TOUCHING);
531
532 testPoint3D << -0.500000000000001, -0.500000000000001, -0.500000000000001;
534 sidelengths3D, center3D, testPoint3D);
535 BOOST_TEST(result == geometry::TOUCHING);
536
537 testPoint3D << 0.5, 0.0, 0.0;
539 sidelengths3D, center3D, testPoint3D);
540 BOOST_TEST(result == geometry::TOUCHING);
541
542 testPoint3D << -0.5, 0.0, 0.0;
544 sidelengths3D, center3D, testPoint3D);
545 BOOST_TEST(result == geometry::TOUCHING);
546
547 testPoint3D << 0.0, 0.5, 0.0;
549 sidelengths3D, center3D, testPoint3D);
550 BOOST_TEST(result == geometry::TOUCHING);
551
552 testPoint3D << 0.0, -0.5, 0.0;
554 sidelengths3D, center3D, testPoint3D);
555 BOOST_TEST(result == geometry::TOUCHING);
556
557 testPoint3D << 0.5, 0.0, -0.5;
559 sidelengths3D, center3D, testPoint3D);
560 BOOST_TEST(result == geometry::TOUCHING);
561
562 testPoint3D << -0.5, 0.0, 0.5;
564 sidelengths3D, center3D, testPoint3D);
565 BOOST_TEST(result == geometry::TOUCHING);
566
567 testPoint3D << 0.0, 0.499999999999999, 0.0;
569 sidelengths3D, center3D, testPoint3D);
570 BOOST_TEST(result == geometry::TOUCHING);
571
572 testPoint3D << 0.0, -0.499999999999999, 0.0;
574 sidelengths3D, center3D, testPoint3D);
575 BOOST_TEST(result == geometry::TOUCHING);
576
577 testPoint3D << 0.0, 0.500000000000001, 0.0;
579 sidelengths3D, center3D, testPoint3D);
580 BOOST_TEST(result == geometry::TOUCHING);
581
582 testPoint3D << 0.0, -0.500000000000001, 0.0;
584 sidelengths3D, center3D, testPoint3D);
585 BOOST_TEST(result == geometry::TOUCHING);
586
587 testPoint3D << 0.499999999999999, 0.0, 0.0;
589 sidelengths3D, center3D, testPoint3D);
590 BOOST_TEST(result == geometry::TOUCHING);
591
592 testPoint3D << -0.499999999999999, 0.0, 0.0;
594 sidelengths3D, center3D, testPoint3D);
595 BOOST_TEST(result == geometry::TOUCHING);
596
597 testPoint3D << 0.500000000000001, 0.0, 0.0;
599 sidelengths3D, center3D, testPoint3D);
600 BOOST_TEST(result == geometry::TOUCHING);
601
602 testPoint3D << -0.500000000000001, 0.0, 0.0;
604 sidelengths3D, center3D, testPoint3D);
605 BOOST_TEST(result == geometry::TOUCHING);
606
607 testPoint3D << 0.0, 0.0, 0.499999999999999;
609 sidelengths3D, center3D, testPoint3D);
610 BOOST_TEST(result == geometry::TOUCHING);
611
612 testPoint3D << 0.0, 0.0, -0.499999999999999;
614 sidelengths3D, center3D, testPoint3D);
615 BOOST_TEST(result == geometry::TOUCHING);
616
617 testPoint3D << 0.0, 0.0, 0.500000000000001;
619 sidelengths3D, center3D, testPoint3D);
620 BOOST_TEST(result == geometry::TOUCHING);
621
622 testPoint3D << 0.0, 0.0, -0.500000000000001;
624 sidelengths3D, center3D, testPoint3D);
625 BOOST_TEST(result == geometry::TOUCHING);
626}
627
628BOOST_AUTO_TEST_SUITE(Convexity)
629
630PRECICE_TEST_SETUP(1_rank)
631BOOST_AUTO_TEST_CASE(ComputeUnitQuadConvexity)
632{
633 PRECICE_TEST();
634 int dim = 3;
635 Eigen::VectorXd coords0(dim);
636 Eigen::VectorXd coords1(dim);
637 Eigen::VectorXd coords2(dim);
638 Eigen::VectorXd coords3(dim);
639 coords0 << 0.0, 0.0, 0.0;
640 coords1 << 1.0, 0.0, 0.0;
641 coords2 << 1.0, 1.0, 0.0;
642 coords3 << 0.0, 1.0, 0.0;
643
644 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
645 auto result = geometry::isConvexQuad(vertexList);
646
647 BOOST_TEST(result.convex);
648 BOOST_TEST(utils::unique_elements(result.vertexOrder));
649 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
650 BOOST_TEST(result.vertexOrder.at(0) == 0);
651 BOOST_TEST(result.vertexOrder.at(1) == 3);
652 BOOST_TEST(result.vertexOrder.at(2) == 2);
653 BOOST_TEST(result.vertexOrder.at(3) == 1);
654}
655
656PRECICE_TEST_SETUP(1_rank)
657BOOST_AUTO_TEST_CASE(ComputeReversedUnitQuadConvexity)
658{
659 PRECICE_TEST();
660 int dim = 3;
661 Eigen::VectorXd coords0(dim);
662 Eigen::VectorXd coords1(dim);
663 Eigen::VectorXd coords2(dim);
664 Eigen::VectorXd coords3(dim);
665 coords0 << 0.0, 1.0, 0.0;
666 coords1 << 1.0, 1.0, 0.0;
667 coords2 << 1.0, 0.0, 0.0;
668 coords3 << 0.0, 0.0, 0.0;
669
670 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
671 auto result = geometry::isConvexQuad(vertexList);
672
673 BOOST_TEST(result.convex);
674 BOOST_TEST(utils::unique_elements(result.vertexOrder));
675 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
676 BOOST_TEST(result.vertexOrder.at(0) == 0);
677 BOOST_TEST(result.vertexOrder.at(1) == 3);
678 BOOST_TEST(result.vertexOrder.at(2) == 2);
679 BOOST_TEST(result.vertexOrder.at(3) == 1);
680}
681
682PRECICE_TEST_SETUP(1_rank)
683BOOST_AUTO_TEST_CASE(ComputeValidQuadConvexity)
684{
685 PRECICE_TEST();
686 int dim = 3;
687 Eigen::VectorXd coords0(dim);
688 Eigen::VectorXd coords1(dim);
689 Eigen::VectorXd coords2(dim);
690 Eigen::VectorXd coords3(dim);
691 coords0 << 0.5, 0.34, 0.0;
692 coords1 << 0.62, 0.32, 0.0;
693 coords2 << 0.6, 0.24, 0.0;
694 coords3 << 0.3, 0.22, 0.0;
695
696 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
697 auto result = geometry::isConvexQuad(vertexList);
698
699 BOOST_TEST(result.convex);
700 BOOST_TEST(utils::unique_elements(result.vertexOrder));
701 BOOST_TEST_MESSAGE(fmt::format("Vertex Order {}", result.vertexOrder));
702 BOOST_TEST(result.vertexOrder.at(0) == 3);
703 BOOST_TEST(result.vertexOrder.at(1) == 2);
704 BOOST_TEST(result.vertexOrder.at(2) == 1);
705 BOOST_TEST(result.vertexOrder.at(3) == 0);
706}
707
708PRECICE_TEST_SETUP(1_rank)
709BOOST_AUTO_TEST_CASE(ComputeValidQuadConvexityWithOffPlane)
710{
711 PRECICE_TEST();
712 int dim = 3;
713 Eigen::VectorXd coords0(dim);
714 Eigen::VectorXd coords1(dim);
715 Eigen::VectorXd coords2(dim);
716 Eigen::VectorXd coords3(dim);
717 coords0 << 0.5, 0.34, 0.0;
718 coords1 << 0.62, 0.32, 0.0;
719 coords2 << 0.6, 0.24, 0.0;
720 coords3 << 0.3, 0.22, 0.5;
721
722 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
723 // This test should result in an error
724 // auto result = geometry::isConvexQuad(vertexList);
725 //
726 // BOOST_TEST(result.convex);
727 // BOOST_TEST(utils::unique_elements(result.vertexOrder));
728 // BOOST_TEST_MESSAGE("Vertex Order" << result.vertexOrder);
729 // BOOST_TEST(result.vertexOrder.at(0) == 3);
730 // BOOST_TEST(result.vertexOrder.at(1) == 2);
731 // BOOST_TEST(result.vertexOrder.at(2) == 1);
732 // BOOST_TEST(result.vertexOrder.at(3) == 0);
733}
734
735PRECICE_TEST_SETUP(1_rank)
736BOOST_AUTO_TEST_CASE(ComputeInvalidUnitQuadConvexity)
737{
738 PRECICE_TEST();
739 int dim = 3;
740 Eigen::VectorXd coords0(dim);
741 Eigen::VectorXd coords1(dim);
742 Eigen::VectorXd coords2(dim);
743 Eigen::VectorXd coords3(dim);
744 coords2 << 1.0, 0.0, 0.0;
745 coords1 << 1.0, 1.0, 0.0;
746 coords0 << 0.5, 0.9, 0.0;
747 coords3 << 0.0, 1.0, 0.0;
748
749 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
750 auto result = geometry::isConvexQuad(vertexList);
751
752 BOOST_TEST(!result.convex);
753}
754
755PRECICE_TEST_SETUP(1_rank)
756BOOST_AUTO_TEST_CASE(ComputeInvalidQuadConvexity)
757{
758 PRECICE_TEST();
759 int dim = 3;
760 Eigen::VectorXd coords0(dim);
761 Eigen::VectorXd coords1(dim);
762 Eigen::VectorXd coords2(dim);
763 Eigen::VectorXd coords3(dim);
764 coords0 << 0.5, 0.34, 0.0;
765 coords1 << 0.62, 0.32, 0.0;
766 coords2 << 0.52, 0.31, 0.0;
767 coords3 << 0.51, 0.22, 0.0;
768
769 auto vertexList = utils::make_array(coords0, coords1, coords2, coords3);
770 auto result = geometry::isConvexQuad(vertexList);
771
772 BOOST_TEST(!result.convex);
773}
774
775BOOST_AUTO_TEST_SUITE_END() // convexity
776
777BOOST_AUTO_TEST_SUITE_END() // geometry
778
BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithSubsteps)
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
#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
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:148
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:195
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:167
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.
bool unique_elements(const Container &c, BinaryPredicate p={})
Definition algorithm.hpp:63
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:50
Main namespace of the precice library.