DGtal 2.1.0
Loading...
Searching...
No Matches
DGtal::QuickHull< TKernel > Struct Template Reference

Aim: Implements the quickhull algorithm by Barber et al. [9], a famous arbitrary dimensional convex hull computation algorithm. It relies on dedicated geometric kernels for computing and comparing facet geometries. More...

#include <DGtal/geometry/tools/QuickHull.h>

Inheritance diagram for DGtal::QuickHull< TKernel >:
[legend]

Data Structures

struct  Facet
 

Public Types

enum  { UNASSIGNED = (Index) -1 }
 Label for points that are not assigned to any facet. More...
 
enum class  Status {
  Uninitialized = 0 , InputInitialized = 1 , SimplexCompleted = 2 , FacetsCompleted = 3 ,
  VerticesCompleted = 4 , AllCompleted = 5 , NotFullDimensional = 10 , InvalidRidge = 11 ,
  InvalidConvexHull = 12
}
 Represents the status of a QuickHull object. More...
 
typedef TKernel Kernel
 
typedef Kernel::CoordinatePoint Point
 
typedef Kernel::CoordinateVector Vector
 
typedef Kernel::CoordinateScalar Scalar
 
typedef Kernel::InternalScalar InternalScalar
 
typedef std::size_t Index
 
typedef std::size_t Size
 
typedef std::vector< IndexIndexRange
 
typedef Kernel::HalfSpace HalfSpace
 
typedef Kernel::CombinatorialPlaneSimplex CombinatorialPlaneSimplex
 
typedef std::pair< Index, IndexRidge
 

Public Member Functions

 BOOST_STATIC_ASSERT ((Point::dimension==Vector::dimension))
 
Standard services (construction, initialization, accessors)
 QuickHull (const Kernel &K_=Kernel(), int dbg=0)
 
Status status () const
 
void clear ()
 Clears the object.
 
Size memory () const
 
Size nbPoints () const
 
Size nbFacets () const
 
Size nbVertices () const
 
Size nbFiniteFacets () const
 
Size nbInfiniteFacets () const
 
Initialization services
template<typename InputPoint >
bool setInput (const std::vector< InputPoint > &input_points, bool remove_duplicates=true)
 
bool setInitialSimplex (const IndexRange &full_splx)
 
Convex hull services
bool computeConvexHull (Status target=Status::VerticesCompleted)
 
bool computeInitialSimplex ()
 
bool computeFacets ()
 
bool computeVertices ()
 
Output services
template<typename OutputPoint >
bool getVertexPositions (std::vector< OutputPoint > &vertex_positions)
 
bool getVertex2Point (IndexRange &vertex_to_point)
 
bool getPoint2Vertex (IndexRange &point_to_vertex)
 
bool getFacetVertices (std::vector< IndexRange > &facet_vertices) const
 
bool getFacetHalfSpaces (std::vector< HalfSpace > &facet_halfspaces)
 
Check hull services
bool check ()
 
bool checkFacets ()
 
bool checkHull ()
 
Interface
void selfDisplay (std::ostream &out) const
 
bool isValid () const
 

Data Fields

public datas
Kernel kernel
 Kernel that is duplicated for computing facet geometries.
 
int debug_level
 debug_level from 0:no to 2
 
std::vector< Pointpoints
 the set of points, indexed as in the array.
 
IndexRange input2comp
 
IndexRange comp2input
 
IndexRange processed_points
 Points already processed (and within the convex hull).
 
std::vector< Indexassignment
 assignment of points to facets
 
std::vector< Facetfacets
 the current set of facets.
 
std::set< Indexdeleted_facets
 set of deleted facets
 
IndexRange p2v
 point index -> vertex index (or UNASSIGNED)
 
IndexRange v2p
 vertex index -> point index
 
Size nb_finite_facets
 Number of finite facets.
 
Size nb_infinite_facets
 Number of infinite facets (!= 0 only for specific kernels)
 
std::vector< double > timings
 Timings of the different phases: 0: init, 1: facets, 2: vertices.
 
std::vector< Sizefacet_counter
 Counts the number of facets with a given number of vertices.
 

Static Public Attributes

static const Size dimension = Point::dimension
 

Protected Attributes

protected datas
Status myStatus
 

protected services

InternalScalar height (const Facet &F, const Point &p) const
 
bool above (const Facet &F, const Point &p) const
 
bool aboveOrOn (const Facet &F, const Point &p) const
 
bool on (const Facet &F, const Point &p) const
 
void cleanFacets ()
 
void renumberInfiniteFacets ()
 
bool processFacet (std::queue< Index > &Q)
 
bool checkFacet (Index f) const
 
Index newFacet ()
 
void deleteFacet (Index f)
 
void makeNeighbors (const Index if1, const Index if2)
 
void unmakeNeighbors (const Index if1, const Index if2)
 
Index mergeParallelFacets (const Index if1, const Index if2)
 
bool areFacetsParallel (const Index if1, const Index if2) const
 
bool areFacetsNeighbor (const Index if1, const Index if2) const
 
Facet makeFacet (const IndexRange &base, Index below) const
 
IndexRange pointsOnRidge (const Ridge &R) const
 
IndexRange orientedFacetPoints (Index f) const
 
void filterVerticesOnFacet (const Index f)
 
IndexRange pickInitialSimplex () const
 
bool computeSimplexConfiguration (const IndexRange &full_simplex)
 
static IndexRange pickIntegers (const Size d, const Size n)
 

Detailed Description

template<typename TKernel>
struct DGtal::QuickHull< TKernel >

Aim: Implements the quickhull algorithm by Barber et al. [9], a famous arbitrary dimensional convex hull computation algorithm. It relies on dedicated geometric kernels for computing and comparing facet geometries.

Description of template class 'QuickHull'

You can use it to build convex hulls of points with integral coordinate (using kernel ConvexHullIntegralKernel), convex hulls with rational coordinates (using kernel ConvexHullRationalKernel) or to build Delaunay triangulations (using kernels DelaunayIntegralKernel and DelaunayRationalKernel).

See also
QuickHull algorithm in arbitrary dimension for convex hull and Delaunay cell complex computation

Below is a complete example that computes the convex hull of points randomly defined in a ball, builds a 3D mesh out of it and output it as an OBJ file.

#include "DGtal/base/Common.h"
#include "DGtal/kernel/PointVector.h"
#include "DGtal/shapes/SurfaceMesh.h"
#include "DGtal/io/writers/SurfaceMeshWriter.h"
#include "QuickHull.h"
using namespace DGtal::Z3i;
int main( int argc, char* argv[] )
{
int nb = argc > 1 ? atoi( argv[ 1 ] ) : 100; // nb points
int R = argc > 2 ? atoi( argv[ 2 ] ) : 10; // x-radius of ellipsoid
// (0) typedefs
// (1) create range of random points in ball
std::vector< Point > V;
const auto R2 = R
for ( int i = 0; i < nb; ) {
Point p( rand() % (2*R+1) - R, rand() % (2*R+1) - R, rand() % (2*R+1) - R );
if ( p.squaredNorm() < R2 ) { V.push_back( p ); i++; }
}
// (2) compute convex hull
hull.setInput( V );
std::cout << "#points=" << hull.nbPoints()
<< " #vertices=" << hull.nbVertices()
<< " #facets=" << hull.nbFacets() << std::endl;
// (3) build mesh
std::vector< RealPoint > positions;
hull.getVertexPositions( positions );
std::vector< std::vector< std::size_t > > facets;
SMesh mesh( positions.cbegin(), positions.cend(), facets.cbegin(), facets.cend() );
// (4) output result as OBJ file
std::ofstream out( "qhull.obj" );
::writeOBJ( out, mesh );
out.close();
return 0;
}
SurfaceMesh< RealPoint, RealVector > SMesh
Z3i this namespace gathers the standard of types for 3D imagery.
Aim: a geometric kernel to compute the convex hull of digital points with integer-only arithmetic.
Aim: Implements the quickhull algorithm by Barber et al. , a famous arbitrary dimensional convex hull...
Definition QuickHull.h:141
bool setInput(const std::vector< InputPoint > &input_points, bool remove_duplicates=true)
Definition QuickHull.h:383
std::vector< Facet > facets
the current set of facets.
Definition QuickHull.h:815
bool getFacetVertices(std::vector< IndexRange > &facet_vertices) const
Definition QuickHull.h:667
Kernel::CoordinatePoint Point
Definition QuickHull.h:143
Size nbVertices() const
Definition QuickHull.h:338
bool getVertexPositions(std::vector< OutputPoint > &vertex_positions)
Definition QuickHull.h:613
Size nbFacets() const
Definition QuickHull.h:329
bool computeConvexHull(Status target=Status::VerticesCompleted)
Definition QuickHull.h:461
Size nbPoints() const
Definition QuickHull.h:324
Aim: An helper class for writing mesh file formats (Waverfront OBJ at this point) and creating a Surf...
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
Definition SurfaceMesh.h:92
int main()
Definition testBits.cpp:56
Note
In opposition with the usual QuickHull implementation, this class uses a kernel that can be chosen in order to provide exact computations. This is the case for lattice points.
In opposition with CGAL 3D convex hull package, or with the arbitrary dimensional dD Triangulation package, this algorithm does not build a simplicial convex hull. Facets may not be trangles or simplices in higher dimensions.
This version is generally more than twice faster than CGAL convex_hull_3 for the usual CGAL kernels Cartesian and Exact_predicates_inexact_constructions_kernel.
However this implementation is not tailored for incremental dynamic convex hull computations.
Template Parameters
TKernelany type of QuickHull kernel, like ConvexHullIntegralKernel.
Examples
geometry/tools/exampleLatticeBallDelaunay2D.cpp, geometry/tools/exampleLatticeBallQuickHull3D.cpp, and geometry/tools/exampleRationalBallQuickHull3D.cpp.

