DGtal 2.1.0
Loading...
Searching...
No Matches
testConvexityHelper.cpp File Reference
#include <iostream>
#include <vector>
#include <algorithm>
#include "DGtal/base/Common.h"
#include "DGtal/kernel/SpaceND.h"
#include "DGtal/geometry/tools/AffineGeometry.h"
#include "DGtal/geometry/volumes/ConvexityHelper.h"
#include "DGtal/shapes/SurfaceMesh.h"
#include "DGtalCatch.h"
Include dependency graph for testConvexityHelper.cpp:

Go to the source code of this file.

Functions

 SCENARIO ("ConvexityHelper< 2 > unit tests", "[convexity_helper][lattice_polytope][2d]")
 
 SCENARIO ("ConvexityHelper< 3 > unit tests", "[convexity_helper][3d]")
 
 SCENARIO ("ConvexityHelper< 3 > triangle tests", "[convexity_helper][triangle][3d]")
 
 SCENARIO ("ConvexityHelper< 3 > degenerated triangle tests", "[convexity_helper][triangle][degenerate][3d]")
 
 SCENARIO ("ConvexityHelper< 3 > open segment tests", "[convexity_helper][open segment][3d]")
 
 SCENARIO ("ConvexityHelper< 3 > open triangle tests", "[convexity_helper][open triangle][3d]")
 
 SCENARIO ("ConvexityHelper< 3 > open triangle unit tests", "[convexity_helper][open triangle][3d]")
 

Detailed Description

This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.

Author
Jacques-Olivier Lachaud (jacqu.nosp@m.es-o.nosp@m.livie.nosp@m.r.la.nosp@m.chaud.nosp@m.@uni.nosp@m.v-sav.nosp@m.oie..nosp@m.fr ) Laboratory of Mathematics (CNRS, UMR 5127), University of Savoie, France
Date
2020/02/01

Functions for testing class ConvexityHelper.

This file is part of the DGtal library.

Definition in file testConvexityHelper.cpp.

Function Documentation

◆ SCENARIO() [1/7]

SCENARIO ( "ConvexityHelper< 2 > unit tests"  ,
""  [convexity_helper][lattice_polytope][2d] 
)

Definition at line 50 of file testConvexityHelper.cpp.

