DGtal
1.5.beta
|
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>
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< Index > | IndexRange |
typedef Kernel::HalfSpace | HalfSpace |
typedef Kernel::CombinatorialPlaneSimplex | CombinatorialPlaneSimplex |
typedef std::pair< Index, Index > | Ridge |
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. More... | |
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. More... | |
int | debug_level |
debug_level from 0:no to 2 More... | |
std::vector< Point > | points |
the set of points, indexed as in the array. More... | |
IndexRange | input2comp |
IndexRange | comp2input |
IndexRange | processed_points |
Points already processed (and within the convex hull). More... | |
std::vector< Index > | assignment |
assignment of points to facets More... | |
std::vector< Facet > | facets |
the current set of facets. More... | |
std::set< Index > | deleted_facets |
set of deleted facets More... | |
IndexRange | p2v |
point index -> vertex index (or UNASSIGNED) More... | |
IndexRange | v2p |
vertex index -> point index More... | |
Size | nb_finite_facets |
Number of finite facets. More... | |
Size | nb_infinite_facets |
Number of infinite facets (!= 0 only for specific kernels) More... | |
std::vector< double > | timings |
Timings of the different phases: 0: init, 1: facets, 2: vertices. More... | |
std::vector< Size > | facet_counter |
Counts the number of facets with a given number of vertices. More... | |
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 | cleanInfiniteFacets () |
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) |
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).
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.
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.TKernel | any type of QuickHull kernel, like ConvexHullIntegralKernel. |
Definition at line 139 of file QuickHull.h.
typedef Kernel::CombinatorialPlaneSimplex DGtal::QuickHull< TKernel >::CombinatorialPlaneSimplex |
Definition at line 151 of file QuickHull.h.
typedef Kernel::HalfSpace DGtal::QuickHull< TKernel >::HalfSpace |
Definition at line 150 of file QuickHull.h.
typedef std::size_t DGtal::QuickHull< TKernel >::Index |
Definition at line 146 of file QuickHull.h.
typedef std::vector< Index > DGtal::QuickHull< TKernel >::IndexRange |
Definition at line 149 of file QuickHull.h.
typedef Kernel::InternalScalar DGtal::QuickHull< TKernel >::InternalScalar |
Definition at line 145 of file QuickHull.h.
typedef TKernel DGtal::QuickHull< TKernel >::Kernel |
Definition at line 141 of file QuickHull.h.
typedef Kernel::CoordinatePoint DGtal::QuickHull< TKernel >::Point |
Definition at line 142 of file QuickHull.h.
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 236 of file QuickHull.h.
typedef Kernel::CoordinateScalar DGtal::QuickHull< TKernel >::Scalar |
Definition at line 144 of file QuickHull.h.
typedef std::size_t DGtal::QuickHull< TKernel >::Size |
Definition at line 147 of file QuickHull.h.
typedef Kernel::CoordinateVector DGtal::QuickHull< TKernel >::Vector |
Definition at line 143 of file QuickHull.h.
anonymous enum |
Label for points that are not assigned to any facet.
Enumerator | |
---|---|
UNASSIGNED |
Definition at line 155 of file QuickHull.h.
|
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 239 of file QuickHull.h.
|
inline |
Default constructor
[in] | K | a kernel for computing facet geometries. |
[in] | dbg | the trace level, from 0 (no) to 3 (very verbose). |
Definition at line 259 of file QuickHull.h.
|
inlineprotected |
F | any valid facet |
p | any point |
Definition at line 856 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.
Referenced by DGtal::QuickHull< TKernel >::checkFacet(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
F | any valid facet |
p | any point |
Definition at line 862 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.
Referenced by DGtal::QuickHull< TKernel >::checkFacet(), and DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
Checks if two facets are neighbors by looking at the points on their boundary.
if1 | a valid facet index |
if2 | a valid facet index |
Definition at line 1292 of file QuickHull.h.
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().
|
inlineprotected |
Checks if two distinct facets are parallel (i.e. should be merged).
if1 | a valid facet index |
if2 | a valid facet index |
Definition at line 1273 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::Facet::H, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::on(), and DGtal::QuickHull< TKernel >::Facet::on_set.
Referenced by DGtal::QuickHull< TKernel >::processFacet().
DGtal::QuickHull< TKernel >::BOOST_STATIC_ASSERT | ( | (Point::dimension==Vector::dimension) | ) |
|
inline |
Global validity check of the convex hull after processing.
Definition at line 715 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::QuickHull< TKernel >::processed_points, DGtal::QuickHull< TKernel >::status(), DGtal::trace, and DGtal::Trace::warning().
|
inlineprotected |
Definition at line 1161 of file QuickHull.h.
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, and DGtal::trace.
Referenced by DGtal::QuickHull< TKernel >::checkFacets(), and DGtal::QuickHull< TKernel >::processFacet().
|
inline |
Definition at line 743 of file QuickHull.h.
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().
|
inline |
Definition at line 762 of file QuickHull.h.
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 >::processed_points, and DGtal::trace.
Referenced by DGtal::QuickHull< TKernel >::check(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
Cleans and renumber the facets so that no one belongs to deleted_facets.
Definition at line 873 of file QuickHull.h.
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().
|
inlineprotected |
Removes infinite facets and renumber the facets in a consecutive way.
Definition at line 901 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::kernel, and DGtal::QuickHull< TKernel >::UNASSIGNED.
|
inline |
Clears the object.
Definition at line 270 of file QuickHull.h.
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 >::processed_points, DGtal::QuickHull< TKernel >::timings, DGtal::QuickHull< TKernel >::Uninitialized, and DGtal::QuickHull< TKernel >::v2p.
Referenced by DGtal::QuickHull< TKernel >::setInput().
|
inline |
Computes the convex hull of the given range of points until a specified target:
[in] | target | the computation target in Status::SimplexCompleted, Status::FacetsCompleted, Status::VerticesCompleted. |
Definition at line 460 of file QuickHull.h.
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::QuickHull< TKernel >::status(), tic(), DGtal::QuickHull< TKernel >::timings, and DGtal::QuickHull< TKernel >::VerticesCompleted.
Referenced by main().
|
inline |
Computes the facets of the convex hull using Quickhull algorithm. If everything went well, the status is Status::FacetsCompleted afterwards.
Definition at line 523 of file QuickHull.h.
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().
|
inline |
Computes the initial full dimensional simplex from the input data.
Definition at line 504 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::myStatus, DGtal::QuickHull< TKernel >::NotFullDimensional, and DGtal::QuickHull< TKernel >::pickInitialSimplex().
Referenced by DGtal::QuickHull< TKernel >::computeConvexHull().
|
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.
Definition at line 1454 of file QuickHull.h.
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 >::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().
|
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.
Definition at line 553 of file QuickHull.h.
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 >::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().
|
inlineprotected |
Deletes the given facet f.
f | a valid facet index |
Definition at line 1218 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::deleted_facets, and DGtal::QuickHull< TKernel >::facets.
Referenced by DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
Filters each vertex on the facet f to keep only the ones that are on or below the neighboring facets
[in] | f | any valid facet index. |
Definition at line 1377 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, and DGtal::QuickHull< TKernel >::on().
|
inline |
This methods return the halfspaces corresponding to each facet. It is useful to build a polytope.
[out] | facet_halfspaces | the range giving for each facet its halfspace as a pair (normal, intercept). |
Definition at line 690 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::H, DGtal::QuickHull< TKernel >::nbFacets(), and DGtal::QuickHull< TKernel >::status().
|
inline |
[out] | facet_vertices | the range giving for each facet the indices of its vertices. |
Definition at line 666 of file QuickHull.h.
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 main().
|
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)
[out] | point_to_vertex | the 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) |
Definition at line 651 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::p2v, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.
|
inline |
Gets for each vertex its index in the input range of points.
[out] | vertex_to_point | the range giving for each vertex its index in the input range of points. |
Definition at line 632 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.
|
inline |
Gets the positions of the convex hull vertices.
OutputPoint | a model of point such that the Point datatype is convertible to it. |
[out] | vertex_positions | the range of vertex positions. |
Definition at line 612 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::kernel, DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.
Referenced by main().
|
inlineprotected |
F | any valid facet |
p | any point |
Definition at line 850 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::Facet::H, and DGtal::QuickHull< TKernel >::kernel.
Referenced by DGtal::QuickHull< TKernel >::processFacet().
|
inline |
Checks the validity/consistency of the object.
Definition at line 1540 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::Uninitialized.
|
inlineprotected |
Builds a facet from a base convex set of at least d different points and a point below.
[in] | base | a range containing at least d distinct points in general position. |
[in] | below | a point below the hyperplane containing simplex. |
Definition at line 1317 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::dimension, and DGtal::QuickHull< TKernel >::kernel.
Referenced by DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), and DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
Makes two distinct facets if1 and if2 as neighbors
if1 | a valid facet index |
if2 | a valid facet index |
Definition at line 1229 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets.
Referenced by DGtal::QuickHull< TKernel >::mergeParallelFacets(), and DGtal::QuickHull< TKernel >::processFacet().
|
inline |
Definition at line 287 of file QuickHull.h.
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 >::processed_points, DGtal::QuickHull< TKernel >::timings, and DGtal::QuickHull< TKernel >::v2p.
|
inlineprotected |
Merge the two facets and return the index of the second one, which is the one deleted.
if1 | a valid facet index |
if2 | a valid facet index |
Definition at line 1250 of file QuickHull.h.
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().
|
inline |
Definition at line 328 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, and DGtal::QuickHull< TKernel >::status().
Referenced by DGtal::QuickHull< TKernel >::getFacetHalfSpaces(), DGtal::QuickHull< TKernel >::getFacetVertices(), main(), and DGtal::QuickHull< TKernel >::selfDisplay().
|
inline |
Definition at line 346 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::nb_finite_facets, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.
|
inline |
Definition at line 355 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::nb_infinite_facets, DGtal::QuickHull< TKernel >::status(), and DGtal::QuickHull< TKernel >::VerticesCompleted.
|
inline |
Definition at line 323 of file QuickHull.h.
Referenced by main(), and DGtal::QuickHull< TKernel >::selfDisplay().
|
inline |
Definition at line 337 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.
Referenced by main(), and DGtal::QuickHull< TKernel >::selfDisplay().
|
inlineprotected |
Definition at line 1208 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets.
Referenced by DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
F | any valid facet |
p | any point |
Definition at line 868 of file QuickHull.h.
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().
|
inlineprotected |
Given a facet index f, return its points oriented consistently with respect to the normal.
f | any valid facet index |
Definition at line 1348 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, and DGtal::QuickHull< TKernel >::Facet::on_set.
Referenced by DGtal::QuickHull< TKernel >::getFacetVertices().
|
inlineprotected |
Definition at line 1397 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::debug_level, DGtal::QuickHull< TKernel >::dimension, DGtal::Trace::info(), DGtal::QuickHull< TKernel >::kernel, max(), DGtal::QuickHull< TKernel >::pickIntegers(), and DGtal::trace.
Referenced by DGtal::QuickHull< TKernel >::computeInitialSimplex().
|
inlinestaticprotected |
{0, 1, ..., n-1}
randomly chosen. [in] | d | the number of returned integers |
[in] | n | the range of possible integers {0, 1, ..., n-1} |
Definition at line 1431 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::pickInitialSimplex().
|
inlineprotected |
[in] | R | a ridge between two facets (a pair of facets). |
Definition at line 1327 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets, and DGtal::QuickHull< TKernel >::Facet::on_set.
Referenced by DGtal::QuickHull< TKernel >::processFacet().
|
inlineprotected |
Process top facet in queue Q as in Quickhull algorithm, and updates the queue accordingly (top facet poped, new facets pushed).
[in,out] | Q | a queue of facet index to process. which is updated by the method. |
Definition at line 969 of file QuickHull.h.
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 >::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().
|
inlineprotected |
Determine infinite facets and renumber them so that finite facets come first and infinite facets come after.
Definition at line 926 of file QuickHull.h.
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().
|
inline |
Writes/Displays the object on an output stream.
out | the output stream where the object is written. |
Definition at line 1519 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::AllCompleted, DGtal::QuickHull< TKernel >::dimension, DGtal::QuickHull< TKernel >::facets, DGtal::QuickHull< TKernel >::FacetsCompleted, DGtal::QuickHull< TKernel >::nbFacets(), DGtal::QuickHull< TKernel >::nbPoints(), DGtal::QuickHull< TKernel >::nbVertices(), DGtal::QuickHull< TKernel >::status(), DGtal::QuickHull< TKernel >::v2p, and DGtal::QuickHull< TKernel >::VerticesCompleted.
|
inline |
Sets the initial full dimensional simplex
full_splx | a dimension+1 -simplex specified as indices in the vector of input point. |
Definition at line 410 of file QuickHull.h.
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 >::status(), and DGtal::trace.
|
inline |
Sets the input data for the QuickHull convex hull algorithm, which is a range of points.
InputPoint | a model of point that is convertible to Point datatype. |
[in] | input_points | the range of input points. |
[in] | remove_duplicates | should be set to 'true' if the input data has duplicates. |
Definition at line 382 of file QuickHull.h.
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, tic(), and DGtal::QuickHull< TKernel >::timings.
Referenced by main().
|
inline |
Definition at line 266 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::myStatus.
Referenced by DGtal::QuickHull< TKernel >::check(), DGtal::QuickHull< TKernel >::computeConvexHull(), DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::getFacetHalfSpaces(), DGtal::QuickHull< TKernel >::getFacetVertices(), DGtal::QuickHull< TKernel >::getPoint2Vertex(), DGtal::QuickHull< TKernel >::getVertex2Point(), DGtal::QuickHull< TKernel >::getVertexPositions(), DGtal::QuickHull< TKernel >::isValid(), DGtal::QuickHull< TKernel >::nbFacets(), DGtal::QuickHull< TKernel >::nbFiniteFacets(), DGtal::QuickHull< TKernel >::nbInfiniteFacets(), DGtal::QuickHull< TKernel >::nbVertices(), DGtal::QuickHull< TKernel >::processFacet(), DGtal::QuickHull< TKernel >::selfDisplay(), and DGtal::QuickHull< TKernel >::setInitialSimplex().
|
inlineprotected |
Makes two distinct facets if1 and if2 no more neighbors
if1 | a valid facet index |
if2 | a valid facet index |
Definition at line 1238 of file QuickHull.h.
References DGtal::QuickHull< TKernel >::facets.
std::vector< Index > DGtal::QuickHull< TKernel >::assignment |
assignment of points to facets
Definition at line 812 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::processFacet().
IndexRange DGtal::QuickHull< TKernel >::comp2input |
the injective mapping between the convex hull point range and the input range.
Definition at line 808 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::setInput().
int DGtal::QuickHull< TKernel >::debug_level |
debug_level from 0:no to 2
Definition at line 800 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::pickInitialSimplex(), and DGtal::QuickHull< TKernel >::processFacet().
std::set< Index > DGtal::QuickHull< TKernel >::deleted_facets |
set of deleted facets
Definition at line 816 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::cleanFacets(), DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::deleteFacet(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::processFacet().
|
static |
Definition at line 152 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::areFacetsNeighbor(), DGtal::QuickHull< TKernel >::checkFacet(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::filterVerticesOnFacet(), DGtal::QuickHull< TKernel >::makeFacet(), DGtal::QuickHull< TKernel >::orientedFacetPoints(), DGtal::QuickHull< TKernel >::pickInitialSimplex(), DGtal::QuickHull< TKernel >::processFacet(), DGtal::QuickHull< TKernel >::selfDisplay(), DGtal::QuickHull< TKernel >::setInitialSimplex(), and DGtal::QuickHull< TKernel >::setInput().
std::vector< Size > DGtal::QuickHull< TKernel >::facet_counter |
Counts the number of facets with a given number of vertices.
Definition at line 828 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::computeVertices().
std::vector< Facet > DGtal::QuickHull< TKernel >::facets |
the current set of facets.
Definition at line 814 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::areFacetsNeighbor(), DGtal::QuickHull< TKernel >::areFacetsParallel(), DGtal::QuickHull< TKernel >::checkFacet(), DGtal::QuickHull< TKernel >::checkFacets(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::cleanFacets(), DGtal::QuickHull< TKernel >::cleanInfiniteFacets(), DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::deleteFacet(), DGtal::QuickHull< TKernel >::filterVerticesOnFacet(), DGtal::QuickHull< TKernel >::getFacetHalfSpaces(), DGtal::QuickHull< TKernel >::makeNeighbors(), DGtal::QuickHull< TKernel >::memory(), DGtal::QuickHull< TKernel >::mergeParallelFacets(), DGtal::QuickHull< TKernel >::nbFacets(), DGtal::QuickHull< TKernel >::newFacet(), DGtal::QuickHull< TKernel >::orientedFacetPoints(), DGtal::QuickHull< TKernel >::pointsOnRidge(), DGtal::QuickHull< TKernel >::processFacet(), DGtal::QuickHull< TKernel >::renumberInfiniteFacets(), DGtal::QuickHull< TKernel >::selfDisplay(), and DGtal::QuickHull< TKernel >::unmakeNeighbors().
IndexRange DGtal::QuickHull< TKernel >::input2comp |
the surjective mapping between the input range and the output range used for convex hull computation.
Definition at line 805 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::setInput().
|
mutable |
Kernel that is duplicated for computing facet geometries.
Definition at line 798 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::above(), DGtal::QuickHull< TKernel >::aboveOrOn(), DGtal::QuickHull< TKernel >::areFacetsParallel(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::cleanInfiniteFacets(), DGtal::QuickHull< TKernel >::getVertexPositions(), DGtal::QuickHull< TKernel >::height(), DGtal::QuickHull< TKernel >::makeFacet(), DGtal::QuickHull< TKernel >::memory(), DGtal::QuickHull< TKernel >::on(), DGtal::QuickHull< TKernel >::pickInitialSimplex(), DGtal::QuickHull< TKernel >::renumberInfiniteFacets(), DGtal::QuickHull< TKernel >::setInitialSimplex(), and DGtal::QuickHull< TKernel >::setInput().
|
protected |
The status of the object: Uninitialized, InputInitialized, SimplexCompleted, FacetsCompleted, VerticesCompleted, InvalidRidge, InvalidConvexHull, NotFullDimensional
Definition at line 839 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeConvexHull(), DGtal::QuickHull< TKernel >::computeFacets(), DGtal::QuickHull< TKernel >::computeInitialSimplex(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::processFacet(), DGtal::QuickHull< TKernel >::setInitialSimplex(), DGtal::QuickHull< TKernel >::setInput(), and DGtal::QuickHull< TKernel >::status().
Size DGtal::QuickHull< TKernel >::nb_finite_facets |
Number of finite facets.
Definition at line 822 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::nbFiniteFacets(), and DGtal::QuickHull< TKernel >::renumberInfiniteFacets().
Size DGtal::QuickHull< TKernel >::nb_infinite_facets |
Number of infinite facets (!= 0 only for specific kernels)
Definition at line 824 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::nbInfiniteFacets(), and DGtal::QuickHull< TKernel >::renumberInfiniteFacets().
IndexRange DGtal::QuickHull< TKernel >::p2v |
point index -> vertex index (or UNASSIGNED)
Definition at line 818 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::getFacetVertices(), DGtal::QuickHull< TKernel >::getPoint2Vertex(), and DGtal::QuickHull< TKernel >::memory().
std::vector< Point > DGtal::QuickHull< TKernel >::points |
the set of points, indexed as in the array.
Definition at line 802 of file QuickHull.h.
IndexRange DGtal::QuickHull< TKernel >::processed_points |
Points already processed (and within the convex hull).
Definition at line 810 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::check(), DGtal::QuickHull< TKernel >::checkHull(), DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeSimplexConfiguration(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::processFacet().
std::vector< double > DGtal::QuickHull< TKernel >::timings |
Timings of the different phases: 0: init, 1: facets, 2: vertices.
Definition at line 826 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeConvexHull(), main(), DGtal::QuickHull< TKernel >::memory(), and DGtal::QuickHull< TKernel >::setInput().
IndexRange DGtal::QuickHull< TKernel >::v2p |
vertex index -> point index
Definition at line 820 of file QuickHull.h.
Referenced by DGtal::QuickHull< TKernel >::clear(), DGtal::QuickHull< TKernel >::computeVertices(), DGtal::QuickHull< TKernel >::getVertex2Point(), DGtal::QuickHull< TKernel >::getVertexPositions(), DGtal::QuickHull< TKernel >::memory(), DGtal::QuickHull< TKernel >::nbVertices(), and DGtal::QuickHull< TKernel >::selfDisplay().