Definition at line 140 of file QuickHull.h.

Member Typedef Documentation

◆ CombinatorialPlaneSimplex

template<typename TKernel >
typedef Kernel::CombinatorialPlaneSimplex DGtal::QuickHull< TKernel >::CombinatorialPlaneSimplex

Definition at line 152 of file QuickHull.h.

◆ HalfSpace

template<typename TKernel >
typedef Kernel::HalfSpace DGtal::QuickHull< TKernel >::HalfSpace

Definition at line 151 of file QuickHull.h.

◆ Index

template<typename TKernel >
typedef std::size_t DGtal::QuickHull< TKernel >::Index

Definition at line 147 of file QuickHull.h.

◆ IndexRange

template<typename TKernel >
typedef std::vector< Index > DGtal::QuickHull< TKernel >::IndexRange

Definition at line 150 of file QuickHull.h.

◆ InternalScalar

template<typename TKernel >
typedef Kernel::InternalScalar DGtal::QuickHull< TKernel >::InternalScalar

Definition at line 146 of file QuickHull.h.

◆ Kernel

template<typename TKernel >
typedef TKernel DGtal::QuickHull< TKernel >::Kernel

Definition at line 142 of file QuickHull.h.

◆ Point

template<typename TKernel >
typedef Kernel::CoordinatePoint DGtal::QuickHull< TKernel >::Point

Definition at line 143 of file QuickHull.h.

◆ Ridge

template<typename TKernel >
typedef std::pair< Index, Index > DGtal::QuickHull< TKernel >::Ridge

A ridge for point p is a pair of facets, such that p is visible from first facet, but not from second facet.

Definition at line 237 of file QuickHull.h.

◆ Scalar

template<typename TKernel >
typedef Kernel::CoordinateScalar DGtal::QuickHull< TKernel >::Scalar

Definition at line 145 of file QuickHull.h.

◆ Size

template<typename TKernel >
typedef std::size_t DGtal::QuickHull< TKernel >::Size

Definition at line 148 of file QuickHull.h.

◆ Vector

template<typename TKernel >
typedef Kernel::CoordinateVector DGtal::QuickHull< TKernel >::Vector

Definition at line 144 of file QuickHull.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename TKernel >
anonymous enum

Label for points that are not assigned to any facet.

Enumerator
UNASSIGNED 

Definition at line 156 of file QuickHull.h.

156{ UNASSIGNED = (Index) -1 };
std::size_t Index
Definition QuickHull.h:147

◆ Status

template<typename TKernel >
enum class DGtal::QuickHull::Status
strong

Represents the status of a QuickHull object.

Enumerator
Uninitialized 

QuickHull is empty and has just been instantiated.

InputInitialized 

A range of input points has been given to QuickHull.

SimplexCompleted 

An initial full-dimensional simplex has been found. QuickHull core algorithm can start.

FacetsCompleted 

All facets of the convex hull are identified.

VerticesCompleted 

All vertices of the convex hull are determined.

AllCompleted 

Same as VerticesCompleted.

NotFullDimensional 

Error: the initial simplex is not full dimensional.

InvalidRidge 

Error: some ridge is not consistent (probably integer overflow).

InvalidConvexHull 

Error: the convex hull does not contain all its points (probably integer overflow).

Definition at line 240 of file QuickHull.h.

240 {
241 Uninitialized = 0,
242 InputInitialized = 1,
243 SimplexCompleted = 2,
244 FacetsCompleted = 3,
246 AllCompleted = 5,
248 InvalidRidge = 11,
250 };
@ InvalidRidge
Error: some ridge is not consistent (probably integer overflow).
@ InvalidConvexHull
Error: the convex hull does not contain all its points (probably integer overflow).
@ SimplexCompleted
An initial full-dimensional simplex has been found. QuickHull core algorithm can start.
@ FacetsCompleted
All facets of the convex hull are identified.
@ NotFullDimensional
Error: the initial simplex is not full dimensional.
@ VerticesCompleted
All vertices of the convex hull are determined.
@ InputInitialized
A range of input points has been given to QuickHull.
@ AllCompleted
Same as VerticesCompleted.
@ Uninitialized
QuickHull is empty and has just been instantiated.

Constructor & Destructor Documentation

◆ QuickHull()

template<typename TKernel >
DGtal::QuickHull< TKernel >::QuickHull ( const Kernel K_ = Kernel(),
int  dbg = 0 
)
inline

Default constructor

Parameters
[in]K_a kernel for computing facet geometries.
[in]dbgthe trace level, from 0 (no) to 3 (very verbose).

Definition at line 260 of file QuickHull.h.

262 {}
Kernel kernel
Kernel that is duplicated for computing facet geometries.
Definition QuickHull.h:799
int debug_level
debug_level from 0:no to 2
Definition QuickHull.h:801

Member Function Documentation

◆ above()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::above ( const Facet F,
const Point p 
) const
inlineprotected

◆ aboveOrOn()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::aboveOrOn ( const Facet F,
const Point p 
) const
inlineprotected
Parameters
Fany valid facet
pany point
Returns
'true' iff p is above F or on F.

Definition at line 863 of file QuickHull.h.

864 { return kernel.aboveOrOn( F.H, p ); }

References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.

Referenced by DGtal::QuickHull< TKernel >::checkFacet(), and DGtal::QuickHull< TKernel >::processFacet().

◆ areFacetsNeighbor()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::areFacetsNeighbor ( const Index  if1,
const Index  if2 
) const
inlineprotected

Checks if two facets are neighbors by looking at the points on their boundary.

Parameters
if1a valid facet index
if2a valid facet index
Returns
'true' if the two facets have enough common points to be direct neighbors.

Definition at line 1268 of file QuickHull.h.

1269 {
1270 const Facet& f1 = facets[ if1 ];
1271 const Facet& f2 = facets[ if2 ];
1272 Index nb = 0;
1273 for ( Index i1 = 0, i2 = 0; i1 < f1.on_set.size() && i2 < f2.on_set.size(); )
1274 {
1275 if ( f1.on_set[ i1 ] == f2.on_set[ i2 ] ) { nb++; i1++; i2++; }
1276 else if ( f1.on_set[ i1 ] < f2.on_set[ i2 ] ) i1++;
1277 else i2++;
1278 }
1279 return nb >= ( dimension - 1 );
1280 }
SMesh::Index Index
static const Size dimension
Definition QuickHull.h:153

References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, and DGtal::QuickHull< TKernel >::Facet::on_set.

Referenced by DGtal::QuickHull< TKernel >::checkFacet(), and DGtal::QuickHull< TKernel >::processFacet().

◆ areFacetsParallel()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::areFacetsParallel ( const Index  if1,
const Index  if2 
) const
inlineprotected

Checks if two distinct facets are parallel (i.e. should be merged).

Parameters
if1a valid facet index
if2a valid facet index
Returns
'true' if the two facets are parallel.

Definition at line 1249 of file QuickHull.h.

1250 {
1251 ASSERT( if1 != if2 );
1252 const Facet& f1 = facets[ if1 ];
1253 const Facet& f2 = facets[ if2 ];
1254 if ( kernel.equal( f1.H, f2.H ) ) return true;
1255 // Need to check if one N is a multiple of the other.
1256 for ( auto&& v : f1.on_set )
1257 if ( ! on( f2, points[ v ] ) ) return false;
1258 return true;
1259 }
bool on(const Facet &F, const Point &p) const
Definition QuickHull.h:869
std::vector< Point > points
the set of points, indexed as in the array.
Definition QuickHull.h:803

References DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::Facet::H, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::on(), DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::QuickHull< TKernel >::points.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ BOOST_STATIC_ASSERT()

template<typename TKernel >
DGtal::QuickHull< TKernel >::BOOST_STATIC_ASSERT ( (Point::dimension==Vector::dimension)  )

◆ check()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::check ( )
inline

Global validity check of the convex hull after processing.

Note
Be careful, this function is much slower than computing the convex hull. It can take a long time since its complexity is \( O(n f ) \), where n is the number of input points and f the number of facets.

Definition at line 716 of file QuickHull.h.

717 {
718 bool ok = true;
720 trace.warning() << "[Quickhull::check] invalid status="
721 << (int)status() << std::endl;
722 ok = false;
723 }
724 if ( processed_points.size() != points.size() ) {
725 trace.warning() << "[Quickhull::check] not all points processed: "
726 << processed_points.size() << " / " << points.size()
727 << std::endl;
728 ok = false;
729 }
730 if ( ! checkHull() ) {
731 trace.warning() << "[Quickhull::check] Check hull is invalid. "
732 << std::endl;
733 ok = false;
734 }
735 if ( ! checkFacets() ) {
736 trace.warning() << "[Quickhull::check] Check facets is invalid. "
737 << std::endl;
738 ok = false;
739 }
740 return ok;
741 }
std::ostream & warning()
Trace trace
IndexRange processed_points
Points already processed (and within the convex hull).
Definition QuickHull.h:811
Status status() const
Definition QuickHull.h:267

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::status(), DGtal::trace, and DGtal::Trace::warning().

◆ checkFacet()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::checkFacet ( Index  f) const
inlineprotected
Returns
true if the facet is valid

Definition at line 1137 of file QuickHull.h.