52{
53 typedef ConvexityHelper< 2 > Helper;
54 typedef Helper::Point Point;
55 typedef ConvexCellComplex< Point > CvxCellComplex;
56 GIVEN( "Given a star { (0,0), (-4,-1), (-3,5), (7,3), (5, -2) } " ) {
57 std::vector<Point> V
58 = { Point(0,0), Point(-4,-1), Point(-3,5), Point(7,3), Point(5, -2) };
59 WHEN( "Computing its lattice polytope" ){
60 const auto P = Helper::computeLatticePolytope( V, false, true );
61 CAPTURE( P );
62 THEN( "The polytope is valid and has 4 non trivial facets" ) {
63 REQUIRE( P.nbHalfSpaces() - 4 == 4 );
64 }
65 THEN( "The polytope is Minkowski summable" ) {
66 REQUIRE( P.canBeSummed() );
67 }
68 THEN( "The polytope contains the input points" ) {
69 REQUIRE( P.isInside( V[ 0 ] ) );
70 REQUIRE( P.isInside( V[ 1 ] ) );
71 REQUIRE( P.isInside( V[ 2 ] ) );
72 REQUIRE( P.isInside( V[ 3 ] ) );
73 REQUIRE( P.isInside( V[ 4 ] ) );
74 }
75 THEN( "The polytope contains 58 points" ) {
76 REQUIRE( P.count() == 58 );
77 }
78 THEN( "The interior of the polytope contains 53 points" ) {
79 REQUIRE( P.countInterior() == 53 );
80 }
81 THEN( "The boundary of the polytope contains 5 points" ) {
82 REQUIRE( P.countBoundary() == 5 );
83 }
84 }
85 }
86 GIVEN( "Given a square with an additional outside vertex " ) {
87 std::vector<Point> V
88 = { Point(-10,-10), Point(10,-10), Point(-10,10), Point(10,10),
89 Point(0,18) };
90 WHEN( "Computing its Delaunay cell complex" ){
91 CvxCellComplex complex;
92 bool ok = Helper::computeDelaunayCellComplex( complex, V, false );
93 CAPTURE( complex );
94 THEN( "The complex has 2 cells, 6 faces, 5 vertices" ) {
95 REQUIRE( ok );
96 REQUIRE( complex.nbCells() == 2 );
97 REQUIRE( complex.nbFaces() == 6 );
98 REQUIRE( complex.nbVertices() == 5 );
99 }
100 THEN( "The faces of cells are finite" ) {
101 bool ok_finite = true;
102 for ( CvxCellComplex::Size c = 0; c < complex.nbCells(); ++c ) {
103 const auto faces = complex.cellFaces( c );
104 for ( auto f : faces )
105 ok_finite = ok_finite && ! complex.isInfinite( complex.faceCell( f ) );
106 }
107 REQUIRE( ok_finite );
108 }
109 THEN( "The opposite of faces of cells are infinite except two" ) {
110 int nb_finite = 0;
111 for ( CvxCellComplex::Size c = 0; c < complex.nbCells(); ++c ) {
112 const auto faces = complex.cellFaces( c );
113 for ( auto f : faces ) {
114 const auto opp_f = complex.opposite( f );
115 nb_finite += complex.isInfinite( complex.faceCell( opp_f ) ) ? 0 : 1;
116 }
117 }
118 REQUIRE( nb_finite == 2 );
119 }
120 }
121 }
122 GIVEN( "Given a degenerated polytope { (0,0), (3,-1), (9,-3), (-6,2) } " ) {
123 std::vector<Point> V
124 = { Point(0,0), Point(3,-1), Point(9,-3), Point(-6,2) };
125 WHEN( "Computing its lattice polytope" ){
126 const auto P = Helper::computeLatticePolytope( V, false, true );
127 CAPTURE( P );
128 THEN( "The polytope is valid and has 2 non trivial facets" ) {
129 REQUIRE( P.nbHalfSpaces() - 4 == 2 );
130 }
131 THEN( "The polytope contains 6 points" ) {
132 REQUIRE( P.count() == 6 );
133 }
134 THEN( "The polytope contains no interior points" ) {
135 REQUIRE( P.countInterior() == 0 );
136 }
137 }
138 WHEN( "Computing the vertices of its convex hull" ){
139 auto X = Helper::computeConvexHullVertices( V, false );
140 std::sort( X.begin(), X.end() );
141 CAPTURE( X );
142 THEN( "The polytope is a segment defined by two points" ) {
143 REQUIRE( X.size() == 2 );
144 REQUIRE( X[ 0 ] == Point(-6, 2) );
145 REQUIRE( X[ 1 ] == Point( 9,-3) );
146 }
147 }
148 }
149 GIVEN( "Given a degenerated simplex { (4,0), (7,2), (-5,-6) } " ) {
150 std::vector<Point> V
151 = { Point(4,0), Point(7,2), Point(-5,-6) };
152 WHEN( "Computing its lattice polytope" ){
153 const auto P = Helper::computeLatticePolytope( V, false, true );
154 CAPTURE( P );
155 THEN( "The polytope is valid and has 2 non trivial facets" ) {
156 REQUIRE( P.nbHalfSpaces() - 4 == 2 );
157 }
158 THEN( "The polytope contains 5 points" ) {
159 REQUIRE( P.count() == 5 );
160 }
161 THEN( "The polytope contains no interior points" ) {
162 REQUIRE( P.countInterior() == 0 );
163 }
164 }
165 WHEN( "Computing the vertices of its convex hull" ){
166 auto X = Helper::computeConvexHullVertices( V, false );
167 std::sort( X.begin(), X.end() );
168 CAPTURE( X );
169 THEN( "The polytope is a segment defined by two points" ) {
170 REQUIRE( X.size() == 2 );
171 REQUIRE( X[ 0 ] == Point(-5,-6) );
172 REQUIRE( X[ 1 ] == Point( 7, 2) );
173 }
174 }
175 }
176}
Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and re...
Aim: Provides a set of functions to facilitate the computation of convex hulls and polytopes,...
MyPointD Point
CAPTURE(thicknessHV)
GIVEN("A cubical complex with random 3-cells")
REQUIRE(domain.isInside(aPoint))

References CAPTURE(), GIVEN(), and REQUIRE().

◆ SCENARIO() [2/7]

SCENARIO ( "ConvexityHelper< 3 > degenerated triangle tests"  ,
""  [convexity_helper][triangle][degenerate][3d] 
)

Definition at line 537 of file testConvexityHelper.cpp.

