File failed to load: https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.3/config/TeX-MML-AM_CHTML/MathJax.js
DGtal 2.0.0
QuickHull.h
1
16
17#pragma once
18
30
31#if defined(QuickHull_RECURSES)
32#error Recursive header files inclusion detected in QuickHull.h
33#else // defined(QuickHull_RECURSES)
35#define QuickHull_RECURSES
36
37#if !defined QuickHull_h
39#define QuickHull_h
40
42// Inclusions
43#include <iostream>
44#include <string>
45#include <vector>
46#include <queue>
47#include <set>
48#include "DGtal/base/Common.h"
49#include "DGtal/base/Clock.h"
50#include "DGtal/geometry/tools/QuickHullKernels.h"
51
52namespace DGtal
53{
55 // template class QuickHull
56
138 template < typename TKernel >
140 {
141 typedef TKernel Kernel;
142 typedef typename Kernel::CoordinatePoint Point;
143 typedef typename Kernel::CoordinateVector Vector;
144 typedef typename Kernel::CoordinateScalar Scalar;
145 typedef typename Kernel::InternalScalar InternalScalar;
146 typedef std::size_t Index;
147 typedef std::size_t Size;
149 typedef std::vector< Index > IndexRange;
150 typedef typename Kernel::HalfSpace HalfSpace;
151 typedef typename Kernel::CombinatorialPlaneSimplex CombinatorialPlaneSimplex;
153
155 enum { UNASSIGNED = (Index) -1 };
156
160 struct Facet {
166
167 Facet() = default;
168 Facet( const Facet& ) = default;
169 Facet( Facet&& ) = default;
170 Facet& operator=( Facet&& ) = default;
171 Facet& operator=( const Facet& ) = default;
172 Facet( const HalfSpace& aH, Index b )
173 : H( aH ), below( b ) {}
174
175 void clear()
176 {
177 H = HalfSpace();
178 neighbors.clear();
179 outside_set.clear();
180 on_set.clear();
182 }
184 {
185 const auto it = std::find( on_set.cbegin(), on_set.cend(), p );
186 if ( it == on_set.cend() ) on_set.push_back( p );
187 }
188 void display( std::ostream& out ) const
189 {
190 const auto N = H.internalNormal();
191 out << "[Facet iN=(" << N[0];
192 for ( Dimension i = 1; i < N.dimension; i++ ) out << "," << N[ i ];
193 out << ") c=" << H.internalIntercept() << " b=" << below << " n={";
194 for ( auto&& n : neighbors ) out << " " << n;
195 out << " } #out=" << outside_set.size();
196 out << " on={";
197 for ( auto&& n : on_set ) out << " " << n;
198 out << " }]" << std::endl;
199 }
200
202 {
203 const auto it = std::find( neighbors.cbegin(), neighbors.cend(), n );
204 if ( it == neighbors.cend() ) neighbors.push_back( n );
205 }
207 {
208 auto it = std::find( neighbors.begin(), neighbors.end(), n );
209 if ( it != neighbors.end() ) {
210 std::swap( *it, neighbors.back() );
211 neighbors.pop_back();
212 }
213 }
214 void swap( Facet& other )
215 {
216 if ( this != &other ) {
217 std::swap( H, other.H );
218 neighbors.swap ( other.neighbors );
219 outside_set.swap( other.outside_set );
220 on_set.swap ( other.on_set );
221 std::swap( below, other.below );
222 }
223 }
225 {
226 Size M;
227 M += neighbors.capacity() * sizeof( Index );
228 M += outside_set.capacity() * sizeof( Index );
229 M += on_set.capacity() * sizeof( Index );
230 return M;
231 }
232 };
233
236 typedef std::pair< Index, Index > Ridge;
237
250
251 // ----------------------- standard services --------------------------
252 public:
255
259 QuickHull( const Kernel& K = Kernel(), int dbg = 0 )
261 {}
262
267 { return myStatus; }
268
270 void clear()
271 {
273 points.clear();
274 processed_points.clear();
275 input2comp.clear();
276 comp2input.clear();
277 assignment.clear();
278 facets.clear();
279 deleted_facets.clear();
280 p2v.clear();
281 v2p.clear();
282 timings.clear();
283 }
284
285
287 Size memory() const
288 {
289 // int debug_level;
290 Size M = sizeof( kernel ) + sizeof( int );
291 // std::vector< Point > points;
292 M += sizeof( std::vector< Point > )
293 + points.capacity() * sizeof( Point );
294 M += sizeof( std::vector< Index > )
295 + processed_points.capacity() * sizeof( Index );
296 M += sizeof( std::vector< Index > )
297 + input2comp.capacity() * sizeof( Index );
298 M += sizeof( std::vector< Index > )
299 + comp2input.capacity() * sizeof( Index );
300 // std::vector< Index > assignment;
301 M += sizeof( std::vector< Index > )
302 + assignment.capacity() * sizeof( Index );
303 // std::vector< Facet > facets;
304 M += sizeof( std::vector< Facet > )
305 + facets.capacity() * sizeof( Facet );
306 for ( const auto& f : facets ) M += f.variableMemory();
307 // std::set< Index > deleted_facets;
308 M += sizeof( std::set< Index > )
309 + deleted_facets.size() * ( sizeof( Index ) + 2*sizeof(Index*) );
310 // IndexRange p2v;
311 M += sizeof( std::vector< Index > )
312 + p2v.capacity() * sizeof( Index );
313 // IndexRange v2p;
314 M += sizeof( std::vector< Index > )
315 + v2p.capacity() * sizeof( Index );
316 // std::vector< double > timings;
317 M += sizeof( std::vector< double > )
318 + timings.capacity() * sizeof( double );
319 return M;
320 }
321
323 Size nbPoints() const
324 { return points.size(); }
325
328 Size nbFacets() const
329 {
332 return facets.size();
333 }
334
338 {
341 return v2p.size();
342 }
343
347 {
350 return nb_finite_facets;
351 }
352
356 {
359 return nb_infinite_facets;
360 }
361
363 // -------------------------- Convex hull services ----------------------------
364 public:
367
381 template < typename InputPoint >
382 bool setInput( const std::vector< InputPoint >& input_points,
383 bool remove_duplicates = true )
384 {
385 Clock tic;
386 tic.startClock();
387 clear();
388 timings.clear();
389 kernel.makeInput( points, input2comp, comp2input,
390 input_points, remove_duplicates );
391 timings.push_back( tic.stopClock() );
392 if ( points.size() <= dimension ) {
394 return false;
395 }
397 return true;
398 }
399
410 bool setInitialSimplex( const IndexRange& full_splx )
411 {
412 if ( status() != Status::InputInitialized ) return false;
413 if ( full_splx.size() != dimension + 1 )
414 {
415 trace.error() << "[QuickHull::setInitialSimplex]"
416 << " not a full dimensional simplex" << std::endl;
418 return false;
419 }
421 for ( Index j = 0; j < dimension; ++j )
422 splx[ j ] = full_splx[ j ];
423 const auto H = kernel.compute( points, splx, full_splx.back() );
424 const auto volume = kernel.volume( H, points[ full_splx.back() ] );
425 if ( volume > 0 )
426 return computeSimplexConfiguration( full_splx );
428 return false;
429 }
430
432 // -------------------------- Convex hull services ----------------------------
433 public:
436
461 {
462 if ( target < Status::InputInitialized || target > Status::AllCompleted )
463 return false;
464 Clock tic;
466 { // Initialization
467 tic.startClock();
468 bool ok1 = computeInitialSimplex();
469 timings.push_back( tic.stopClock() );
470 if ( ! ok1 ) return false;
471 if ( status() == target ) return true;
472 }
474 { // Computes facets
475 tic.startClock();
476 bool ok2 = computeFacets();
477 timings.push_back( tic.stopClock() );
478 if ( ! ok2 ) return false;
479 if ( status() == target ) return true;
480 }
482 { // Computes vertices
483 tic.startClock();
484 bool ok3 = computeVertices();
485 timings.push_back( tic.stopClock() );
486 if ( ! ok3 ) return false;
487 if ( status() == target ) return true;
488 }
489 if ( target == Status::AllCompleted
491 { // for now, Status::VerticesCompleted and
492 // Status::AllCompleted are the same.
494 return true;
495 }
496 return false;
497 }
498
505 {
506 const auto full_simplex = pickInitialSimplex();
507 if ( full_simplex.empty() ) {
509 return false;
510 }
511 return computeSimplexConfiguration( full_simplex );
512 }
513
524 {
525 if ( status() != Status::SimplexCompleted ) return false;
526 std::queue< Index > Q;
527 for ( Index fi = 0; fi < facets.size(); ++fi )
528 Q.push( fi );
529 Index n = 0;
530 while ( processFacet( Q ) ) {
531 if ( debug_level >= 1 )
532 trace.info() << "---- Iteration " << n++ << " #Q=" << Q.size() << std::endl;
533 }
534 cleanFacets();
535 if ( debug_level >= 2 ) {
536 trace.info() << ".... #facets=" << facets.size()
537 << " #deleted=" << deleted_facets.size() << std::endl;
538 }
540 return true;
541 }
542
554 {
555 static const int MAX_NB_VPF = 10 * dimension;
556 if ( status() != Status::FacetsCompleted ) return false;
557
558 // Renumber infinite facets in case of Delaunay triangulation computation.
560
561 // Builds the maps v2p: vertex -> point, and p2v : point -> vertex.
562 facet_counter = IndexRange( MAX_NB_VPF, 0 );
563 v2p.clear();
564 p2v.resize( points.size() );
565 std::vector< IndexRange > p2f( points.size() );
566 for ( Index f = 0; f < facets.size(); ++f ) {
567 for ( auto&& p : facets[ f ].on_set ) p2f[ p ].push_back( f );
568 }
569
570 // vertices belong to at least d facets
571 Index v = 0;
572 for ( Index p = 0; p < points.size(); ++p ) {
573 const auto nbf = p2f[ p ].size();
574 facet_counter[ std::min( (int) nbf, MAX_NB_VPF-1 ) ] += 1;
575 if ( nbf >= dimension ) {
576 v2p.push_back( p );
577 p2v[ p ] = v++;
578 }
579 else p2v[ p ] = UNASSIGNED;
580 }
581
582 // Display debug informations
583 if ( debug_level >= 1 ) {
584 trace.info() << "#vertices=" << v2p.size() << " #facets=" << facets.size()
585 << std::endl;
586 trace.info() << "#inc_facets/point= ";
587 for ( auto n : facet_counter ) trace.info() << n << " ";
588 trace.info() << std::endl;
589 }
591 return true;
592 }
593
594
596 // -------------------------- Output services ----------------------------
597 public:
600
611 template < typename OutputPoint >
612 bool getVertexPositions( std::vector< OutputPoint >& vertex_positions )
613 {
614 vertex_positions.clear();
616 && status() <= Status::AllCompleted ) ) return false;
617 vertex_positions.resize( v2p.size() );
618 for ( Index i = 0; i < v2p.size(); i++ ) {
619 kernel.convertPointTo( points[ v2p[ i ] ], vertex_positions[ i ] );
620 }
621 return true;
622 }
623
632 bool getVertex2Point( IndexRange& vertex_to_point )
633 {
634 vertex_to_point.clear();
636 && status() <= Status::AllCompleted ) ) return false;
637 vertex_to_point = v2p;
638 return true;
639 }
640
651 bool getPoint2Vertex( IndexRange& point_to_vertex )
652 {
653 point_to_vertex.clear();
655 && status() <= Status::AllCompleted ) ) return false;
656 point_to_vertex = p2v;
657 return true;
658 }
659
666 bool getFacetVertices( std::vector< IndexRange >& facet_vertices ) const
667 {
668 facet_vertices.clear();
670 && status() <= Status::AllCompleted ) ) return false;
671 facet_vertices.reserve( nbFacets() );
672 for ( Index f = 0; f < nbFacets(); ++f ) {
673 IndexRange ofacet = orientedFacetPoints( f );
674 for ( auto& v : ofacet ) v = p2v[ v ];
675 facet_vertices.push_back( ofacet );
676 }
677 return true;
678 }
679
690 bool getFacetHalfSpaces( std::vector< HalfSpace >& facet_halfspaces )
691 {
692 facet_halfspaces.clear();
693 if ( ! ( status() >= Status::FacetsCompleted
694 && status() <= Status::AllCompleted ) ) return false;
695 facet_halfspaces.reserve( nbFacets() );
696 for ( Index f = 0; f < nbFacets(); ++f ) {
697 facet_halfspaces.push_back( facets[ f ].H );
698 }
699 return true;
700 }
701
702
704 // -------------------------- Check hull services ----------------------------
705 public:
708
715 bool check()
716 {
717 bool ok = true;
719 trace.warning() << "[Quickhull::check] invalid status="
720 << (int)status() << std::endl;
721 ok = false;
722 }
723 if ( processed_points.size() != points.size() ) {
724 trace.warning() << "[Quickhull::check] not all points processed: "
725 << processed_points.size() << " / " << points.size()
726 << std::endl;
727 ok = false;
728 }
729 if ( ! checkHull() ) {
730 trace.warning() << "[Quickhull::check] Check hull is invalid. "
731 << std::endl;
732 ok = false;
733 }
734 if ( ! checkFacets() ) {
735 trace.warning() << "[Quickhull::check] Check facets is invalid. "
736 << std::endl;
737 ok = false;
738 }
739 return ok;
740 }
741
744 {
745 Size nb = 0;
746 Size nbok = 0;
747 for ( Index f = 0; f < facets.size(); ++f )
748 if ( deleted_facets.count( f ) == 0 ) {
749 bool ok = checkFacet( f );
750 nbok += ok ? 1 : 0;
751 nb += 1;
752 }
753 return nb == nbok;
754 }
755
763 {
764 Index nbok = 0;
765 // for ( Index v = 0; v < points.size(); ++v ) {
766 for ( auto v : processed_points ) {
767 bool ok = true;
768 for ( Index f = 0; f < facets.size(); ++f )
769 if ( deleted_facets.count( f ) == 0 ) {
770 if ( above( facets[ f ], points[ v ] ) ) {
771 ok = false;
772 trace.error() << "- bad vertex " << v << " " << points[ v ]
773 << " dist="
774 << ( kernel.height( facets[ f ].H, points[ v ] ) )
775 << std::endl;
776 break;
777 }
778 }
779 nbok += ok ? 1 : 0;
780 }
781 if ( debug_level >= 2 ) {
782 trace.info() << nbok << "/"
783 << processed_points.size() << " vertices inside convex hull."
784 << std::endl;
785 }
787 return nbok == processed_points.size();
788 }
789
791 // ------------------------ public datas --------------------------
792 public:
795
796 public:
798 mutable Kernel kernel;
802 std::vector< Point > points;
812 std::vector< Index > assignment;
814 std::vector< Facet > facets;
816 std::set< Index > deleted_facets;
826 std::vector< double > timings;
828 std::vector< Size > facet_counter;
829
831 // ------------------------ protected datas --------------------------
832 protected:
835
840
842 // --------------------- protected services --------------------------
843 protected:
846
850 InternalScalar height( const Facet& F, const Point& p ) const
851 { return kernel.height( F.H, p ); }
852
856 bool above( const Facet& F, const Point& p ) const
857 { return kernel.above( F.H, p ); }
858
862 bool aboveOrOn( const Facet& F, const Point& p ) const
863 { return kernel.aboveOrOn( F.H, p ); }
864
868 bool on( const Facet& F, const Point& p ) const
869 { return kernel.on( F.H, p ); }
870
874 {
875 if ( deleted_facets.empty() ) return;
876 IndexRange renumbering( facets.size() );
877 Index i = 0;
878 Index j = 0;
879 for ( auto& l : renumbering ) {
880 if ( ! deleted_facets.count( j ) ) l = i++;
881 else l = UNASSIGNED;
882 j++;
883 }
884 const Index nf = facets.size() - deleted_facets.size();
885 deleted_facets.clear();
886 for ( Index f = 0; f < facets.size(); f++ )
887 if ( ( renumbering[ f ] != UNASSIGNED ) && ( f != renumbering[ f ] ) )
888 facets[ renumbering[ f ] ] = facets[ f ];
889 facets.resize( nf );
890 for ( auto& F : facets ) {
891 for ( auto& N : F.neighbors ) {
892 if ( renumbering[ N ] == UNASSIGNED )
893 trace.error() << "Invalid deleted neighboring facet." << std::endl;
894 else N = renumbering[ N ];
895 }
896 }
897 }
898
902 {
903 nb_finite_facets = facets.size();
905 if ( ! kernel.hasInfiniteFacets() ) return;
906 IndexRange renumbering( facets.size() );
907 Index i = 0;
908 Index k = facets.size();
909 Index j = 0;
910 for ( auto& l : renumbering ) {
911 if ( ! kernel.isHalfSpaceFacetInfinite( facets[ j ].H ) ) l = i++;
912 else l = --k;
913 j++;
914 }
915 if ( i != k )
916 trace.error() << "[Quickhull::renumberInfiniteFacets]"
917 << " Error renumbering infinite facets "
918 << " up finite=" << i << " low infinite=" << k << std::endl;
919 std::vector< Facet > new_facets( facets.size() );
920 for ( Index f = 0; f < facets.size(); f++ )
921 new_facets[ renumbering[ f ] ].swap( facets[ f ] );
922 facets.swap( new_facets );
923 for ( auto& F : facets ) {
924 for ( auto& N : F.neighbors ) {
925 N = renumbering[ N ];
926 }
927 }
928 // Assign correct number of facets.
931 }
932
933
944 bool processFacet( std::queue< Index >& Q )
945 {
946 // If Q empty, we are done
947 if ( Q.empty() ) return false;
948 Index F = Q.front();
949 Q.pop();
950 // If F is already deleted, proceed to next in queue.
951 if ( deleted_facets.count( F ) ) return true;
952 // Take car of current facet.
953 const Facet& facet = facets[ F ];
954 if ( debug_level >= 3 ) {
955 trace.info() << "---------------------------------------------" << std::endl;
956 trace.info() << "---- ACTIVE FACETS---------------------------" << std::endl;
957 bool ok = true;
958 for ( Index i = 0; i < facets.size(); i++ )
959 if ( ! deleted_facets.count( i ) ) {
960 trace.info() << "- facet " << i << " ";
961 facets[ i ].display( trace.info() );
962 ok = ok && checkFacet( i );
963 }
964 if ( ! ok ) { // should not happen.
966 return false;
967 }
968 }
969 if ( debug_level >= 2 ) {
970 trace.info() << "---------------------------------------------" << std::endl;
971 trace.info() << "Processing facet " << F << " ";
972 facet.display( trace.info() );
973 }
974 if ( facet.outside_set.empty() ) return true;
975 // Selects furthest vertex
976 Index furthest_v = facet.outside_set[ 0 ];
977 auto furthest_h = height( facet, points[ furthest_v ] );
978 for ( Index v = 1; v < facet.outside_set.size(); v++ ) {
979 auto h = height( facet, points[ facet.outside_set[ v ] ] );
980 if ( h > furthest_h ) {
981 furthest_h = h;
982 furthest_v = facet.outside_set[ v ];
983 }
984 }
985 const Point& p = points[ furthest_v ];
986 // Extracts Visible facets V and Horizon Ridges H
987 std::vector< Index > V; // visible facets
988 std::set< Index > M; // marked facets (are in E or were in E)
989 std::queue< Index > E; // queue to extract visible facets
990 std::vector< Ridge > H; // visible facets
991 E.push ( F );
992 M.insert( F );
993 while ( ! E.empty() ) {
994 Index G = E.front(); E.pop();
995 V.push_back( G );
996 for ( auto& N : facets[ G ].neighbors ) {
997 if ( aboveOrOn( facets[ N ], p ) ) {
998 if ( M.count( N ) ) continue;
999 E.push( N );
1000 } else {
1001 H.push_back( { G, N } );
1002 }
1003 M.insert( N );
1004 }
1005 } // while ( ! E.empty() )
1006 if ( debug_level >= 1 ) {
1007 trace.info() << "#Visible=" << V.size() << " #Horizon=" << H.size()
1008 << " furthest_v=" << furthest_v << std::endl;
1009 }
1010 // Create new facets
1011 IndexRange new_facets;
1012 // For each ridge R in H
1013 for ( Index i = 0; i < H.size(); i++ )
1014 {
1015 // Create a new facet from R and p
1016 IndexRange ridge = pointsOnRidge( H[ i ] );
1017 if ( debug_level >= 3 ) {
1018 trace.info() << "Ridge (" << H[i].first << "," << H[i].second << ") = {";
1019 for ( auto&& r : ridge ) trace.info() << " " << r;
1020 trace.info() << " } furthest_v=" << furthest_v << std::endl;
1021 }
1022 IndexRange base( 1 + ridge.size() );
1023 Index j = 0;
1024 base[ j++ ] = furthest_v;
1025 for ( auto&& v : ridge ) base[ j++ ] = v;
1026 if ( j < dimension ) {
1027 trace.error() << "Bad ridge between " << std::endl
1028 << "- facet " << H[i].first << " ";
1029 facets[ H[i].first ].display( trace.error() );
1030 trace.error() << "- facet " << H[i].second << " ";
1031 facets[ H[i].second ].display( trace.error() );
1032 }
1033 Index nf = newFacet();
1034 new_facets.push_back( nf );
1035 facets[ nf ] = makeFacet( base, facets[ H[i].first ].below );
1036 facets[ nf ].on_set = IndexRange { base.cbegin(), base.cend() };
1037 std::sort( facets[ nf ].on_set.begin(), facets[ nf ].on_set.end() );
1038 makeNeighbors( nf, H[ i ].second );
1039 if ( debug_level >= 3 ) {
1040 trace.info() << "* New facet " << nf << " ";
1041 facets[ nf ].display( trace.info() );
1042 }
1043 // Checks that the facet is not parallel to another in the Horizon
1044 for ( Index k = 0; k < new_facets.size() - 1; k++ )
1045 if ( areFacetsParallel( new_facets[ k ], nf ) ) {
1046 if ( debug_level >= 1 ) {
1047 trace.info() << "Facets " << new_facets[ k ] << " and " << nf
1048 << " are parallel => merge." << std::endl;
1049 }
1050 mergeParallelFacets( new_facets[ k ], nf );
1051 new_facets.pop_back();
1052 deleteFacet( nf );
1053 if ( debug_level >= 3 ) {
1054 facets[ new_facets[ k ] ].display( trace.info() );
1055 }
1056 }
1057 }
1058 // For each new facet
1059 for ( Index i = 0; i < new_facets.size(); i++ )
1060 { // link the new facet to its neighbors
1061 for ( Index j = i + 1; j < new_facets.size(); j++ )
1062 {
1063 const Index nfi = new_facets[ i ];
1064 const Index nfj = new_facets[ j ];
1065 if ( areFacetsNeighbor( nfi, nfj ) )
1066 makeNeighbors( nfi, nfj );
1067 }
1068 }
1069 // Extracts all outside points from visible facets V
1070 IndexRange outside_pts;
1071 for ( auto&& vf : V ) {
1072 for ( auto&& v : facets[ vf ].outside_set ) {
1073 if ( v != furthest_v ) {
1074 outside_pts.push_back( v );
1075 assignment[ v ] = UNASSIGNED;
1076 }
1077 }
1078 }
1079 // For each new facet F'
1080 for ( Index i = 0; i < new_facets.size(); i++ ) {
1081 Facet& Fp = facets[ new_facets[ i ] ];
1082 Index max_j = outside_pts.size();
1083 for ( Index j = 0; j < max_j; ) {
1084 const Index v = outside_pts[ j ];
1085 if ( above( Fp, points[ v ] ) ) {
1086 Fp.outside_set.push_back( v );
1087 assignment[ v ] = new_facets[ i ];
1088 outside_pts[ j ] = outside_pts.back();
1089 outside_pts.pop_back();
1090 max_j--;
1091 } else j++;
1092 }
1093 if ( debug_level >= 3 ) {
1094 trace.info() << "- New facet " << new_facets[ i ] << " ";
1095 Fp.display( trace.info() );
1096 }
1097 }
1098 // Update processed points
1099 processed_points.push_back( furthest_v );
1100 for ( auto v : outside_pts ) processed_points.push_back( v );
1101
1102 // Delete the facets in V
1103 for ( auto&& v : V ) {
1104 if ( debug_level >= 2 ) {
1105 trace.info() << "Delete facet " << v << " ";
1106 facets[ v ].display( trace.info() );
1107 }
1108 deleteFacet( v );
1109 }
1110
1111 // Add new facets to queue
1112 for ( Index i = 0; i < new_facets.size(); i++ )
1113 Q.push( new_facets[ i ] );
1114 if ( debug_level >= 1 ) {
1115 trace.info() << "#facets=" << facets.size()
1116 << " #deleted=" << deleted_facets.size() << std::endl;
1117 }
1118
1119 // Checks that everything is ok.
1120 if ( debug_level >= 1 ) {
1121 trace.info() << "[CHECK INVARIANT] " << processed_points.size()
1122 << " / " << points.size() << " points processed." << std::endl;
1123 bool okh = checkHull();
1124 if ( ! okh )
1125 trace.error() << "[computeFacet] Invalid convex hull" << std::endl;
1126 bool okf = checkFacets();
1127 if ( ! okf )
1128 trace.error() << "[computeFacet] Invalid facets" << std::endl;
1129 if ( ! ( okh && okf ) ) myStatus = Status::InvalidConvexHull;
1130 }
1131
1132 return status() == Status::SimplexCompleted;
1133 }
1134
1136 bool checkFacet( Index f ) const
1137 {
1138 const Facet& F = facets[ f ];
1139 bool ok = F.on_set.size() >= dimension;
1140 for ( auto v : F.on_set )
1141 if ( ! on( F, points[ v ] ) ) {
1142 trace.error() << "[QuickHull::checkFacet( " << f
1143 << ") Invalid 'on' vertex " << v << std::endl;
1144 ok = false;
1145 }
1146 if ( F.neighbors.size() < dimension ) {
1147 trace.error() << "[QuickHull::checkFacet( " << f
1148 << ") Not enough neighbors " << F.neighbors.size() << std::endl;
1149 ok = false;
1150 }
1151 for ( auto nf : F.neighbors )
1152 if ( ! areFacetsNeighbor( f, nf ) ) {
1153 trace.error() << "[QuickHull::checkFacet( " << f
1154 << ") Invalid neighbor " << nf << std::endl;
1155 ok = false;
1156 }
1157 if ( aboveOrOn( F, points[ F.below ] ) ) {
1158 trace.error() << "[QuickHull::checkFacet( " << f
1159 << ") Bad below point " << F.below << std::endl;
1160 ok = false;
1161 }
1162 for ( auto ov : F.outside_set )
1163 if ( ! above( F, points[ ov ] ) ) {
1164 trace.error() << "[QuickHull::checkFacet( " << f
1165 << ") Bad outside point " << ov << std::endl;
1166 ok = false;
1167 }
1168 for ( auto && v : F.on_set ) {
1169 Size n = 0;
1170 for ( auto&& N : facets[ f ].neighbors )
1171 if ( on( facets[ N ], points[ v ] ) ) n += 1;
1172 if ( n < dimension-1 ) {
1173 trace.error() << "[QuickHull::checkFacet( " << f << ") 'on' point " << v
1174 << " is a vertex of " << n << " facets "
1175 << "(should be >= " << dimension-1 << ")" << std::endl;
1176 ok = false;
1177 }
1178 }
1179 return ok;
1180 }
1181
1184 {
1185 // SLightly faster to postpone deletion of intermediate facets.
1186 const Index f = facets.size();
1187 facets.push_back( Facet() );
1188 return f;
1189 }
1190
1194 {
1195 for ( auto n : facets[ f ].neighbors )
1196 facets[ n ].subNeighbor( f );
1197 deleted_facets.insert( f );
1198 facets[ f ].clear();
1199 }
1200
1204 void makeNeighbors( const Index if1, const Index if2 )
1205 {
1206 facets[ if1 ].addNeighbor( if2 );
1207 facets[ if2 ].addNeighbor( if1 );
1208 }
1209
1213 void unmakeNeighbors( const Index if1, const Index if2 )
1214 {
1215 facets[ if1 ].subNeighbor( if2 );
1216 facets[ if2 ].subNeighbor( if1 );
1217 }
1218
1225 Index mergeParallelFacets( const Index if1, const Index if2 )
1226 {
1227 Facet& f1 = facets[ if1 ];
1228 Facet& f2 = facets[ if2 ];
1229 std::copy( f2.outside_set.cbegin(), f2.outside_set.cend(),
1230 std::back_inserter( f1.outside_set ) );
1231 IndexRange merge_idx;
1232 std::set_union( f1.on_set.cbegin(), f1.on_set.cend(),
1233 f2.on_set.cbegin(), f2.on_set.cend(),
1234 std::back_inserter( merge_idx ) );
1235 f1.on_set.swap( merge_idx );
1236 for ( auto && nf2 : f2.neighbors ) {
1237 if ( nf2 == if1 ) continue;
1238 facets[ nf2 ].subNeighbor( if2 );
1239 makeNeighbors( if1, nf2 );
1240 }
1241 return if2;
1242 }
1243
1248 bool areFacetsParallel( const Index if1, const Index if2 ) const
1249 {
1250 ASSERT( if1 != if2 );
1251 const Facet& f1 = facets[ if1 ];
1252 const Facet& f2 = facets[ if2 ];
1253 if ( kernel.equal( f1.H, f2.H ) ) return true;
1254 // Need to check if one N is a multiple of the other.
1255 for ( auto&& v : f1.on_set )
1256 if ( ! on( f2, points[ v ] ) ) return false;
1257 return true;
1258 }
1259
1267 bool areFacetsNeighbor( const Index if1, const Index if2 ) const
1268 {
1269 const Facet& f1 = facets[ if1 ];
1270 const Facet& f2 = facets[ if2 ];
1271 Index nb = 0;
1272 for ( Index i1 = 0, i2 = 0; i1 < f1.on_set.size() && i2 < f2.on_set.size(); )
1273 {
1274 if ( f1.on_set[ i1 ] == f2.on_set[ i2 ] ) { nb++; i1++; i2++; }
1275 else if ( f1.on_set[ i1 ] < f2.on_set[ i2 ] ) i1++;
1276 else i2++;
1277 }
1278 return nb >= ( dimension - 1 );
1279 }
1280
1292 Facet makeFacet( const IndexRange& base, Index below ) const
1293 {
1295 for ( Size i = 0; i < dimension; i++ ) simplex[ i ] = base[ i ];
1296 auto plane = kernel.compute( points, simplex, below );
1297 return Facet( plane, below );
1298 }
1299
1303 {
1304 // Slightly faster to look at points instead at true geometry.
1305 IndexRange result;
1306 const Facet& f1 = facets[ R.first ];
1307 const Facet& f2 = facets[ R.second ];
1308 std::set_intersection( f1.on_set.cbegin(), f1.on_set.cend(),
1309 f2.on_set.cbegin(), f2.on_set.cend(),
1310 std::back_inserter( result ) );
1311 return result;
1312 }
1313
1324 {
1325 const Facet& F = facets[ f ];
1326 IndexRange result = F.on_set;
1327 // Sort a facet such that its points, taken in order, have
1328 // always the same orientation of the facet. More precisely,
1329 // facets span a `dimension-1` vector space. There are thus
1330 // dimension-2 fixed points, and the last ones (at least two)
1331 // may be reordered.
1333 for ( Dimension k = 0; k < dimension-2; k++ )
1334 splx[ k ] = result[ k ];
1335 std::sort( result.begin()+dimension-2, result.end(),
1336 [&] ( Index i, Index j )
1337 {
1338 splx[ dimension-2 ] = i;
1339 splx[ dimension-1 ] = j;
1340 const auto H = kernel.compute( points, splx );
1341 const auto orient = kernel.dot( F.H, H );
1342 return orient > 0;
1343 } );
1344 return result;
1345 }
1346
1347
1353 {
1354 auto & on_set = facets[ f ].on_set;
1355 for ( Index i = 0; i < on_set.size(); )
1356 {
1357 Index v = on_set[ i ];
1358 Size n = 0;
1359 for ( auto&& N : facets[ f ].neighbors )
1360 if ( on( facets[ N ], points[ v ] ) ) n += 1;
1361 if ( n >= dimension-1 ) i++;
1362 else {
1363 on_set[ i ] = on_set.back();
1364 on_set.pop_back();
1365 }
1366 }
1367 std::sort( on_set.begin(), on_set.end() );
1368 }
1369
1373 {
1374 const Size nb = points.size();
1375 if ( nb < dimension + 1 ) return IndexRange();
1376 IndexRange best = pickIntegers( dimension + 1, nb );
1378 for ( Index j = 0; j < dimension; ++j ) splx[ j ] = best[ j ];
1379 const auto first_H = kernel.compute( points, splx, best.back() );
1380 auto best_volume = kernel.volume ( first_H, points[ best.back() ] );
1381 const Size nbtries = std::min( (Size) 10, 1 + nb / 10 );
1382 const Size max_nbtries = std::max( (Size) 10, 2 * nb );
1383 for ( Size i = 0; i < max_nbtries; i++ )
1384 {
1385 IndexRange tmp = pickIntegers( dimension + 1, nb );
1386 for ( Index j = 0; j < dimension; ++j ) splx[ j ] = tmp[ j ];
1387 const auto tmp_H = kernel.compute( points, splx, tmp.back() );
1388 const auto tmp_volume = kernel.volume ( tmp_H, points[ tmp.back() ] );
1389 if ( best_volume < tmp_volume ) {
1390 if ( debug_level >= 1 ) {
1391 trace.info() << "(" << i << ")"
1392 << " new_volume = " << tmp_volume
1393 << " > " << best_volume << std::endl;
1394 }
1395 best = tmp;
1396 best_volume = tmp_volume;
1397 }
1398 if ( i >= nbtries && best_volume > 0 ) return best;
1399 }
1400 return IndexRange();
1401 }
1402
1406 static IndexRange pickIntegers( const Size d, const Size n )
1407 {
1408 IndexRange result( d );
1409 bool distinct = false;
1410 while ( ! distinct )
1411 {
1412 distinct = true;
1413 for ( Index i = 0; i < d; i++ ) result[ i ] = rand() % n;
1414 std::sort( result.begin(), result.end() );
1415 for ( Index i = 1; distinct && i < d; i++ )
1416 distinct = result[ i-1 ] != result[ i ];
1417 }
1418 return result;
1419 }
1420
1429 bool computeSimplexConfiguration( const IndexRange& full_simplex )
1430 {
1431 assignment = std::vector< Index >( points.size(), UNASSIGNED );
1432 facets.resize( dimension + 1 );
1433 deleted_facets.clear();
1434 for ( Index j = 0; j < full_simplex.size(); ++j )
1435 {
1436 IndexRange lsimplex( dimension );
1437 IndexRange isimplex( dimension );
1438 Index s = 0;
1439 for ( Index i = 0; i <= dimension; i++ )
1440 if ( i != j ) {
1441 lsimplex[ s ] = i;
1442 isimplex[ s ] = full_simplex[ i ];
1443 s++;
1444 }
1445 facets[ j ] = makeFacet( isimplex, full_simplex[ j ] );
1446 facets[ j ].neighbors = lsimplex;
1447 for ( auto&& v : isimplex ) facets[ j ].on_set.push_back( v );
1448 std::sort( facets[ j ].on_set.begin(), facets[ j ].on_set.end() );
1449 }
1450 // List of unassigned vertices
1451 for ( Index fi = 0; fi < facets.size(); ++fi ) {
1452 Facet& f = facets[ fi ];
1453 for ( Index v = 0; v < points.size(); v++ )
1454 {
1455 if ( assignment[ v ] == UNASSIGNED && above( f, points[ v ] ) ) {
1456 f.outside_set.push_back( v );
1457 assignment[ v ] = fi;
1458 }
1459 }
1460 }
1461 for ( Index v = 0; v < points.size(); v++ )
1462 if ( assignment[ v ] == UNASSIGNED )
1463 processed_points.push_back( v );
1464
1465 // Display some information
1466 if ( debug_level >= 2 ) {
1467 for ( auto&& f : facets ) f.display( trace.info() );
1468 }
1470 if ( debug_level >= 1 ) {
1471 trace.info() << "[CHECK INVARIANT] " << processed_points.size()
1472 << " / " << points.size() << " points processed." << std::endl;
1473 bool okh = checkHull();
1474 if ( ! okh )
1475 trace.error() << "[computeInitialSimplex] Invalid convex hull" << std::endl;
1476 bool okf = checkFacets();
1477 if ( ! okf )
1478 trace.error() << "[computeInitialSimplex] Invalid facets" << std::endl;
1479 if ( ! ( okh && okf ) ) myStatus = Status::InvalidConvexHull;
1480 }
1481 return status() == Status::SimplexCompleted;
1482 }
1483
1485 // ----------------------- Interface --------------------------------------
1486 public:
1489
1494 void selfDisplay ( std::ostream & out ) const
1495 {
1496 out << "[QuickHull " << dimension << "D"
1497 // << " status=" << status()
1498 << " #P=" << nbPoints();
1500 out << " #F=" << nbFacets();
1502 out << " #V=" << nbVertices();
1503 out << "]";
1505 && nbFacets() == 24 ) {
1506 for ( auto f : facets ) f.display( out );
1507 for ( auto v : v2p ) out << points[ v2p[ v ] ] << std::endl;
1508 }
1509 }
1510
1515 bool isValid() const
1516 {
1517 return status() >= Status::Uninitialized
1519 }
1520
1521
1522 };
1523
1531 template < typename TKernel >
1532 std::ostream&
1533 operator<< ( std::ostream & out,
1534 const QuickHull< TKernel > & object )
1535 {
1536 object.selfDisplay( out );
1537 return out;
1538 }
1539
1540} // namespace DGtal
1541
1543// Includes inline functions.
1544// //
1546
1547#endif // !defined QuickHull_h
1548
1549#undef QuickHull_RECURSES
1550#endif // else defined(QuickHull_RECURSES)
1551
void tic()
Starts timer.
DGtal is the top-level namespace which contains all DGtal functions and types.
void swap(UnorderedSetByBlock< Key, TSplitter, Hash, KeyEqual, UnorderedMapAllocator > &s1, UnorderedSetByBlock< Key, TSplitter, Hash, KeyEqual, UnorderedMapAllocator > &s2) noexcept
std::ostream & operator<<(std::ostream &out, const ClosedIntegerHalfPlane< TSpace > &object)
DGtal::uint32_t Dimension
Definition Common.h:119
Trace trace
void subNeighbor(Index n)
Definition QuickHull.h:206
Facet & operator=(Facet &&)=default
void addPointOn(Index p)
Definition QuickHull.h:183
Facet(const HalfSpace &aH, Index b)
Definition QuickHull.h:172
void display(std::ostream &out) const
Definition QuickHull.h:188
Facet(Facet &&)=default
void addNeighbor(Index n)
Definition QuickHull.h:201
Index below
index of point that is below this facet
Definition QuickHull.h:165
Facet(const Facet &)=default
HalfSpace H
the facet geometry
Definition QuickHull.h:161
IndexRange outside_set
outside set, i.e. points above this facet
Definition QuickHull.h:163
void swap(Facet &other)
Definition QuickHull.h:214
IndexRange neighbors
neighbor facets
Definition QuickHull.h:162
Facet & operator=(const Facet &)=default
IndexRange on_set
on set, i.e. points on this facet, sorted
Definition QuickHull.h:164
Size variableMemory() const
Definition QuickHull.h:224
Aim: Implements the quickhull algorithm by Barber et al. barber1996, a famous arbitrary dimensional c...
Definition QuickHull.h:140
bool computeInitialSimplex()
Definition QuickHull.h:504
bool computeFacets()
Definition QuickHull.h:523
void unmakeNeighbors(const Index if1, const Index if2)
Definition QuickHull.h:1213
bool setInput(const std::vector< InputPoint > &input_points, bool remove_duplicates=true)
Definition QuickHull.h:382
Kernel::CoordinateVector Vector
Definition QuickHull.h:143
bool areFacetsNeighbor(const Index if1, const Index if2) const
Definition QuickHull.h:1267
QuickHull(const Kernel &K=Kernel(), int dbg=0)
Definition QuickHull.h:259
bool setInitialSimplex(const IndexRange &full_splx)
Definition QuickHull.h:410
Kernel::InternalScalar InternalScalar
Definition QuickHull.h:145
std::vector< Index > IndexRange
Definition QuickHull.h:149
bool on(const Facet &F, const Point &p) const
Definition QuickHull.h:868
void filterVerticesOnFacet(const Index f)
Definition QuickHull.h:1352
bool checkFacet(Index f) const
Definition QuickHull.h:1136
BOOST_STATIC_ASSERT((Point::dimension==Vector::dimension))
Facet makeFacet(const IndexRange &base, Index below) const
Definition QuickHull.h:1292
IndexRange pointsOnRidge(const Ridge &R) const
Definition QuickHull.h:1302
Kernel::CombinatorialPlaneSimplex CombinatorialPlaneSimplex
Definition QuickHull.h:151
bool computeSimplexConfiguration(const IndexRange &full_simplex)
Definition QuickHull.h:1429
std::vector< Facet > facets
Definition QuickHull.h:814
bool getFacetVertices(std::vector< IndexRange > &facet_vertices) const
Definition QuickHull.h:666
std::size_t Index
Definition QuickHull.h:146
std::pair< Index, Index > Ridge
Definition QuickHull.h:236
Kernel::CoordinatePoint Point
Definition QuickHull.h:142
std::vector< Point > points
Definition QuickHull.h:802
std::set< Index > deleted_facets
Definition QuickHull.h:816
void clear()
Clears the object.
Definition QuickHull.h:270
Size nbFiniteFacets() const
Definition QuickHull.h:346
bool areFacetsParallel(const Index if1, const Index if2) const
Definition QuickHull.h:1248
bool computeVertices()
Definition QuickHull.h:553
bool getVertex2Point(IndexRange &vertex_to_point)
Definition QuickHull.h:632
IndexRange pickInitialSimplex() const
Definition QuickHull.h:1372
bool isValid() const
Definition QuickHull.h:1515
IndexRange orientedFacetPoints(Index f) const
Definition QuickHull.h:1323
std::vector< Size > facet_counter
Definition QuickHull.h:828
Status status() const
Definition QuickHull.h:266
bool processFacet(std::queue< Index > &Q)
Definition QuickHull.h:944
static IndexRange pickIntegers(const Size d, const Size n)
Definition QuickHull.h:1406
Size nbVertices() const
Definition QuickHull.h:337
std::size_t Size
Definition QuickHull.h:147
void deleteFacet(Index f)
Definition QuickHull.h:1193
bool getVertexPositions(std::vector< OutputPoint > &vertex_positions)
Definition QuickHull.h:612
std::vector< Index > assignment
Definition QuickHull.h:812
static const Size dimension
Definition QuickHull.h:152
Status
Represents the status of a QuickHull object.
Definition QuickHull.h:239
@ InvalidRidge
Error: some ridge is not consistent (probably integer overflow).
Definition QuickHull.h:247
@ InvalidConvexHull
Error: the convex hull does not contain all its points (probably integer overflow).
Definition QuickHull.h:248
@ SimplexCompleted
An initial full-dimensional simplex has been found. QuickHull core algorithm can start.
Definition QuickHull.h:242
@ FacetsCompleted
All facets of the convex hull are identified.
Definition QuickHull.h:243
@ NotFullDimensional
Error: the initial simplex is not full dimensional.
Definition QuickHull.h:246
@ VerticesCompleted
All vertices of the convex hull are determined.
Definition QuickHull.h:244
@ InputInitialized
A range of input points has been given to QuickHull.
Definition QuickHull.h:241
@ AllCompleted
Same as VerticesCompleted.
Definition QuickHull.h:245
@ Uninitialized
QuickHull is empty and has just been instantiated.
Definition QuickHull.h:240
Size nbInfiniteFacets() const
Definition QuickHull.h:355
InternalScalar height(const Facet &F, const Point &p) const
Definition QuickHull.h:850
bool getFacetHalfSpaces(std::vector< HalfSpace > &facet_halfspaces)
Definition QuickHull.h:690
Index mergeParallelFacets(const Index if1, const Index if2)
Definition QuickHull.h:1225
void makeNeighbors(const Index if1, const Index if2)
Definition QuickHull.h:1204
void selfDisplay(std::ostream &out) const
Definition QuickHull.h:1494
bool aboveOrOn(const Facet &F, const Point &p) const
Definition QuickHull.h:862
Size memory() const
Definition QuickHull.h:287
Size nbFacets() const
Definition QuickHull.h:328
bool getPoint2Vertex(IndexRange &point_to_vertex)
Definition QuickHull.h:651
void renumberInfiniteFacets()
Definition QuickHull.h:901
bool computeConvexHull(Status target=Status::VerticesCompleted)
Definition QuickHull.h:460
Kernel::HalfSpace HalfSpace
Definition QuickHull.h:150
Kernel::CoordinateScalar Scalar
Definition QuickHull.h:144
std::vector< double > timings
Definition QuickHull.h:826
bool above(const Facet &F, const Point &p) const
Definition QuickHull.h:856
Size nbPoints() const
Definition QuickHull.h:323
KSpace K