1138 {
1139 const Facet& F = facets[ f ];
1140 bool ok = F.on_set.size() >= dimension;
1141 for ( auto v : F.on_set )
1142 if ( ! on( F, points[ v ] ) ) {
1143 trace.error() << "[QuickHull::checkFacet( " << f
1144 << ") Invalid 'on' vertex " << v << std::endl;
1145 ok = false;
1146 }
1147 if ( F.neighbors.size() < dimension ) {
1148 trace.error() << "[QuickHull::checkFacet( " << f
1149 << ") Not enough neighbors " << F.neighbors.size() << std::endl;
1150 ok = false;
1151 }
1152 for ( auto nf : F.neighbors )
1153 if ( ! areFacetsNeighbor( f, nf ) ) {
1154 trace.error() << "[QuickHull::checkFacet( " << f
1155 << ") Invalid neighbor " << nf << std::endl;
1156 ok = false;
1157 }
1158 if ( aboveOrOn( F, points[ F.below ] ) ) {
1159 trace.error() << "[QuickHull::checkFacet( " << f
1160 << ") Bad below point " << F.below << std::endl;
1161 ok = false;
1162 }
1163 for ( auto ov : F.outside_set )
1164 if ( ! above( F, points[ ov ] ) ) {
1165 trace.error() << "[QuickHull::checkFacet( " << f
1166 << ") Bad outside point " << ov << std::endl;
1167 ok = false;
1168 }
1169 for ( auto && v : F.on_set ) {
1170 Size n = 0;
1171 for ( auto&& N : facets[ f ].neighbors )
1172 if ( on( facets[ N ], points[ v ] ) ) n += 1;
1173 if ( n < dimension-1 ) {
1174 trace.error() << "[QuickHull::checkFacet( " << f << ") 'on' point " << v
1175 << " is a vertex of " << n << " facets "
1176 << "(should be >= " << dimension-1 << ")" << std::endl;
1177 ok = false;
1178 }
1179 }
1180 return ok;
1181 }
std::ostream & error()
bool areFacetsNeighbor(const Index if1, const Index if2) const
Definition QuickHull.h:1268
bool aboveOrOn(const Facet &F, const Point &p) const
Definition QuickHull.h:863
bool above(const Facet &F, const Point &p) const
Definition QuickHull.h:857
HalfEdgeDataStructure::Size Size

References DGtal::QuickHull< TKernel >::above(), DGtal::QuickHull< TKernel >::aboveOrOn(), DGtal::QuickHull< TKernel >::areFacetsNeighbor(), DGtal::QuickHull< TKernel >::Facet::below, DGtal::QuickHull< TKernel >::dimension, DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::Facet::neighbors, DGtal::QuickHull< TKernel >::on(), DGtal::QuickHull< TKernel >::Facet::on_set, DGtal::QuickHull< TKernel >::Facet::outside_set, DGtal::QuickHull< TKernel >::points, and DGtal::trace.

Referenced by DGtal::QuickHull< TKernel >::checkFacets(), and DGtal::QuickHull< TKernel >::processFacet().

◆ checkFacets()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::checkFacets ( )
inline
Returns
'true' if all facets are consistent with their neighors

Definition at line 744 of file QuickHull.h.

745 {
746 Size nb = 0;
747 Size nbok = 0;
748 for ( Index f = 0; f < facets.size(); ++f )
749 if ( deleted_facets.count( f ) == 0 ) {
750 bool ok = checkFacet( f );
751 nbok += ok ? 1 : 0;
752 nb += 1;
753 }
754 return nb == nbok;
755 }
bool checkFacet(Index f) const
Definition QuickHull.h:1137
std::set< Index > deleted_facets
set of deleted facets
Definition QuickHull.h:817

References DGtal::QuickHull< TKernel >::checkFacet(), DGtal::QuickHull< TKernel >::deleted_facets, and DGtal::QuickHull< TKernel >::facets.