539{
540 typedef ConvexityHelper< 3 > Helper;
541 typedef Helper::Point Point;
542 typedef Helper::Vector Vector;
543 GIVEN( "Given degenerated triplets of points" ) {
544 Helper::LatticePolytope P;
545 Helper::LatticePolytope T;
546 Point a, b, c;
547 Vector n;
548 int nb_total = 0;
549 int nb_ok = 0;
550 int nbi_ok = 0;
551 int nb_P = 0, nb_T = 0, nbi_P = 0, nbi_T = 0;
552 for ( int i = 0; i < 20; i++ )
553 {
554 do {
555 a = Point( rand() % 10, rand() % 10, rand() % 10 );
556 b = Point( rand() % 10, rand() % 10, rand() % 10 );
557 c = Point( rand() % 10, rand() % 10, rand() % 10 );
558 n = (b-a).crossProduct(c-a);
559 } while ( n != Vector::zero );
560 std::vector<Point> V = { a, b, c };
561 P = Helper::computeLatticePolytope( V, true, false );
562 T = Helper::compute3DTriangle( a, b, c );
563 nb_P = P.count();
564 nb_T = T.count();
565 nbi_P = P.countInterior();
566 nbi_T = T.countInterior();
567 nb_total += 1;
568 nb_ok += ( nb_P == nb_T ) ? 1 : 0;
569 nbi_ok += ( nbi_P == 0 && nbi_T == 0 ) ? 1 : 0;
570 if ( ( nb_ok != nb_total ) || ( nbi_ok != nb_total ) ) break;
571 }
572 CAPTURE( a ); CAPTURE( b ); CAPTURE( c ); CAPTURE( n );
573 CAPTURE( P ); CAPTURE( T );
574 CAPTURE( nb_P ); CAPTURE( nb_T ); CAPTURE( nbi_P ); CAPTURE( nbi_T );
575 WHEN( "Computing their tightiest polytope or triangle" ) {
576 THEN( "They have the same number of inside points" ) {
577 REQUIRE( nb_ok == nb_total );
578 }
579 THEN( "They do not have interior points" ) {
580 REQUIRE( nbi_ok == nb_total );
581 }
582 }
583 }
584
585 GIVEN( "Given degenerated triplets of points" ) {
586 typedef Helper::LatticePolytope::UnitSegment UnitSegment;
587 Helper::LatticePolytope P;
588 Helper::LatticePolytope T;
589 Point a, b, c;
590 Vector n;
591 int nb_total = 0;
592 int nb_ok = 0;
593 int nbi_ok = 0;
594 int nb_P = 0, nb_T = 0, nbi_P = 0, nbi_T = 0;
595 for ( int i = 0; i < 20; i++ )
596 {
597 do {
598 a = Point( rand() % 10, rand() % 10, rand() % 10 );
599 b = Point( rand() % 10, rand() % 10, rand() % 10 );
600 c = Point( rand() % 10, rand() % 10, rand() % 10 );
601 n = (b-a).crossProduct(c-a);
602 } while ( n != Vector::zero );
603 std::vector<Point> V = { a, b, c };
604 P = Helper::computeLatticePolytope( V, true, true );
605 T = Helper::compute3DTriangle( a, b, c, true );
606 for ( Dimension k = 0; k < 3; k++ )
607 {
608 P += UnitSegment( k );
609 T += UnitSegment( k );
610 }
611 nb_P = P.count();
612 nb_T = T.count();
613 nbi_P = P.countInterior();
614 nbi_T = T.countInterior();
615 nb_total += 1;
616 nb_ok += ( nb_P == nb_T ) ? 1 : 0;
617 nbi_ok += ( nbi_P == nbi_T ) ? 1 : 0;
618 if ( ( nb_ok != nb_total ) || ( nbi_ok != nb_total ) ) break;
619 }
620 CAPTURE( a ); CAPTURE( b ); CAPTURE( c ); CAPTURE( n );
621 CAPTURE( P ); CAPTURE( T );
622 CAPTURE( nb_P ); CAPTURE( nb_T ); CAPTURE( nbi_P ); CAPTURE( nbi_T );
623 WHEN( "Computing their tightiest polytope or triangle, dilated by a cube" ) {
624 THEN( "They have the same number of inside points" ) {
625 REQUIRE( nb_ok == nb_total );
626 }
627 THEN( "They have the same number of interior points" ) {
628 REQUIRE( nbi_ok == nb_total );
629 }
630 }
631 }
632}
DigitalPlane::Point Vector
auto crossProduct(PointVector< 3, LeftEuclideanRing, LeftContainer > const &lhs, PointVector< 3, RightEuclideanRing, RightContainer > const &rhs) -> decltype(DGtal::constructFromArithmeticConversion(lhs, rhs))
Cross product of two 3D Points/Vectors.
DGtal::uint32_t Dimension
Definition Common.h:119

References CAPTURE(), DGtal::crossProduct(), GIVEN(), and REQUIRE().

◆ SCENARIO() [3/7]

SCENARIO ( "ConvexityHelper< 3 > open segment tests"  ,
""  [convexity_helper][open segment][3d] 
)

Definition at line 635 of file testConvexityHelper.cpp.

637{
638 typedef ConvexityHelper< 3 > Helper;
639 typedef Helper::Point Point;
640// typedef Helper::Vector Vector;
641 typedef Helper::LatticePolytope::UnitSegment UnitSegment;
642
643 auto nb_open_segment_smaller_than_segment = 0;
644 auto nb_dilated_open_segment_smaller_than_dilated_segment = 0;
645 auto nb_dilated_open_segment_greater_than_interior_dilated_segment = 0;
646 const int N = 100;
647 const int K = 10;
648 for ( auto n = 0; n < N; n++ )
649 {
650 Point a( rand() % K, rand() % K, rand() % K );
651 Point b( rand() % K, rand() % K, rand() % K );
652 if ( a == b ) b[ rand() % 3 ] += 1;
653 Helper::LatticePolytope CS = Helper::computeSegment( a, b );
654 Helper::LatticePolytope OS = Helper::computeOpenSegment( a, b );
655 {
656 CAPTURE( CS );
657 CAPTURE( OS );
658 auto cs_in = CS.count();
659 auto os_in = OS.count();
660 REQUIRE( os_in < cs_in );
661 nb_open_segment_smaller_than_segment += ( os_in < cs_in ) ? 1 : 0;
662 }
663 for ( Dimension k = 0; k < 3; k++ )
664 {
665 CS += UnitSegment( k );
666 OS += UnitSegment( k );
667 }
668 {
669 CAPTURE( CS );
670 CAPTURE( OS );
671 auto cs_in = CS.count();
672 auto cs_int = CS.countInterior();
673 auto os_in = OS.count();
674 REQUIRE( os_in < cs_in );
675 REQUIRE( cs_int <= os_in );
676 nb_dilated_open_segment_smaller_than_dilated_segment
677 += ( os_in < cs_in ) ? 1 : 0;
678 nb_dilated_open_segment_greater_than_interior_dilated_segment
679 += ( cs_int <= os_in ) ? 1 : 0;
680 }
681 }
682 WHEN( "Computing open segments" ) {
683 THEN( "They contain less lattice points than closed segments" ) {
684 REQUIRE( nb_open_segment_smaller_than_segment == N );
685 }
686 THEN( "When dilated, they contain less lattice points than dilated closed segments" ) {
687 REQUIRE( nb_dilated_open_segment_smaller_than_dilated_segment == N );
688 }
689 THEN( "When dilated, they contain more lattice points than the interior of dilated closed segments" ) {
690 REQUIRE( nb_dilated_open_segment_greater_than_interior_dilated_segment == N );
691 }
692 }
693}
KSpace K

References CAPTURE(), K, and REQUIRE().

◆ SCENARIO() [4/7]

SCENARIO ( "ConvexityHelper< 3 > open triangle tests"  ,
""  [convexity_helper][open triangle][3d] 
)

Definition at line 696 of file testConvexityHelper.cpp.

698{
699 typedef ConvexityHelper< 3 > Helper;
700 typedef Helper::Point Point;
701// typedef Helper::Vector Vector;
702 typedef Helper::LatticePolytope::UnitSegment UnitSegment;
703
704 auto nb_open_triangle_smaller_than_triangle = 0;
705 auto nb_dilated_open_triangle_smaller_than_dilated_triangle = 0;
706 auto nb_dilated_open_triangle_greater_than_interior_dilated_triangle = 0;
707 const int N = 100;
708 const int K = 10;
709 for ( auto n = 0; n < N; n++ )
710 {
711 Point a, b, c;
712 Dimension d = 0;
713 do {
714 a = Point( rand() % K, rand() % K, rand() % K );
715 b = Point( rand() % K, rand() % K, rand() % K );
716 c = Point( rand() % K, rand() % K, rand() % K );
717 std::vector< Point > X = { a, b, c };
719 } while ( d != 2 );
720 CAPTURE( a, b, c );
721 Helper::LatticePolytope CS = Helper::compute3DTriangle( a, b, c, true );
722 Helper::LatticePolytope OS = Helper::compute3DOpenTriangle( a, b, c, true );
723 {
724 CAPTURE( CS );
725 CAPTURE( OS );
726 auto cs_in = CS.count();
727 auto os_in = OS.count();
728 REQUIRE( os_in < cs_in );
729 nb_open_triangle_smaller_than_triangle += ( os_in < cs_in ) ? 1 : 0;
730 }
731 for ( Dimension k = 0; k < 3; k++ )
732 {
733 CS += UnitSegment( k );
734 OS += UnitSegment( k );
735 }
736 {
737 CAPTURE( CS );
738 CAPTURE( OS );
739 auto cs_in = CS.count();
740 auto cs_int = CS.countInterior();
741 auto os_in = OS.count();
742 REQUIRE( os_in < cs_in );
743 REQUIRE( cs_int <= os_in );
744 nb_dilated_open_triangle_smaller_than_dilated_triangle
745 += ( os_in < cs_in ) ? 1 : 0;
746 nb_dilated_open_triangle_greater_than_interior_dilated_triangle
747 += ( cs_int <= os_in ) ? 1 : 0;
748 }
749 }
750 WHEN( "Computing open triangles" ) {
751 THEN( "They contain less lattice points than closed triangles" ) {
752 REQUIRE( nb_open_triangle_smaller_than_triangle == N );
753 }
754 THEN( "When dilated, they contain less lattice points than dilated closed triangles" ) {
755 REQUIRE( nb_dilated_open_triangle_smaller_than_dilated_triangle == N );
756 }
757 THEN( "When dilated, they contain more lattice points than the interior of dilated closed triangles" ) {
758 REQUIRE( nb_dilated_open_triangle_greater_than_interior_dilated_triangle == N );
759 }
760 }
761}
static DGtal::int64_t affineDimension(const std::vector< TInputPoint > &X, const double tolerance=1e-12)

References DGtal::AffineGeometry< TPoint >::affineDimension(), CAPTURE(), K, and REQUIRE().

◆ SCENARIO() [5/7]

SCENARIO ( "ConvexityHelper< 3 > open triangle unit tests"  ,
""  [convexity_helper][open triangle][3d] 
)

Definition at line 763 of file testConvexityHelper.cpp.

765{
766 typedef ConvexityHelper< 3 > Helper;
767 typedef Helper::Point Point;
768// typedef Helper::Vector Vector;
769 typedef Helper::LatticePolytope::UnitSegment UnitSegment;
770
771 Point a( 0, 0, 0 );
772 Point b( 2, 1, 1 );
773 Point c( 2, 2, 1 );
774 CAPTURE( a, b, c );
775 Helper::LatticePolytope CS = Helper::compute3DTriangle( a, b, c, true );
776 Helper::LatticePolytope OS = Helper::compute3DOpenTriangle( a, b, c, true );
777 for ( Dimension k = 0; k < 3; k++ )
778 {
779 CS += UnitSegment( k );
780 OS += UnitSegment( k );
781 }
782 std::vector< Point > PCS, POS;
783 CS.getPoints( PCS );
784 OS.getPoints( POS );
785 CAPTURE( PCS, POS );
786 REQUIRE( PCS.size() == 21 );
787 REQUIRE( POS.size() == 3 );
788}

References CAPTURE(), and REQUIRE().

◆ SCENARIO() [6/7]

SCENARIO ( "ConvexityHelper< 3 > triangle tests"  ,
""  [convexity_helper][triangle][3d] 
)

Definition at line 439 of file testConvexityHelper.cpp.

441{
442 typedef ConvexityHelper< 3 > Helper;
443 typedef Helper::Point Point;
444 typedef Helper::Vector Vector;
445 GIVEN( "Given non degenerated triplets of points" ) {
446 Helper::LatticePolytope P;
447 Helper::LatticePolytope T;
448 Point a, b, c;
449 Vector n;
450 int nb_total = 0;
451 int nb_ok = 0;
452 int nbi_ok = 0;
453 int nb_P = 0, nb_T = 0, nbi_P = 0, nbi_T = 0;
454 for ( int i = 0; i < 20; i++ )
455 {
456 do {
457 a = Point( rand() % 10, rand() % 10, rand() % 10 );
458 b = Point( rand() % 10, rand() % 10, rand() % 10 );
459 c = Point( rand() % 10, rand() % 10, rand() % 10 );
460 n = (b-a).crossProduct(c-a);
461 } while ( n == Vector::zero );
462 std::vector<Point> V = { a, b, c };
463 P = Helper::computeLatticePolytope( V, false, false );
464 T = Helper::compute3DTriangle( a, b, c );
465 nb_P = P.count();
466 nb_T = T.count();
467 nbi_P = P.countInterior();
468 nbi_T = T.countInterior();
469 nb_total += 1;
470 nb_ok += ( nb_P == nb_T ) ? 1 : 0;
471 nbi_ok += ( nbi_P == 0 && nbi_T == 0 ) ? 1 : 0;
472 if ( ( nb_ok != nb_total ) || ( nbi_ok != nb_total ) ) break;
473 }
474 CAPTURE( a ); CAPTURE( b ); CAPTURE( c ); CAPTURE( n );
475 CAPTURE( P ); CAPTURE( T );
476 CAPTURE( nb_P ); CAPTURE( nb_T ); CAPTURE( nbi_P ); CAPTURE( nbi_T );
477 WHEN( "Computing their tightiest polytope or triangle" ) {
478 THEN( "They have the same number of inside points" ) {
479 REQUIRE( nb_ok == nb_total );
480 }
481 THEN( "They do not have interior points" ) {
482 REQUIRE( nbi_ok == nb_total );
483 }
484 }
485 }
486
487 GIVEN( "Given non degenerated triplets of points" ) {
488 typedef Helper::LatticePolytope::UnitSegment UnitSegment;
489 Helper::LatticePolytope P;
490 Helper::LatticePolytope T;
491 Point a, b, c;
492 Vector n;
493 int nb_total = 0;
494 int nb_ok = 0;
495 int nbi_ok = 0;
496 int nb_P = 0, nb_T = 0, nbi_P = 0, nbi_T = 0;
497 for ( int i = 0; i < 20; i++ )
498 {
499 do {
500 a = Point( rand() % 10, rand() % 10, rand() % 10 );
501 b = Point( rand() % 10, rand() % 10, rand() % 10 );
502 c = Point( rand() % 10, rand() % 10, rand() % 10 );
503 n = (b-a).crossProduct(c-a);
504 } while ( n == Vector::zero );
505 std::vector<Point> V = { a, b, c };
506 P = Helper::computeLatticePolytope( V, false, true );
507 T = Helper::compute3DTriangle( a, b, c, true );
508 for ( Dimension k = 0; k < 3; k++ )
509 {
510 P += UnitSegment( k );
511 T += UnitSegment( k );
512 }
513 nb_P = P.count();
514 nb_T = T.count();
515 nbi_P = P.countInterior();
516 nbi_T = T.countInterior();
517 nb_total += 1;
518 nb_ok += ( nb_P == nb_T ) ? 1 : 0;
519 nbi_ok += ( nbi_P == nbi_T ) ? 1 : 0;
520 if ( ( nb_ok != nb_total ) || ( nbi_ok != nb_total ) ) break;
521 }
522 CAPTURE( a ); CAPTURE( b ); CAPTURE( c ); CAPTURE( n );
523 CAPTURE( P ); CAPTURE( T );
524 CAPTURE( nb_P ); CAPTURE( nb_T ); CAPTURE( nbi_P ); CAPTURE( nbi_T );
525 WHEN( "Computing their tightiest polytope or triangle, dilated by a cube" ) {
526 THEN( "They have the same number of inside points" ) {
527 REQUIRE( nb_ok == nb_total );
528 }
529 THEN( "They have the same number of interior points" ) {
530 REQUIRE( nbi_ok == nb_total );
531 }
532 }
533 }
534}

References CAPTURE(), DGtal::crossProduct(), GIVEN(), and REQUIRE().

◆ SCENARIO() [7/7]

SCENARIO ( "ConvexityHelper< 3 > unit tests"  ,
""  [convexity_helper][3d] 
)

Definition at line 182 of file testConvexityHelper.cpp.

184{
185 typedef ConvexityHelper< 3 > Helper;
186 typedef Helper::Point Point;
187 typedef Helper::RealPoint RealPoint;
188 typedef Helper::RealVector RealVector;
190 typedef PolygonalSurface< Point > LatticePolySurf;
191 typedef ConvexCellComplex< Point > CvxCellComplex;
192 GIVEN( "Given an octahedron star { (0,0,0), (-2,0,0), (2,0,0), (0,-2,0), (0,2,0), (0,0,-2), (0,0,2) } " ) {
193 std::vector<Point> V
194 = { Point(0,0,0), Point(-2,0,0), Point(2,0,0), Point(0,-2,0), Point(0,2,0),
195 Point(0,0,-2), Point(0,0,2) };
196 WHEN( "Computing its lattice polytope" ){
197 const auto P = Helper::computeLatticePolytope( V, false, true );
198 CAPTURE( P );
199 THEN( "The polytope is valid and has 8 non trivial facets plus 12 edge constraints" ) {
200 REQUIRE( P.nbHalfSpaces() - 6 == 20 );
201 }
202 THEN( "The polytope is Minkowski summable" ) {
203 REQUIRE( P.canBeSummed() );
204 }
205 THEN( "The polytope contains the input points" ) {
206 REQUIRE( P.isInside( V[ 0 ] ) );
207 REQUIRE( P.isInside( V[ 1 ] ) );
208 REQUIRE( P.isInside( V[ 2 ] ) );
209 REQUIRE( P.isInside( V[ 3 ] ) );
210 REQUIRE( P.isInside( V[ 4 ] ) );
211 REQUIRE( P.isInside( V[ 5 ] ) );
212 REQUIRE( P.isInside( V[ 6 ] ) );
213 }
214 THEN( "The polytope contains 25 points" ) {
215 REQUIRE( P.count() == 25 );
216 }
217 THEN( "The interior of the polytope contains 7 points" ) {
218 REQUIRE( P.countInterior() == 7 );
219 }
220 THEN( "The boundary of the polytope contains 18 points" ) {
221 REQUIRE( P.countBoundary() == 18 );
222 }
223 }
224 WHEN( "Computing the boundary of its convex hull as a SurfaceMesh" ){
225 SMesh smesh;
226 bool ok = Helper::computeConvexHullBoundary( smesh, V, false );
227 CAPTURE( smesh );
228 THEN( "The surface mesh is valid and has 6 vertices, 12 edges and 8 faces" ) {
229 REQUIRE( ok );
230 REQUIRE( smesh.nbVertices() == 6 );
231 REQUIRE( smesh.nbEdges() == 12 );
232 REQUIRE( smesh.nbFaces() == 8 );
233 }
234 THEN( "The surface mesh has the topology of a sphere" ) {
235 REQUIRE( smesh.Euler() == 2 );
236 REQUIRE( smesh.computeManifoldBoundaryEdges().size() == 0 );
237 REQUIRE( smesh.computeNonManifoldEdges().size() == 0 );
238 }
239 }
240 WHEN( "Computing the boundary of its convex hull as a lattice PolygonalSurface" ){
241 LatticePolySurf lpsurf;
242 bool ok = Helper::computeConvexHullBoundary( lpsurf, V, false );
243 CAPTURE( lpsurf );
244 THEN( "The polygonal surface is valid and has 6 vertices, 12 edges and 8 faces" ) {
245 REQUIRE( ok );
246 REQUIRE( lpsurf.nbVertices() == 6 );
247 REQUIRE( lpsurf.nbEdges() == 12 );
248 REQUIRE( lpsurf.nbFaces() == 8 );
249 REQUIRE( lpsurf.nbArcs() == 24 );
250 }
251 THEN( "The polygonal surface has the topology of a sphere and no boundary" ) {
252 REQUIRE( lpsurf.Euler() == 2 );
253 REQUIRE( lpsurf.allBoundaryArcs().size() == 0 );
254 REQUIRE( lpsurf.allBoundaryVertices().size() == 0 );
255 }
256 }
257 WHEN( "Computing its convex hull as a ConvexCellComplex" ){
258 CvxCellComplex complex;
259 bool ok = Helper::computeConvexHullCellComplex( complex, V, false );
260 CAPTURE( complex );
261 THEN( "The convex cell complex is valid and has 6 vertices, 8 faces and 1 finite cell" ) {
262 REQUIRE( ok );
263 REQUIRE( complex.nbVertices() == 6 );
264 REQUIRE( complex.nbFaces() == 8 );
265 REQUIRE( complex.nbCells() == 1 );
266 }
267 }
268 WHEN( "Computing the vertices of its convex hull" ){
269 const auto X = Helper::computeConvexHullVertices( V, false );
270 CAPTURE( X );
271 THEN( "The polytope has 6 vertices" ) {
272 REQUIRE( X.size() == 6 );
273 }
274 }
275 }
276 GIVEN( "Given a cube with an additional outside vertex " ) {
277 std::vector<Point> V
278 = { Point(-10,-10,-10), Point(10,-10,-10), Point(-10,10,-10), Point(10,10,-10),
279 Point(-10,-10,10), Point(10,-10,10), Point(-10,10,10), Point(10,10,10),
280 Point(0,0,18) };
281 WHEN( "Computing its Delaunay cell complex" ){
282 CvxCellComplex complex;
283 bool ok = Helper::computeDelaunayCellComplex( complex, V, false );
284 CAPTURE( complex );
285 THEN( "The complex has 2 cells, 10 faces, 9 vertices" ) {
286 REQUIRE( ok );
287 REQUIRE( complex.nbCells() == 2 );
288 REQUIRE( complex.nbFaces() == 10 );
289 REQUIRE( complex.nbVertices() == 9 );
290 }
291 THEN( "The faces of cells are finite" ) {
292 bool ok_finite = true;
293 for ( CvxCellComplex::Size c = 0; c < complex.nbCells(); ++c ) {
294 const auto faces = complex.cellFaces( c );
295 for ( auto f : faces )
296 ok_finite = ok_finite && ! complex.isInfinite( complex.faceCell( f ) );
297 }
298 REQUIRE( ok_finite );
299 }
300 THEN( "The opposite of faces of cells are infinite except two" ) {
301 int nb_finite = 0;
302 for ( CvxCellComplex::Size c = 0; c < complex.nbCells(); ++c ) {
303 const auto faces = complex.cellFaces( c );
304 for ( auto f : faces ) {
305 const auto opp_f = complex.opposite( f );
306 nb_finite += complex.isInfinite( complex.faceCell( opp_f ) ) ? 0 : 1;
307 }
308 }
309 REQUIRE( nb_finite == 2 );
310 }
311 }
312 WHEN( "Computing the vertices of its convex hull" ){
313 const auto X = Helper::computeConvexHullVertices( V, false );
314 CAPTURE( X );
315 THEN( "The polytope has 9 vertices" ) {
316 REQUIRE( X.size() == 9 );
317 }
318 }
319 }
320 GIVEN( "Given a degenerated 1d polytope { (0,0,1), (3,-1,2), (9,-3,4), (-6,2,-1) } " ) {
321 std::vector<Point> V
322 = { Point(0,0,1), Point(3,-1,2), Point(9,-3,4), Point(-6,2,-1) };
323 WHEN( "Computing its lattice polytope" ){
324 const auto P = Helper::computeLatticePolytope( V, false, true );
325 CAPTURE( P );
326 THEN( "The polytope is valid and has 6 non trivial facets" ) {
327 REQUIRE( P.nbHalfSpaces() - 6 == 6 );
328 }
329 THEN( "The polytope contains 6 points" ) {
330 REQUIRE( P.count() == 6 );
331 }
332 THEN( "The polytope contains no interior points" ) {
333 REQUIRE( P.countInterior() == 0 );
334 }
335 }
336 WHEN( "Computing the vertices of its convex hull" ){
337 auto X = Helper::computeConvexHullVertices( V, false );
338 std::sort( X.begin(), X.end() );
339 CAPTURE( X );
340 THEN( "The polytope is a segment defined by two points" ) {
341 REQUIRE( X.size() == 2 );
342 REQUIRE( X[ 0 ] == Point(-6, 2,-1) );
343 REQUIRE( X[ 1 ] == Point( 9,-3, 4) );
344 }
345 }
346 }
347 GIVEN( "Given a degenerated 1d simplex { (1,0,-1), Point(4,-1,-2), Point(10,-3,-4) } " ) {
348 std::vector<Point> V
349 = { Point(1,0,-1), Point(4,-1,-2), Point(10,-3,-4) };
350 WHEN( "Computing its lattice polytope" ){
351 const auto P = Helper::computeLatticePolytope( V, false, true );
352 CAPTURE( P );
353 THEN( "The polytope is valid and has 6 non trivial facets" ) {
354 REQUIRE( P.nbHalfSpaces() - 6 == 6 );
355 }
356 THEN( "The polytope contains 4 points" ) {
357 REQUIRE( P.count() == 4 );
358 }
359 THEN( "The polytope contains no interior points" ) {
360 REQUIRE( P.countInterior() == 0 );
361 }
362 }
363 WHEN( "Computing the vertices of its convex hull" ){
364 auto X = Helper::computeConvexHullVertices( V, false );
365 std::sort( X.begin(), X.end() );
366 CAPTURE( X );
367 THEN( "The polytope is a segment defined by two points" ) {
368 REQUIRE( X.size() == 2 );
369 REQUIRE( X[ 0 ] == Point( 1, 0,-1) );
370 REQUIRE( X[ 1 ] == Point(10,-3,-4) );
371 }
372 }
373 }
374 GIVEN( "Given a degenerated 2d polytope { (2,1,0), (1,0,1), (1,2,1), (0,1,2), (0,3,2) } " ) {
375 std::vector<Point> V
376 = { Point(2,1,0), Point(1,0,1), Point(1,2,1), Point(0,1,2), Point(0,3,2) };
377 WHEN( "Computing its lattice polytope" ){
378 const auto P = Helper::computeLatticePolytope( V, false, true );
379 CAPTURE( P );
380 THEN( "The polytope is valid and has more than 6 non trivial facets" ) {
381 REQUIRE( P.nbHalfSpaces() - 6 == 6 );
382 }
383 THEN( "The polytope contains 7 points" ) {
384 REQUIRE( P.count() == 7 );
385 }
386 THEN( "The polytope contains no interior points" ) {
387 REQUIRE( P.countInterior() == 0 );
388 }
389 }
390 WHEN( "Computing the vertices of its convex hull" ){
391 auto X = Helper::computeConvexHullVertices( V, false );
392 std::sort( X.begin(), X.end() );
393 CAPTURE( X );
394 THEN( "The polytope is a quad" ) {
395 REQUIRE( X.size() == 4 );
396 REQUIRE( X[ 0 ] == Point( 0, 1, 2) );
397 REQUIRE( X[ 1 ] == Point( 0, 3, 2) );
398 REQUIRE( X[ 2 ] == Point( 1, 0, 1) );
399 REQUIRE( X[ 3 ] == Point( 2, 1, 0) );
400 }
401 }
402 }
403 GIVEN( "Given a degenerated 2d simplex { (2,1,0), (1,0,1), (1,5,1), (0,3,2) } " ) {
404 std::vector<Point> V
405 = { Point(2,1,0), Point(1,0,1), Point(1,5,1), Point(0,3,2) };
406 WHEN( "Computing its lattice polytope" ){
407 const auto P = Helper::computeLatticePolytope( V, false, true );
408 CAPTURE( P );
409 THEN( "The polytope is valid and has more than 6 non trivial facets" ) {
410 REQUIRE( P.nbHalfSpaces() - 6 == 6 );
411 }
412 THEN( "The polytope contains 8 points" ) {
413 REQUIRE( P.count() == 8 );
414 }
415 THEN( "The polytope contains no interior points" ) {
416 REQUIRE( P.countInterior() == 0 );
417 }
418 }
419 WHEN( "Computing the vertices of its convex hull" ){
420 auto X = Helper::computeConvexHullVertices( V, false );
421 std::sort( X.begin(), X.end() );
422 CAPTURE( X );
423 THEN( "The polytope is a quad" ) {
424 REQUIRE( X.size() == 4 );
425 REQUIRE( X[ 0 ] == Point( 0, 3, 2) );
426 REQUIRE( X[ 1 ] == Point( 1, 0, 1) );
427 REQUIRE( X[ 2 ] == Point( 1, 5, 1) );
428 REQUIRE( X[ 3 ] == Point( 2, 1, 0) );
429 }
430 }
431 }
432}
Aim: Implements basic operations that will be used in Point and Vector classes.
Aim: Represents a polygon mesh, i.e. a 2-dimensional combinatorial surface whose faces are (topologic...
SurfaceMesh< RealPoint, RealVector > SMesh
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
Definition SurfaceMesh.h:92
long Euler() const
Size nbFaces() const
Size nbEdges() const
Edges computeManifoldBoundaryEdges() const
Size nbVertices() const
Edges computeNonManifoldEdges() const
PointVector< 3, double > RealPoint

References CAPTURE(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::computeManifoldBoundaryEdges(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::computeNonManifoldEdges(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::Euler(), GIVEN(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbEdges(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbFaces(), DGtal::SurfaceMesh< TRealPoint, TRealVector >::nbVertices(), and REQUIRE().