Referenced by DGtal::QuickHull< TKernel >::check(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().

◆ checkHull()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::checkHull ( )
inline
Returns
'true' iff the current set of facets encloses all the points
Note
Be careful, this function is much slower than computing the convex hull. It can take a long time since its complexity is \( O(n f ) \), where n is the number of input points and f the number of facets.

Definition at line 763 of file QuickHull.h.

764 {
765 Index nbok = 0;
766 // for ( Index v = 0; v < points.size(); ++v ) {
767 for ( auto v : processed_points ) {
768 bool ok = true;
769 for ( Index f = 0; f < facets.size(); ++f )
770 if ( deleted_facets.count( f ) == 0 ) {
771 if ( above( facets[ f ], points[ v ] ) ) {
772 ok = false;
773 trace.error() << "- bad vertex " << v << " " << points[ v ]
774 << " dist="
775 << ( kernel.height( facets[ f ].H, points[ v ] ) )
776 << std::endl;
777 break;
778 }
779 }
780 nbok += ok ? 1 : 0;
781 }
782 if ( debug_level >= 2 ) {
783 trace.info() << nbok << "/"
784 << processed_points.size() << " vertices inside convex hull."
785 << std::endl;
786 }
788 return nbok == processed_points.size();
789 }
std::ostream & info()

References DGtal::QuickHull< TKernel >::above(), DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::H, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::InvalidConvexHull, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::processed_points, and DGtal::trace.

Referenced by DGtal::QuickHull< TKernel >::check(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().

◆ cleanFacets()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::cleanFacets ( )
inlineprotected

Cleans and renumber the facets so that no one belongs to deleted_facets.

Definition at line 874 of file QuickHull.h.

875 {
876 if ( deleted_facets.empty() ) return;
877 IndexRange renumbering( facets.size() );
878 Index i = 0;
879 Index j = 0;
880 for ( auto& l : renumbering ) {
881 if ( ! deleted_facets.count( j ) ) l = i++;
882 else l = UNASSIGNED;
883 j++;
884 }
885 const Index nf = facets.size() - deleted_facets.size();
886 deleted_facets.clear();
887 for ( Index f = 0; f < facets.size(); f++ )
888 if ( ( renumbering[ f ] != UNASSIGNED ) && ( f != renumbering[ f ] ) )
889 facets[ renumbering[ f ] ] = facets[ f ];
890 facets.resize( nf );
891 for ( auto& F : facets ) {
892 for ( auto& N : F.neighbors ) {
893 if ( renumbering[ N ] == UNASSIGNED )
894 trace.error() << "Invalid deleted neighboring facet." << std::endl;
895 else N = renumbering[ N ];
896 }
897 }
898 }
std::vector< Index > IndexRange
Definition QuickHull.h:150

References DGtal::QuickHull< TKernel >::deleted_facets, DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::trace, and DGtal::QuickHull< TKernel >::UNASSIGNED.

Referenced by DGtal::QuickHull< TKernel >::computeFacets().

◆ clear()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::clear ( )
inline

Clears the object.

Definition at line 271 of file QuickHull.h.

272 {
274 points.clear();
275 processed_points.clear();
276 input2comp.clear();
277 comp2input.clear();
278 assignment.clear();
279 facets.clear();
280 deleted_facets.clear();
281 p2v.clear();
282 v2p.clear();
283 timings.clear();
284 }
IndexRange v2p
vertex index -> point index
Definition QuickHull.h:821
IndexRange input2comp
Definition QuickHull.h:806
std::vector< Index > assignment
assignment of points to facets
Definition QuickHull.h:813
IndexRange p2v
point index -> vertex index (or UNASSIGNED)
Definition QuickHull.h:819
IndexRange comp2input
Definition QuickHull.h:809
std::vector< double > timings
Timings of the different phases: 0: init, 1: facets, 2: vertices.
Definition QuickHull.h:827

References DGtal::QuickHull< TKernel >::assignment, DGtal::QuickHull< TKernel >::comp2input, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::input2comp, DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::timings, DGtal::QuickHull< TKernel >::Uninitialized, and DGtal::QuickHull< TKernel >::v2p.

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::clear(), DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::compute(), and DGtal::QuickHull< TKernel >::setInput().

◆ computeConvexHull()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::computeConvexHull ( Status  target = Status::VerticesCompleted)
inline

Computes the convex hull of the given range of points until a specified target:

  • Status::SimplexCompleted: an initial full dimensional simplex is determined.
  • Status::FacetsCompleted: all the facets of the convex hull are extracted (core of quickhull). You can stop here if you only need an H-representation of the output polytope.
  • Status::VerticesCompleted: all the vertices of the hull are determined and are consistently ordered on each facet. You need to stop here if you need a V-representation of the output polytope.
Precondition
status() must be at least Status::InputInitialized
Parameters
[in]targetthe computation target in Status::SimplexCompleted, Status::FacetsCompleted, Status::VerticesCompleted.
Returns
'true' if the computation target has been successfully achieved, 'false' if the achieved status is not the one specified.
Examples
geometry/tools/exampleLatticeBallQuickHull3D.cpp.

Definition at line 461 of file QuickHull.h.

462 {
463 if ( target < Status::InputInitialized || target > Status::AllCompleted )
464 return false;
465 Clock tic;
467 { // Initialization
468 tic.startClock();
469 bool ok1 = computeInitialSimplex();
470 timings.push_back( tic.stopClock() );
471 if ( ! ok1 ) return false;
472 if ( status() == target ) return true;
473 }
475 { // Computes facets
476 tic.startClock();
477 bool ok2 = computeFacets();
478 timings.push_back( tic.stopClock() );
479 if ( ! ok2 ) return false;
480 if ( status() == target ) return true;
481 }
483 { // Computes vertices
484 tic.startClock();
485 bool ok3 = computeVertices();
486 timings.push_back( tic.stopClock() );
487 if ( ! ok3 ) return false;
488 if ( status() == target ) return true;
489 }
490 if ( target == Status::AllCompleted
492 { // for now, Status::VerticesCompleted and
493 // Status::AllCompleted are the same.
495 return true;
496 }
497 return false;
498 }
void tic()
Starts timer.
bool computeInitialSimplex()
Definition QuickHull.h:505
bool computeFacets()
Definition QuickHull.h:524
bool computeVertices()
Definition QuickHull.h:554

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeInitialSimplex(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::QuickHull< TKernel >::InputInitialized, DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::SimplexCompleted, DGtal::Clock::startClock(), DGtal::QuickHull< TKernel >::status(), tic(), DGtal::QuickHull< TKernel >::timings, and DGtal::QuickHull< TKernel >::VerticesCompleted.

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::compute(), and main().

◆ computeFacets()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::computeFacets ( )
inline

Computes the facets of the convex hull using Quickhull algorithm. If everything went well, the status is Status::FacetsCompleted afterwards.

Precondition
the status shoud be Status::SimplexCompleted (computeInitialSimplex should have been called).
Returns
'true' except if the status is not Initialized when called.

Definition at line 524 of file QuickHull.h.

525 {
526 if ( status() != Status::SimplexCompleted ) return false;
527 std::queue< Index > Q;
528 for ( Index fi = 0; fi < facets.size(); ++fi )
529 Q.push( fi );
530 Index n = 0;
531 while ( processFacet( Q ) ) {
532 if ( debug_level >= 1 )
533 trace.info() << "---- Iteration " << n++ << " #Q=" << Q.size() << std::endl;
534 }
535 cleanFacets();
536 if ( debug_level >= 2 ) {
537 trace.info() << ".... #facets=" << facets.size()
538 << " #deleted=" << deleted_facets.size() << std::endl;
539 }
541 return true;
542 }
bool processFacet(std::queue< Index > &Q)
Definition QuickHull.h:945

References DGtal::QuickHull< TKernel >::cleanFacets(), DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::processFacet(), DGtal::QuickHull< TKernel >::SimplexCompleted, DGtal::QuickHull< TKernel >::status(), and DGtal::trace.

Referenced by DGtal::QuickHull< TKernel >::computeConvexHull().

◆ computeInitialSimplex()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::computeInitialSimplex ( )
inline

Computes the initial full dimensional simplex from the input data.

Returns
'true' iff the input data contains d+1 points in general position, the object has then the status SimplexCompleted, otherwise returns 'false' and the status is NotFullDimensional.

Definition at line 505 of file QuickHull.h.

506 {
507 const auto full_simplex = pickInitialSimplex();
508 if ( full_simplex.empty() ) {
510 return false;
511 }
512 return computeSimplexConfiguration( full_simplex );
513 }
bool computeSimplexConfiguration(const IndexRange &full_simplex)
Definition QuickHull.h:1453
IndexRange pickInitialSimplex() const
Definition QuickHull.h:1382

References DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::NotFullDimensional, and DGtal::QuickHull< TKernel >::pickInitialSimplex().

Referenced by DGtal::QuickHull< TKernel >::computeConvexHull().

◆ computeSimplexConfiguration()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::computeSimplexConfiguration ( const IndexRange full_simplex)
inlineprotected

Computes the initial configuration induced by the given full dimensional simplex. Follows Status::InputInitialized and and terminate with Status::SimplexCompleted if everything went well.

Returns
'true' iff the input data contains d+1 points in general position, the object has then the status SimplexCompleted, otherwise returns 'false' and the status is NotFullDimensional.

Definition at line 1453 of file QuickHull.h.

1454 {
1455 assignment = std::vector< Index >( points.size(), UNASSIGNED );
1456 facets.resize( dimension + 1 );
1457 deleted_facets.clear();
1458 for ( Index j = 0; j < full_simplex.size(); ++j )
1459 {
1460 IndexRange lsimplex( dimension );
1461 IndexRange isimplex( dimension );
1462 Index s = 0;
1463 for ( Index i = 0; i <= dimension; i++ )
1464 if ( i != j ) {
1465 lsimplex[ s ] = i;
1466 isimplex[ s ] = full_simplex[ i ];
1467 s++;
1468 }
1469 facets[ j ] = makeFacet( isimplex, full_simplex[ j ] );
1470 facets[ j ].neighbors = lsimplex;
1471 for ( auto&& v : isimplex ) facets[ j ].on_set.push_back( v );
1472 std::sort( facets[ j ].on_set.begin(), facets[ j ].on_set.end() );
1473 }
1474 // List of unassigned vertices
1475 for ( Index fi = 0; fi < facets.size(); ++fi ) {
1476 Facet& f = facets[ fi ];
1477 for ( Index v = 0; v < points.size(); v++ )
1478 {
1479 if ( assignment[ v ] == UNASSIGNED && above( f, points[ v ] ) ) {
1480 f.outside_set.push_back( v );
1481 assignment[ v ] = fi;
1482 }
1483 }
1484 }
1485 for ( Index v = 0; v < points.size(); v++ )
1486 if ( assignment[ v ] == UNASSIGNED )
1487 processed_points.push_back( v );
1488
1489 // Display some information
1490 if ( debug_level >= 2 ) {
1491 for ( auto&& f : facets ) f.display( trace.info() );
1492 }
1494 if ( debug_level >= 1 ) {
1495 trace.info() << "[CHECK INVARIANT] " << processed_points.size()
1496 << " / " << points.size() << " points processed." << std::endl;
1497 bool okh = checkHull();
1498 if ( ! okh )
1499 trace.error() << "[computeInitialSimplex] Invalid convex hull" << std::endl;
1500 bool okf = checkFacets();
1501 if ( ! okf )
1502 trace.error() << "[computeInitialSimplex] Invalid facets" << std::endl;
1503 if ( ! ( okh && okf ) ) myStatus = Status::InvalidConvexHull;
1504 }
1505 return status() == Status::SimplexCompleted;
1506 }
Facet makeFacet(const IndexRange &base, Index below) const
Definition QuickHull.h:1293
void display(ostream &out, const AContainer &C)

References DGtal::QuickHull< TKernel >::above(), DGtal::QuickHull< TKernel >::assignment, DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::QuickHull< TKernel >::dimension, DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::InvalidConvexHull, DGtal::QuickHull< TKernel >::makeFacet(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::Facet::outside_set, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::SimplexCompleted, DGtal::QuickHull< TKernel >::status(), DGtal::trace, and DGtal::QuickHull< TKernel >::UNASSIGNED.

Referenced by DGtal::QuickHull< TKernel >::computeInitialSimplex(), and DGtal::QuickHull< TKernel >::setInitialSimplex().

◆ computeVertices()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::computeVertices ( )
inline

Computes the vertices of the convex hull once the facets have been computed. It computes for each facet its vertices and reorder them so that, taken in order, their orientation matches the orientation of the facet. If everything went well, the status is Status::VerticesCompleted afterwards.

Precondition
the status shoud be Status::FacetsCompleted (computeFacets should have been called).
Returns
'true' except if the status is not FacetsCompleted when called.

Definition at line 554 of file QuickHull.h.

555 {
556 static const int MAX_NB_VPF = 10 * dimension;
557 if ( status() != Status::FacetsCompleted ) return false;
558
559 // Renumber infinite facets in case of Delaunay triangulation computation.
561
562 // Builds the maps v2p: vertex -> point, and p2v : point -> vertex.
563 facet_counter = IndexRange( MAX_NB_VPF, 0 );
564 v2p.clear();
565 p2v.resize( points.size() );
566 std::vector< IndexRange > p2f( points.size() );
567 for ( Index f = 0; f < facets.size(); ++f ) {
568 for ( auto&& p : facets[ f ].on_set ) p2f[ p ].push_back( f );
569 }
570
571 // vertices belong to at least d facets
572 Index v = 0;
573 for ( Index p = 0; p < points.size(); ++p ) {
574 const auto nbf = p2f[ p ].size();
575 facet_counter[ std::min( (int) nbf, MAX_NB_VPF-1 ) ] += 1;
576 if ( nbf >= dimension ) {
577 v2p.push_back( p );
578 p2v[ p ] = v++;
579 }
580 else p2v[ p ] = UNASSIGNED;
581 }
582
583 // Display debug informations
584 if ( debug_level >= 1 ) {
585 trace.info() << "#vertices=" << v2p.size() << " #facets=" << facets.size()
586 << std::endl;
587 trace.info() << "#inc_facets/point= ";
588 for ( auto n : facet_counter ) trace.info() << n << " ";
589 trace.info() << std::endl;
590 }
592 return true;
593 }
std::vector< Size > facet_counter
Counts the number of facets with a given number of vertices.
Definition QuickHull.h:829
void renumberInfiniteFacets()
Definition QuickHull.h:902

References DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facet_counter, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::renumberInfiniteFacets(), DGtal::QuickHull< TKernel >::status(), DGtal::trace, DGtal::QuickHull< TKernel >::UNASSIGNED, DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.

Referenced by DGtal::QuickHull< TKernel >::computeConvexHull().

◆ deleteFacet()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::deleteFacet ( Index  f)
inlineprotected

Deletes the given facet f.

Parameters
fa valid facet index

Definition at line 1194 of file QuickHull.h.

1195 {
1196 for ( auto n : facets[ f ].neighbors )
1197 facets[ n ].subNeighbor( f );
1198 deleted_facets.insert( f );
1199 facets[ f ].clear();
1200 }

References DGtal::QuickHull< TKernel >::deleted_facets, and DGtal::QuickHull< TKernel >::facets.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ filterVerticesOnFacet()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::filterVerticesOnFacet ( const Index  f)
inlineprotected

Filters each vertex on the facet f to keep only the ones that are on or below the neighboring facets

Note
intended for debugging purposes.
Parameters
[in]fany valid facet index.

Definition at line 1362 of file QuickHull.h.

1363 {
1364 auto & on_set = facets[ f ].on_set;
1365 for ( Index i = 0; i < on_set.size(); )
1366 {
1367 Index v = on_set[ i ];
1368 Size n = 0;
1369 for ( auto&& N : facets[ f ].neighbors )
1370 if ( on( facets[ N ], points[ v ] ) ) n += 1;
1371 if ( n >= dimension-1 ) i++;
1372 else {
1373 on_set[ i ] = on_set.back();
1374 on_set.pop_back();
1375 }
1376 }
1377 std::sort( on_set.begin(), on_set.end() );
1378 }

References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::on(), and DGtal::QuickHull< TKernel >::points.

◆ getFacetHalfSpaces()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::getFacetHalfSpaces ( std::vector< HalfSpace > &  facet_halfspaces)
inline

This methods return the halfspaces corresponding to each facet. It is useful to build a polytope.

Parameters
[out]facet_halfspacesthe range giving for each facet its halfspace as a pair (normal, intercept).
Returns
'true' if the convex hull was computed before and status() was Status::FacetsCompleted, Status::VerticesCompleted or Status::AllCompleted, 'false' otherwise.

Definition at line 691 of file QuickHull.h.

692 {
693 facet_halfspaces.clear();
694 if ( ! ( status() >= Status::FacetsCompleted
695 && status() <= Status::AllCompleted ) ) return false;
696 facet_halfspaces.reserve( nbFacets() );
697 for ( Index f = 0; f < nbFacets(); ++f ) {
698 facet_halfspaces.push_back( facets[ f ].H );
699 }
700 return true;
701 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::H, DGtal::QuickHull< TKernel >::nbFacets(), and DGtal::QuickHull< TKernel >::status().

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::makePolytope().

◆ getFacetVertices()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::getFacetVertices ( std::vector< IndexRange > &  facet_vertices) const
inline
Parameters
[out]facet_verticesthe range giving for each facet the indices of its vertices.
Returns
'true' if the convex hull was computed before and status() was Status::VerticesCompleted or Status::AllCompleted, 'false' otherwise.
Examples
geometry/tools/exampleLatticeBallQuickHull3D.cpp.

Definition at line 667 of file QuickHull.h.

668 {
669 facet_vertices.clear();
671 && status() <= Status::AllCompleted ) ) return false;
672 facet_vertices.reserve( nbFacets() );
673 for ( Index f = 0; f < nbFacets(); ++f ) {
674 IndexRange ofacet = orientedFacetPoints( f );
675 for ( auto& v : ofacet ) v = p2v[ v ];
676 facet_vertices.push_back( ofacet );
677 }
678 return true;
679 }
IndexRange orientedFacetPoints(Index f) const
Definition QuickHull.h:1324

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::nbFacets(), DGtal::QuickHull< TKernel >::orientedFacetPoints(), DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::compute(), and main().

◆ getPoint2Vertex()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::getPoint2Vertex ( IndexRange point_to_vertex)
inline

Gets for each point its index in the output range of vertices (or UNASSIGNED if the point is not part of the convex hull)

Parameters
[out]point_to_vertexthe range giving for each point its index in the output range of vertices (or UNASSIGNED if the point is not part of the convex hull)
Returns
'true' if the convex hull was computed before and status() was Status::VerticesCompleted or Status::AllCompleted, 'false' otherwise.

Definition at line 652 of file QuickHull.h.

653 {
654 point_to_vertex.clear();
656 && status() <= Status::AllCompleted ) ) return false;
657 point_to_vertex = p2v;
658 return true;
659 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.

◆ getVertex2Point()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::getVertex2Point ( IndexRange vertex_to_point)
inline

Gets for each vertex its index in the input range of points.

Parameters
[out]vertex_to_pointthe range giving for each vertex its index in the input range of points.
Returns
'true' if the convex hull was computed before and status() was Status::VerticesCompleted or Status::AllCompleted, 'false' otherwise.

Definition at line 633 of file QuickHull.h.

634 {
635 vertex_to_point.clear();
637 && status() <= Status::AllCompleted ) ) return false;
638 vertex_to_point = v2p;
639 return true;
640 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::compute().

◆ getVertexPositions()

template<typename TKernel >
template<typename OutputPoint >
bool DGtal::QuickHull< TKernel >::getVertexPositions ( std::vector< OutputPoint > &  vertex_positions)
inline

Gets the positions of the convex hull vertices.

Template Parameters
OutputPointa model of point such that the Point datatype is convertible to it.
Parameters
[out]vertex_positionsthe range of vertex positions.
Returns
'true' if the convex hull was computed before and status() was Status::VerticesCompleted or Status::AllCompleted, 'false' otherwise.
Examples
geometry/tools/exampleLatticeBallQuickHull3D.cpp.

Definition at line 613 of file QuickHull.h.

614 {
615 vertex_positions.clear();
617 && status() <= Status::AllCompleted ) ) return false;
618 vertex_positions.resize( v2p.size() );
619 for ( Index i = 0; i < v2p.size(); i++ ) {
620 kernel.convertPointTo( points[ v2p[ i ] ], vertex_positions[ i ] );
621 }
622 return true;
623 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.

Referenced by main().

◆ height()

template<typename TKernel >
InternalScalar DGtal::QuickHull< TKernel >::height ( const Facet F,
const Point p 
) const
inlineprotected
Parameters
Fany valid facet
pany point
Returns
the height of p wrt F (0: on, >0: above ).

Definition at line 851 of file QuickHull.h.

852 { return kernel.height( F.H, p ); }

References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ isValid()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::isValid ( ) const
inline

Checks the validity/consistency of the object.

Returns
'true' if the object is valid, 'false' otherwise.

Definition at line 1539 of file QuickHull.h.

1540 {
1541 return status() >= Status::Uninitialized
1543 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::Uninitialized.

◆ makeFacet()

template<typename TKernel >
Facet DGtal::QuickHull< TKernel >::makeFacet ( const IndexRange base,
Index  below 
) const
inlineprotected

Builds a facet from a base convex set of at least d different points and a point below.

Parameters
[in]basea range containing at least d distinct points in general position.
[in]belowa point below the hyperplane containing simplex.
Returns
a facet object of normal and c parameter such as the points of simplex lies on the facet hyperplane while point below lies under.

Definition at line 1293 of file QuickHull.h.

1294 {
1296 for ( Size i = 0; i < dimension; i++ ) simplex[ i ] = base[ i ];
1297 auto plane = kernel.compute( points, simplex, below );
1298 return Facet( plane, below );
1299 }
Kernel::CombinatorialPlaneSimplex CombinatorialPlaneSimplex
Definition QuickHull.h:152

References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::kernel, and DGtal::QuickHull< TKernel >::points.

Referenced by DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().

◆ makeNeighbors()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::makeNeighbors ( const Index  if1,
const Index  if2 
)
inlineprotected

Makes two distinct facets if1 and if2 as neighbors

Parameters
if1a valid facet index
if2a valid facet index

Definition at line 1205 of file QuickHull.h.

1206 {
1207 facets[ if1 ].addNeighbor( if2 );
1208 facets[ if2 ].addNeighbor( if1 );
1209 }

References DGtal::QuickHull< TKernel >::facets.

Referenced by DGtal::QuickHull< TKernel >::mergeParallelFacets(), and DGtal::QuickHull< TKernel >::processFacet().

◆ memory()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::memory ( ) const
inline
Returns
an estimation of the current memory occupation of this object.

Definition at line 288 of file QuickHull.h.

289 {
290 // int debug_level;
291 Size M = sizeof( kernel ) + sizeof( int );
292 // std::vector< Point > points;
293 M += sizeof( std::vector< Point > )
294 + points.capacity() * sizeof( Point );
295 M += sizeof( std::vector< Index > )
296 + processed_points.capacity() * sizeof( Index );
297 M += sizeof( std::vector< Index > )
298 + input2comp.capacity() * sizeof( Index );
299 M += sizeof( std::vector< Index > )
300 + comp2input.capacity() * sizeof( Index );
301 // std::vector< Index > assignment;
302 M += sizeof( std::vector< Index > )
303 + assignment.capacity() * sizeof( Index );
304 // std::vector< Facet > facets;
305 M += sizeof( std::vector< Facet > )
306 + facets.capacity() * sizeof( Facet );
307 for ( const auto& f : facets ) M += f.variableMemory();
308 // std::set< Index > deleted_facets;
309 M += sizeof( std::set< Index > )
310 + deleted_facets.size() * ( sizeof( Index ) + 2*sizeof(Index*) );
311 // IndexRange p2v;
312 M += sizeof( std::vector< Index > )
313 + p2v.capacity() * sizeof( Index );
314 // IndexRange v2p;
315 M += sizeof( std::vector< Index > )
316 + v2p.capacity() * sizeof( Index );
317 // std::vector< double > timings;
318 M += sizeof( std::vector< double > )
319 + timings.capacity() * sizeof( double );
320 return M;
321 }

References DGtal::QuickHull< TKernel >::assignment, DGtal::QuickHull< TKernel >::comp2input, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::input2comp, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::timings, and DGtal::QuickHull< TKernel >::v2p.

◆ mergeParallelFacets()

template<typename TKernel >
Index DGtal::QuickHull< TKernel >::mergeParallelFacets ( const Index  if1,
const Index  if2 
)
inlineprotected

Merge the two facets and return the index of the second one, which is the one deleted.

Parameters
if1a valid facet index
if2a valid facet index
Precondition
the two facets should be distinct and parallel

Definition at line 1226 of file QuickHull.h.

1227 {
1228 Facet& f1 = facets[ if1 ];
1229 Facet& f2 = facets[ if2 ];
1230 std::copy( f2.outside_set.cbegin(), f2.outside_set.cend(),
1231 std::back_inserter( f1.outside_set ) );
1232 IndexRange merge_idx;
1233 std::set_union( f1.on_set.cbegin(), f1.on_set.cend(),
1234 f2.on_set.cbegin(), f2.on_set.cend(),
1235 std::back_inserter( merge_idx ) );
1236 f1.on_set.swap( merge_idx );
1237 for ( auto && nf2 : f2.neighbors ) {
1238 if ( nf2 == if1 ) continue;
1239 facets[ nf2 ].subNeighbor( if2 );
1240 makeNeighbors( if1, nf2 );
1241 }
1242 return if2;
1243 }
void makeNeighbors(const Index if1, const Index if2)
Definition QuickHull.h:1205

References DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::makeNeighbors(), DGtal::QuickHull< TKernel >::Facet::neighbors, DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::QuickHull< TKernel >::Facet::outside_set.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ nbFacets()

◆ nbFiniteFacets()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nbFiniteFacets ( ) const
inline
Returns
the number of finite facets of the convex hull.
Precondition
status() >= Status::VerticesCompleted

Definition at line 347 of file QuickHull.h.

348 {
351 return nb_finite_facets;
352 }
Size nb_finite_facets
Number of finite facets.
Definition QuickHull.h:823

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::nb_finite_facets, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.

◆ nbInfiniteFacets()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nbInfiniteFacets ( ) const
inline
Returns
the number of infinite facets of the convex hull.
Precondition
status() >= Status::VerticesCompleted

Definition at line 356 of file QuickHull.h.

357 {
360 return nb_infinite_facets;
361 }
Size nb_infinite_facets
Number of infinite facets (!= 0 only for specific kernels)
Definition QuickHull.h:825

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::nb_infinite_facets, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.

◆ nbPoints()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nbPoints ( ) const
inline
Returns
the number of points used as input.
Examples
geometry/tools/exampleLatticeBallQuickHull3D.cpp.

Definition at line 324 of file QuickHull.h.

325 { return points.size(); }

References DGtal::QuickHull< TKernel >::points.

Referenced by main(), and DGtal::QuickHull< TKernel >::selfDisplay().

◆ nbVertices()

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nbVertices ( ) const
inline

◆ newFacet()

template<typename TKernel >
Index DGtal::QuickHull< TKernel >::newFacet ( )
inlineprotected
Returns
an unused facet index

Definition at line 1184 of file QuickHull.h.

1185 {
1186 // SLightly faster to postpone deletion of intermediate facets.
1187 const Index f = facets.size();
1188 facets.push_back( Facet() );
1189 return f;
1190 }

References DGtal::QuickHull< TKernel >::facets.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ on()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::on ( const Facet F,
const Point p 
) const
inlineprotected
Parameters
Fany valid facet
pany point
Returns
'true' iff p lies on F.

Definition at line 869 of file QuickHull.h.

870 { return kernel.on( F.H, p ); }

References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.

Referenced by DGtal::QuickHull< TKernel >::areFacetsParallel(), DGtal::QuickHull< TKernel >::checkFacet(), and DGtal::QuickHull< TKernel >::filterVerticesOnFacet().

◆ orientedFacetPoints()

template<typename TKernel >
IndexRange DGtal::QuickHull< TKernel >::orientedFacetPoints ( Index  f) const
inlineprotected

Given a facet index f, return its points oriented consistently with respect to the normal.

Parameters
fany valid facet index
Returns
the range of point indices that support this facet, in a consistent ordering.
Note
the order of points is consistent if, picking any d-simplex in order in this range, their associated half-spaces have all the same orientation.

Definition at line 1324 of file QuickHull.h.

1325 {
1326 const Facet& F = facets[ f ];
1327 IndexRange result = F.on_set;
1328 // Sort a facet such that its points, taken in order, have
1329 // always the same orientation of the facet. More precisely,
1330 // facets span a `dimension-1` vector space. There are thus
1331 // dimension-2 fixed points, and the last ones (at least two)
1332 // may be reordered.
1334 for ( Dimension k = 0; k < dimension-2; k++ )
1335 splx[ k ] = result[ k ];
1336 // std::cout << "Orienting face " << f << " (";
1337 // for ( auto i : result ) std::cout << " " << i;
1338 std::sort( result.begin()+dimension-2, result.end(),
1339 [&] ( Index i, Index j )
1340 {
1341 splx[ dimension-2 ] = i;
1342 splx[ dimension-1 ] = j;
1343 const auto H = kernel.compute( points, splx );
1344 const auto orient = kernel.dot( F.H, H );
1345 return orient > 0;
1346 } );
1347 for ( Dimension k = 0; k < dimension; k++ )
1348 splx[ k ] = result[ k ];
1349 // const auto H = kernel.compute( points, splx );
1350 // const auto orient = kernel.dot( F.H, H );
1351 // std::cout << " ) => (";
1352 // for ( auto i : result ) std::cout << " " << i;
1353 // std::cout << " ) N=" << H.internalNormal() << " dot=" << orient << "\n";
1354 return result;
1355 }
DGtal::uint32_t Dimension
Definition Common.h:119

References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, and DGtal::QuickHull< TKernel >::Facet::on_set.

Referenced by DGtal::QuickHull< TKernel >::getFacetVertices().

◆ pickInitialSimplex()

template<typename TKernel >
IndexRange DGtal::QuickHull< TKernel >::pickInitialSimplex ( ) const
inlineprotected
Returns
a full dimensional simplex as a vector of d + 1 distinct indices of input points, or an empty vector if none was found.

Definition at line 1382 of file QuickHull.h.

1383 {
1384 const Size nb = points.size();
1385 if ( nb < dimension + 1 ) return IndexRange();
1386 IndexRange best = pickIntegers( dimension + 1, nb );
1388 for ( Index j = 0; j < dimension; ++j ) splx[ j ] = best[ j ];
1389 const auto first_H = kernel.compute( points, splx, best.back() );
1390 auto best_volume = kernel.volume ( first_H, points[ best.back() ] );
1391 // Randomized approach to find full dimensional simplex.
1392 //
1393 // Let a be the proportion of full dimensional simplices among
1394 // all simplices of the input points. Let p be the desired
1395 // probability to find a full dimensional simplex after t tries.
1396 // Then t >= log(1-p)/log(1-a)
1397 // For a=0.25, p=99% we found 16 tries.
1398 // For a=0.20, p=99% we found 20 tries.
1399 const Size nbtries = std::min( (Size) 10, 1 + nb / 10 );
1400 // const Size max_nbtries = std::max( (Size) 10, 2 * nb );
1401 const Size max_nbtries = 20;
1402 for ( Size i = 0; i < max_nbtries; i++ )
1403 {
1404 IndexRange tmp = pickIntegers( dimension + 1, nb );
1405 for ( Index j = 0; j < dimension; ++j ) splx[ j ] = tmp[ j ];
1406 const auto tmp_H = kernel.compute( points, splx, tmp.back() );
1407 const auto tmp_volume = kernel.volume ( tmp_H, points[ tmp.back() ] );
1408 if ( best_volume < tmp_volume ) {
1409 if ( debug_level >= 1 ) {
1410 trace.info() << "(" << i << ")"
1411 << " new_volume = " << tmp_volume
1412 << " > " << best_volume << std::endl;
1413 }
1414 best = tmp;
1415 best_volume = tmp_volume;
1416 }
1417 if ( i >= nbtries && best_volume > 0 )
1418 return best;
1419 }
1420 // If not found, we adopt a deterministic algorithm based on Gauss reduction.
1422 if ( debug_level >= 1 )
1423 trace.info() << "[QuickHull::pickInitialSimplex] #affine subset = " << best.size() << std::endl;
1424 return ( best.size() == (dimension+1) ) ? best : IndexRange();
1425 }
static std::vector< Size > affineSubset(const std::vector< TInputPoint > &X, const double tolerance=1e-12)
static IndexRange pickIntegers(const Size d, const Size n)
Definition QuickHull.h:1430

References DGtal::AffineGeometry< TPoint >::affineSubset(), DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::dimension, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::pickIntegers(), DGtal::QuickHull< TKernel >::points, and DGtal::trace.

Referenced by DGtal::QuickHull< TKernel >::computeInitialSimplex().

◆ pickIntegers()

template<typename TKernel >
static IndexRange DGtal::QuickHull< TKernel >::pickIntegers ( const Size  d,
const Size  n 
)
inlinestaticprotected
Returns
a vector of d distinct integers in {0, 1, ..., n-1} randomly chosen.
Parameters
[in]dthe number of returned integers
[in]nthe range of possible integers {0, 1, ..., n-1}

Definition at line 1430 of file QuickHull.h.

1431 {
1432 IndexRange result( d );
1433 bool distinct = false;
1434 while ( ! distinct )
1435 {
1436 distinct = true;
1437 for ( Index i = 0; i < d; i++ ) result[ i ] = rand() % n;
1438 std::sort( result.begin(), result.end() );
1439 for ( Index i = 1; distinct && i < d; i++ )
1440 distinct = result[ i-1 ] != result[ i ];
1441 }
1442 return result;
1443 }

Referenced by DGtal::QuickHull< TKernel >::pickInitialSimplex().

◆ pointsOnRidge()

template<typename TKernel >
IndexRange DGtal::QuickHull< TKernel >::pointsOnRidge ( const Ridge R) const
inlineprotected
Parameters
[in]Ra ridge between two facets (a pair of facets).
Returns
the points lying on a ridge.

Definition at line 1303 of file QuickHull.h.

1304 {
1305 // Slightly faster to look at points instead at true geometry.
1306 IndexRange result;
1307 const Facet& f1 = facets[ R.first ];
1308 const Facet& f2 = facets[ R.second ];
1309 std::set_intersection( f1.on_set.cbegin(), f1.on_set.cend(),
1310 f2.on_set.cbegin(), f2.on_set.cend(),
1311 std::back_inserter( result ) );
1312 return result;
1313 }

References DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::Facet::on_set, and DGtal::R.

Referenced by DGtal::QuickHull< TKernel >::processFacet().

◆ processFacet()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::processFacet ( std::queue< Index > &  Q)
inlineprotected

Process top facet in queue Q as in Quickhull algorithm, and updates the queue accordingly (top facet poped, new facets pushed).

Parameters
[in,out]Qa queue of facet index to process. which is updated by the method.
Returns
'true' if there is still work to do, 'false' when finished
Note
Core of Quickhull algorithm.

Definition at line 945 of file QuickHull.h.

946 {
947 // If Q empty, we are done
948 if ( Q.empty() ) return false;
949 Index F = Q.front();
950 Q.pop();
951 // If F is already deleted, proceed to next in queue.
952 if ( deleted_facets.count( F ) ) return true;
953 // Take car of current facet.
954 const Facet& facet = facets[ F ];
955 if ( debug_level >= 3 ) {
956 trace.info() << "---------------------------------------------" << std::endl;
957 trace.info() << "---- ACTIVE FACETS---------------------------" << std::endl;
958 bool ok = true;
959 for ( Index i = 0; i < facets.size(); i++ )
960 if ( ! deleted_facets.count( i ) ) {
961 trace.info() << "- facet " << i << " ";
962 facets[ i ].display( trace.info() );
963 ok = ok && checkFacet( i );
964 }
965 if ( ! ok ) { // should not happen.
967 return false;
968 }
969 }
970 if ( debug_level >= 2 ) {
971 trace.info() << "---------------------------------------------" << std::endl;
972 trace.info() << "Processing facet " << F << " ";
973 facet.display( trace.info() );
974 }
975 if ( facet.outside_set.empty() ) return true;
976 // Selects furthest vertex
977 Index furthest_v = facet.outside_set[ 0 ];
978 auto furthest_h = height( facet, points[ furthest_v ] );
979 for ( Index v = 1; v < facet.outside_set.size(); v++ ) {
980 auto h = height( facet, points[ facet.outside_set[ v ] ] );
981 if ( h > furthest_h ) {
982 furthest_h = h;
983 furthest_v = facet.outside_set[ v ];
984 }
985 }
986 const Point& p = points[ furthest_v ];
987 // Extracts Visible facets V and Horizon Ridges H
988 std::vector< Index > V; // visible facets
989 std::set< Index > M; // marked facets (are in E or were in E)
990 std::queue< Index > E; // queue to extract visible facets
991 std::vector< Ridge > H; // visible facets
992 E.push ( F );
993 M.insert( F );
994 while ( ! E.empty() ) {
995 Index G = E.front(); E.pop();
996 V.push_back( G );
997 for ( auto& N : facets[ G ].neighbors ) {
998 if ( aboveOrOn( facets[ N ], p ) ) {
999 if ( M.count( N ) ) continue;
1000 E.push( N );
1001 } else {
1002 H.push_back( { G, N } );
1003 }
1004 M.insert( N );
1005 }
1006 } // while ( ! E.empty() )
1007 if ( debug_level >= 1 ) {
1008 trace.info() << "#Visible=" << V.size() << " #Horizon=" << H.size()
1009 << " furthest_v=" << furthest_v << std::endl;
1010 }
1011 // Create new facets
1012 IndexRange new_facets;
1013 // For each ridge R in H
1014 for ( Index i = 0; i < H.size(); i++ )
1015 {
1016 // Create a new facet from R and p
1017 IndexRange ridge = pointsOnRidge( H[ i ] );
1018 if ( debug_level >= 3 ) {
1019 trace.info() << "Ridge (" << H[i].first << "," << H[i].second << ") = {";
1020 for ( auto&& r : ridge ) trace.info() << " " << r;
1021 trace.info() << " } furthest_v=" << furthest_v << std::endl;
1022 }
1023 IndexRange base( 1 + ridge.size() );
1024 Index j = 0;
1025 base[ j++ ] = furthest_v;
1026 for ( auto&& v : ridge ) base[ j++ ] = v;
1027 if ( j < dimension ) {
1028 trace.error() << "Bad ridge between " << std::endl
1029 << "- facet " << H[i].first << " ";
1030 facets[ H[i].first ].display( trace.error() );
1031 trace.error() << "- facet " << H[i].second << " ";
1032 facets[ H[i].second ].display( trace.error() );
1033 }
1034 Index nf = newFacet();
1035 new_facets.push_back( nf );
1036 facets[ nf ] = makeFacet( base, facets[ H[i].first ].below );
1037 facets[ nf ].on_set = IndexRange { base.cbegin(), base.cend() };
1038 std::sort( facets[ nf ].on_set.begin(), facets[ nf ].on_set.end() );
1039 makeNeighbors( nf, H[ i ].second );
1040 if ( debug_level >= 3 ) {
1041 trace.info() << "* New facet " << nf << " ";
1042 facets[ nf ].display( trace.info() );
1043 }
1044 // Checks that the facet is not parallel to another in the Horizon
1045 for ( Index k = 0; k < new_facets.size() - 1; k++ )
1046 if ( areFacetsParallel( new_facets[ k ], nf ) ) {
1047 if ( debug_level >= 1 ) {
1048 trace.info() << "Facets " << new_facets[ k ] << " and " << nf
1049 << " are parallel => merge." << std::endl;
1050 }
1051 mergeParallelFacets( new_facets[ k ], nf );
1052 new_facets.pop_back();
1053 deleteFacet( nf );
1054 if ( debug_level >= 3 ) {
1055 facets[ new_facets[ k ] ].display( trace.info() );
1056 }
1057 }
1058 }
1059 // For each new facet
1060 for ( Index i = 0; i < new_facets.size(); i++ )
1061 { // link the new facet to its neighbors
1062 for ( Index j = i + 1; j < new_facets.size(); j++ )
1063 {
1064 const Index nfi = new_facets[ i ];
1065 const Index nfj = new_facets[ j ];
1066 if ( areFacetsNeighbor( nfi, nfj ) )
1067 makeNeighbors( nfi, nfj );
1068 }
1069 }
1070 // Extracts all outside points from visible facets V
1071 IndexRange outside_pts;
1072 for ( auto&& vf : V ) {
1073 for ( auto&& v : facets[ vf ].outside_set ) {
1074 if ( v != furthest_v ) {
1075 outside_pts.push_back( v );
1076 assignment[ v ] = UNASSIGNED;
1077 }
1078 }
1079 }
1080 // For each new facet F'
1081 for ( Index i = 0; i < new_facets.size(); i++ ) {
1082 Facet& Fp = facets[ new_facets[ i ] ];
1083 Index max_j = outside_pts.size();
1084 for ( Index j = 0; j < max_j; ) {
1085 const Index v = outside_pts[ j ];
1086 if ( above( Fp, points[ v ] ) ) {
1087 Fp.outside_set.push_back( v );
1088 assignment[ v ] = new_facets[ i ];
1089 outside_pts[ j ] = outside_pts.back();
1090 outside_pts.pop_back();
1091 max_j--;
1092 } else j++;
1093 }
1094 if ( debug_level >= 3 ) {
1095 trace.info() << "- New facet " << new_facets[ i ] << " ";
1096 Fp.display( trace.info() );
1097 }
1098 }
1099 // Update processed points
1100 processed_points.push_back( furthest_v );
1101 for ( auto v : outside_pts ) processed_points.push_back( v );
1102
1103 // Delete the facets in V
1104 for ( auto&& v : V ) {
1105 if ( debug_level >= 2 ) {
1106 trace.info() << "Delete facet " << v << " ";
1107 facets[ v ].display( trace.info() );
1108 }
1109 deleteFacet( v );
1110 }
1111
1112 // Add new facets to queue
1113 for ( Index i = 0; i < new_facets.size(); i++ )
1114 Q.push( new_facets[ i ] );
1115 if ( debug_level >= 1 ) {
1116 trace.info() << "#facets=" << facets.size()
1117 << " #deleted=" << deleted_facets.size() << std::endl;
1118 }
1119
1120 // Checks that everything is ok.
1121 if ( debug_level >= 1 ) {
1122 trace.info() << "[CHECK INVARIANT] " << processed_points.size()
1123 << " / " << points.size() << " points processed." << std::endl;
1124 bool okh = checkHull();
1125 if ( ! okh )
1126 trace.error() << "[computeFacet] Invalid convex hull" << std::endl;
1127 bool okf = checkFacets();
1128 if ( ! okf )
1129 trace.error() << "[computeFacet] Invalid facets" << std::endl;
1130 if ( ! ( okh && okf ) ) myStatus = Status::InvalidConvexHull;
1131 }
1132
1133 return status() == Status::SimplexCompleted;
1134 }
IndexRange pointsOnRidge(const Ridge &R) const
Definition QuickHull.h:1303
bool areFacetsParallel(const Index if1, const Index if2) const
Definition QuickHull.h:1249
void deleteFacet(Index f)
Definition QuickHull.h:1194
InternalScalar height(const Facet &F, const Point &p) const
Definition QuickHull.h:851
Index mergeParallelFacets(const Index if1, const Index if2)
Definition QuickHull.h:1226

References DGtal::QuickHull< TKernel >::above(), DGtal::QuickHull< TKernel >::aboveOrOn(), DGtal::QuickHull< TKernel >::areFacetsNeighbor(), DGtal::QuickHull< TKernel >::areFacetsParallel(), DGtal::QuickHull< TKernel >::assignment, DGtal::QuickHull< TKernel >::checkFacet(), DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::deleted_facets, DGtal::QuickHull< TKernel >::deleteFacet(), DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::Facet::display(), DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::H, DGtal::QuickHull< TKernel >::height(), DGtal::Trace::info(), DGtal::QuickHull< TKernel >::InvalidConvexHull, DGtal::QuickHull< TKernel >::makeFacet(), DGtal::QuickHull< TKernel >::makeNeighbors(), DGtal::QuickHull< TKernel >::mergeParallelFacets(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::newFacet(), DGtal::QuickHull< TKernel >::Facet::outside_set, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::pointsOnRidge(), DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::SimplexCompleted, DGtal::QuickHull< TKernel >::status(), DGtal::trace, and DGtal::QuickHull< TKernel >::UNASSIGNED.

Referenced by DGtal::QuickHull< TKernel >::computeFacets().

◆ renumberInfiniteFacets()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::renumberInfiniteFacets ( )
inlineprotected

Determine infinite facets and renumber them so that finite facets come first and infinite facets come after.

Definition at line 902 of file QuickHull.h.

903 {
904 nb_finite_facets = facets.size();
906 if ( ! kernel.hasInfiniteFacets() ) return;
907 IndexRange renumbering( facets.size() );
908 Index i = 0;
909 Index k = facets.size();
910 Index j = 0;
911 for ( auto& l : renumbering ) {
912 if ( ! kernel.isHalfSpaceFacetInfinite( facets[ j ].H ) ) l = i++;
913 else l = --k;
914 j++;
915 }
916 if ( i != k )
917 trace.error() << "[Quickhull::renumberInfiniteFacets]"
918 << " Error renumbering infinite facets "
919 << " up finite=" << i << " low infinite=" << k << std::endl;
920 std::vector< Facet > new_facets( facets.size() );
921 for ( Index f = 0; f < facets.size(); f++ )
922 new_facets[ renumbering[ f ] ].swap( facets[ f ] );
923 facets.swap( new_facets );
924 for ( auto& F : facets ) {
925 for ( auto& N : F.neighbors ) {
926 N = renumbering[ N ];
927 }
928 }
929 // Assign correct number of facets.
932 }
void swap(UnorderedSetByBlock< Key, TSplitter, Hash, KeyEqual, UnorderedMapAllocator > &s1, UnorderedSetByBlock< Key, TSplitter, Hash, KeyEqual, UnorderedMapAllocator > &s2) noexcept

References DGtal::Trace::error(), DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::nb_finite_facets, DGtal::QuickHull< TKernel >::nb_infinite_facets, DGtal::swap(), and DGtal::trace.

Referenced by DGtal::QuickHull< TKernel >::computeVertices().

◆ selfDisplay()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::selfDisplay ( std::ostream &  out) const
inline

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.

Definition at line 1518 of file QuickHull.h.

1519 {
1520 out << "[QuickHull " << dimension << "D"
1521 // << " status=" << status()
1522 << " #P=" << nbPoints();
1524 out << " #F=" << nbFacets();
1526 out << " #V=" << nbVertices();
1527 out << "]";
1528 // if ( status() >= Status::FacetsCompleted && status() <= Status::AllCompleted
1529 // && nbFacets() == 24 ) {
1530 // for ( auto f : facets ) f.display( out );
1531 // for ( auto v : v2p ) out << points[ v2p[ v ] ] << std::endl;
1532 // }
1533 }

References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::QuickHull< TKernel >::nbFacets(), DGtal::QuickHull< TKernel >::nbPoints(), DGtal::QuickHull< TKernel >::nbVertices(), DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.

◆ setInitialSimplex()

template<typename TKernel >
bool DGtal::QuickHull< TKernel >::setInitialSimplex ( const IndexRange full_splx)
inline

Sets the initial full dimensional simplex

Precondition
status() must be Status::InputInitialized
Parameters
full_splxa dimension+1-simplex specified as indices in the vector of input point.
Returns
'true' iff this initial simplex is full dimensional, the object has then the status SimplexCompleted, otherwise returns 'false' and the status is NotFullDimensional.

Definition at line 411 of file QuickHull.h.

412 {
413 if ( status() != Status::InputInitialized ) return false;
414 if ( full_splx.size() != dimension + 1 )
415 {
416 trace.error() << "[QuickHull::setInitialSimplex]"
417 << " not a full dimensional simplex" << std::endl;
419 return false;
420 }
422 for ( Index j = 0; j < dimension; ++j )
423 splx[ j ] = full_splx[ j ];
424 const auto H = kernel.compute( points, splx, full_splx.back() );
425 const auto volume = kernel.volume( H, points[ full_splx.back() ] );
426 if ( volume > 0 )
427 return computeSimplexConfiguration( full_splx );
429 return false;
430 }

References DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::dimension, DGtal::Trace::error(), DGtal::H, DGtal::QuickHull< TKernel >::InputInitialized, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::NotFullDimensional, DGtal::QuickHull< TKernel >::points, DGtal::QuickHull< TKernel >::status(), and DGtal::trace.

◆ setInput()

template<typename TKernel >
template<typename InputPoint >
bool DGtal::QuickHull< TKernel >::setInput ( const std::vector< InputPoint > &  input_points,
bool  remove_duplicates = true 
)
inline

Sets the input data for the QuickHull convex hull algorithm, which is a range of points.

Template Parameters
InputPointa model of point that is convertible to Point datatype.
Parameters
[in]input_pointsthe range of input points.
[in]remove_duplicatesshould be set to 'true' if the input data has duplicates.
Returns
'true' if the object is successfully initialized, status must be Status::InputInitialized, 'false' otherwise.
Examples
geometry/tools/exampleLatticeBallQuickHull3D.cpp.

Definition at line 383 of file QuickHull.h.

385 {
386 Clock tic;
387 tic.startClock();
388 clear();
389 timings.clear();
390 kernel.makeInput( points, input2comp, comp2input,
391 input_points, remove_duplicates );
392 timings.push_back( tic.stopClock() );
393 if ( points.size() <= dimension ) {
395 return false;
396 }
398 return true;
399 }
void clear()
Clears the object.
Definition QuickHull.h:271

References DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::comp2input, DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::input2comp, DGtal::QuickHull< TKernel >::InputInitialized, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::NotFullDimensional, DGtal::QuickHull< TKernel >::points, DGtal::Clock::startClock(), tic(), and DGtal::QuickHull< TKernel >::timings.

Referenced by DGtal::detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K >::compute(), and main().

◆ status()

◆ unmakeNeighbors()

template<typename TKernel >
void DGtal::QuickHull< TKernel >::unmakeNeighbors ( const Index  if1,
const Index  if2 
)
inlineprotected

Makes two distinct facets if1 and if2 no more neighbors

Parameters
if1a valid facet index
if2a valid facet index

Definition at line 1214 of file QuickHull.h.

1215 {
1216 facets[ if1 ].subNeighbor( if2 );
1217 facets[ if2 ].subNeighbor( if1 );
1218 }

References DGtal::QuickHull< TKernel >::facets.

Field Documentation

◆ assignment

template<typename TKernel >
std::vector< Index > DGtal::QuickHull< TKernel >::assignment

◆ comp2input

template<typename TKernel >
IndexRange DGtal::QuickHull< TKernel >::comp2input

the injective mapping between the convex hull point range and the input range.

Definition at line 809 of file QuickHull.h.

Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::setInput().

◆ debug_level

◆ deleted_facets

◆ dimension

◆ facet_counter

template<typename TKernel >
std::vector< Size > DGtal::QuickHull< TKernel >::facet_counter

Counts the number of facets with a given number of vertices.

Definition at line 829 of file QuickHull.h.

Referenced by DGtal::QuickHull< TKernel >::computeVertices().

◆ facets

◆ input2comp

template<typename TKernel >
IndexRange DGtal::QuickHull< TKernel >::input2comp

the surjective mapping between the input range and the output range used for convex hull computation.

Definition at line 806 of file QuickHull.h.

Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::setInput().

◆ kernel

◆ myStatus

◆ nb_finite_facets

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nb_finite_facets

Number of finite facets.

Definition at line 823 of file QuickHull.h.

Referenced by DGtal::QuickHull< TKernel >::nbFiniteFacets(), and DGtal::QuickHull< TKernel >::renumberInfiniteFacets().

◆ nb_infinite_facets

template<typename TKernel >
Size DGtal::QuickHull< TKernel >::nb_infinite_facets

Number of infinite facets (!= 0 only for specific kernels)

Definition at line 825 of file QuickHull.h.

Referenced by DGtal::QuickHull< TKernel >::nbInfiniteFacets(), and DGtal::QuickHull< TKernel >::renumberInfiniteFacets().

◆ p2v

◆ points

◆ processed_points

◆ timings

template<typename TKernel >
std::vector< double > DGtal::QuickHull< TKernel >::timings

◆ v2p


The documentation for this struct was generated from the following file: