DGtal  1.5.beta
DGtal::functions::Hull2D Namespace Reference

Hull2D namespace gathers useful functions to compute and return the convex hull or the alpha-shape of a range of 2D points. More...

Enumerations

enum  ThicknessDefinition { HorizontalVerticalThickness , EuclideanThickness }
 The 2 thickness definitions. More...
 

Functions

template<typename Stack , typename Point , typename Predicate >
void updateHullWithStack (Stack &aStack, const Point &aNewPoint, const Predicate &aPredicate)
 Procedure that updates the hull when an extra point aNewPoint is considered: while the last three consecutive points, ie. aNewPoint and the last two points of the container are not oriented such that the predicate aPredicate returns 'true',
the last point of the container is removed. More...
 
template<typename Stack , typename Point , typename Predicate >
void updateHullWithAdaptedStack (Stack aStack, const Point &aNewPoint, const Predicate &aPredicate)
 Procedure that calls Hull2D::updateHullWithStack on a copy of the stack object used to retrieved the hull vertices. Useful when the first argument is a stack adapter returned by a function: it must be copied before being passed by reference in Hull2D::updateHullWithStack. More...
 
template<typename Stack , typename ForwardIterator , typename Predicate >
void buildHullWithStack (Stack &aStack, const ForwardIterator &itb, const ForwardIterator &ite, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time. This technique is called Sklansky's scan, Graham's scan or 3-coins algorithm.
It works for all WEVP [Toussaint and Avis, 1982 : [118]].
More...
 
template<typename Stack , typename ForwardIterator , typename Predicate >
void buildHullWithAdaptedStack (Stack aStack, const ForwardIterator &itb, const ForwardIterator &ite, const Predicate &aPredicate)
 Procedure that calls Hull2D::buildHullWithStack on a copy of the stack object used to retrieved the hull vertices. Useful when the first argument is a stack adapter returned by a function: it must be copied before being passed by reference in Hull2D::buildHullWithStack. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void openGrahamScan (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void closedGrahamScanFromVertex (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void closedGrahamScanFromAnyPoint (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate , typename PolarComparator >
void grahamConvexHullAlgorithm (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate, PolarComparator &aPolarComparator)
 Procedure that retrieves the vertices of the convex hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Graham's algorithm [Graham, 1972 : [63]]. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate , typename PolarComparator >
void grahamConvexHullAlgorithm (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the convex hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Graham's algorithm [Graham, 1972 : [63]]. More...
 
template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void andrewConvexHullAlgorithm (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, const Predicate &aPredicate)
 Procedure that retrieves the vertices of the hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known monotone-chain algorithm due to [Andrew, 1979 : [6]]. More...
 
template<typename ForwardIterator >
double computeHullThickness (const ForwardIterator &itb, const ForwardIterator &ite, const ThicknessDefinition &def)
 Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical or Euclidean distances). It takes as input the vertices of the hull given by the range [itbn, ite). The procedure applies the classic rotating caliper to recover all anti-podal pairs. More...
 
template<typename ForwardIterator , typename TInputPoint >
double computeHullThickness (const ForwardIterator &itb, const ForwardIterator &ite, const ThicknessDefinition &def, TInputPoint &antipodalEdgeP, TInputPoint &antipodalEdgeQ, TInputPoint &antipodalVertexR)
 Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical or Euclidean distances). It takes as input the vertices of the hull given by the range [itbn, ite). The procedure applies the classic rotating caliper to recover all anti-podal pairs. More...
 
template<typename TInputPoint >
double getAngle (const TInputPoint &a, const TInputPoint &b, const TInputPoint &c, const TInputPoint &d)
 
template<typename TInputPoint >
bool isCoLinearOpp (const TInputPoint &a, const TInputPoint &b, const TInputPoint &c, const TInputPoint &d)
 
template<typename TInputPoint >
double getThicknessAntipodalPair (const TInputPoint &p, const TInputPoint &q, const TInputPoint &r, const ThicknessDefinition &def)
 
template<typename TInputPoint >
double computeHProjDistance (const TInputPoint &a, const TInputPoint &b, const TInputPoint &c, bool &isInside)
 
template<typename TInputPoint >
double computeVProjDistance (const TInputPoint &a, const TInputPoint &b, const TInputPoint &c, bool &isInside)
 
template<typename TInputPoint >
double computeEuclideanDistance (const TInputPoint &a, const TInputPoint &b, const TInputPoint &c, bool &isInside)
 
template<typename ForwardIterator , typename OutputIterator , typename Functor >
void melkmanConvexHullAlgorithm (const ForwardIterator &itb, const ForwardIterator &ite, OutputIterator res, Functor &aFunctor)
 Procedure that retrieves the vertices of the hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Melkman algorithm [Melkman, 1979 : [90]]. More...
 

Variables

constexpr double angleTolerance = 1e-6
 

Detailed Description

Hull2D namespace gathers useful functions to compute and return the convex hull or the alpha-shape of a range of 2D points.

Enumeration Type Documentation

◆ ThicknessDefinition

The 2 thickness definitions.

Enumerator
HorizontalVerticalThickness 
EuclideanThickness 

Definition at line 78 of file Hull2DHelpers.h.

Function Documentation

◆ andrewConvexHullAlgorithm()

template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void DGtal::functions::Hull2D::andrewConvexHullAlgorithm ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known monotone-chain algorithm due to [Andrew, 1979 : [6]].

  • first, points are sorted along the horizontal axis.
  • then, the lower and upper convex hull are computed by a simple Graham scan.
    See also
    Hull2D::openGrahamScan
    Postcondition
    The first point of the resulting list of extremal points follows the one with minimal x-coordinate and y-coordinate. Orientation depends on the predicate.
    Parameters
    itbbegin iterator
    iteend iterator
    resoutput iterator used to export the retrieved points
    aPredicateany ternary predicate
    Template Parameters
    ForwardIteratora model of forward and readable iterator
    OutputIteratora model of incrementable and writable iterator
    Predicatea model of ternary predicate

◆ buildHullWithAdaptedStack()

template<typename Stack , typename ForwardIterator , typename Predicate >
void DGtal::functions::Hull2D::buildHullWithAdaptedStack ( Stack  aStack,
const ForwardIterator &  itb,
const ForwardIterator &  ite,
const Predicate &  aPredicate 
)

Procedure that calls Hull2D::buildHullWithStack on a copy of the stack object used to retrieved the hull vertices. Useful when the first argument is a stack adapter returned by a function: it must be copied before being passed by reference in Hull2D::buildHullWithStack.

Parameters
aStackstack copy
itbbegin iterator
iteend iterator
aPredicatepredicate
Template Parameters
Stacka model of CStack
ForwardIteratora model of forward and readable iterator
Predicatea model of ternary predicate
See also
Hull2D::buildHullWithStack

◆ buildHullWithStack()

template<typename Stack , typename ForwardIterator , typename Predicate >
void DGtal::functions::Hull2D::buildHullWithStack ( Stack &  aStack,
const ForwardIterator &  itb,
const ForwardIterator &  ite,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time. This technique is called Sklansky's scan, Graham's scan or 3-coins algorithm.
It works for all WEVP [Toussaint and Avis, 1982 : [118]].

Parameters
aStackreference to the stack of retrieved vertices
itbbegin iterator
iteend iterator
aPredicatepredicate
Template Parameters
Stacka model of CStack
ForwardIteratora model of forward and readable iterator
Predicatea model of ternary predicate
See also
Hull2D::updateHullWithStack

◆ closedGrahamScanFromAnyPoint()

template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void DGtal::functions::Hull2D::closedGrahamScanFromAnyPoint ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time.

See also
Hull2D::buildHullWithStack Hull2D::closedGrahamScanFromVertex

NB: We do not assume that the starting point of the polygon is an extremal point like in Hull2D::closedGrahamScanFromVertex

Parameters
itbbegin iterator
iteend iterator
resoutput iterator used to export the retrieved points
aPredicateany ternary predicate
Template Parameters
ForwardIteratora model of forward and readable iterator
OutputIteratora model of incremental and writable iterator
Predicatea model of ternary predicate

◆ closedGrahamScanFromVertex()

template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void DGtal::functions::Hull2D::closedGrahamScanFromVertex ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time.

See also
Hull2D::buildHullWithStack
Parameters
itbbegin iterator
iteend iterator
resoutput iterator used to export the retrieved points
aPredicateany ternary predicate
Precondition
we assume that the starting point of the polygon is an extremal point.
Template Parameters
ForwardIteratora model of forward and readable iterator
OutputIteratora model of incremental and writable iterator
Predicatea model of ternary predicate

◆ computeEuclideanDistance()

template<typename TInputPoint >
double DGtal::functions::Hull2D::computeEuclideanDistance ( const TInputPoint &  a,
const TInputPoint &  b,
const TInputPoint &  c,
bool &  isInside 
)

Computes the euclidean distance a point c according to the segment [a, b]. (i.e the distance between c and its projected point on [a,b].

Parameters
[in]aone point of the segment.
[in]ba second point of the segment.
[in]cthe point for which the vertical distance is computed.
[out]isInsideindicates if the projected point is inside the segment or not.

◆ computeHProjDistance()

template<typename TInputPoint >
double DGtal::functions::Hull2D::computeHProjDistance ( const TInputPoint &  a,
const TInputPoint &  b,
const TInputPoint &  c,
bool &  isInside 
)

Computes the horizontal distance a point according to the segment [ a , b ]. (i.e the horizontal projection distance of on [ a , b ]).

Note
if the segment [a, b] is horizontal (i.e a [1]==b[1]) then an infinite value (std::numerics<double>::max()) is returned.
Parameters
[in]aone point of the segment.
[in]ba second point of the segment.
[in]cthe point for which the horizontal distance is computed.
[out]isInsideindicates if the projected point is inside the segment or not.

◆ computeHullThickness() [1/2]

template<typename ForwardIterator >
double DGtal::functions::Hull2D::computeHullThickness ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
const ThicknessDefinition def 
)

Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical or Euclidean distances). It takes as input the vertices of the hull given by the range [itbn, ite). The procedure applies the classic rotating caliper to recover all anti-podal pairs.

Typical use:

typedef PointVector<2,DGtal::int32_t> Point;
typedef InHalfPlaneBySimple3x3Matrix<Point, DGtal::int32_t> Functor;
ch.add(Point(0,0));
ch.add(Point(11,1));
ch.add(Point(12,3));
ch.add(Point(8,3));
ch.add(Point(4,5));
ch.add(Point(2,6));
ch.add(Point(1,4));
double th = computeHullThickness(ch.begin(), ch.end(),
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of...
double computeHullThickness(const ForwardIterator &itb, const ForwardIterator &ite, const ThicknessDefinition &def)
Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical ...
MyPointD Point
Definition: testClone2.cpp:383
DGtal::MelkmanConvexHull< Point, Functor > ch
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
Parameters
[in]itbbegin iterator on the convex hull points.
[in]iteend iterator on the convex hull points.
[in]defdefinition of the thickness used in the estimation (i.e HorizontalVerticalThickness or EuclideanThickness)
Note
If the convex hull contains 0, 1 or 2 points the thickness of 0 is returned.
Warning
The convex hull should be oriented in counter clockwise else it will return wrong result.
Examples
geometry/tools/exampleConvexHull2D.cpp.

Referenced by convexHull(), and testConvexHullCompThickness().

◆ computeHullThickness() [2/2]

template<typename ForwardIterator , typename TInputPoint >
double DGtal::functions::Hull2D::computeHullThickness ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
const ThicknessDefinition def,
TInputPoint &  antipodalEdgeP,
TInputPoint &  antipodalEdgeQ,
TInputPoint &  antipodalVertexR 
)

Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical or Euclidean distances). It takes as input the vertices of the hull given by the range [itbn, ite). The procedure applies the classic rotating caliper to recover all anti-podal pairs.

Typical use:

typedef PointVector<2,DGtal::int32_t> Point;
typedef InHalfPlaneBySimple3x3Matrix<Point, DGtal::int32_t> Functor;
ch.add(Point(0,0));
ch.add(Point(11,1));
ch.add(Point(12,3));
ch.add(Point(8,3));
ch.add(Point(4,5));
ch.add(Point(2,6));
ch.add(Point(1,4));
Point p, q, s;
double th = computeHullThickness(ch.begin(), ch.end(),
Parameters
[in]itbbegin iterator on the convex hull points.
[in]iteend iterator on the convex hull points.
[in]defdefinition of the thickness used in the estimation (i.e HorizontalVerticalThickness or EuclideanThickness)
[out]antipodalEdgePone point of the antipodal edge associated to the minimal value of convex hull thickness.
[out]antipodalEdgeQone point of the antipodal edge associated to the minimal value of convex hull thickness.
[out]antipodalVertexRthe vertex of the antipodal pair associated to the minimal value of convex hull thickness.
Note
If the convex hull contains 0, 1 or 2 points the thickness of 0 is returned and the antipodal points are updated with the first points (if they exist).
Warning
The convex hull should be oriented in counter clockwise else it will return wrong result.

◆ computeVProjDistance()

template<typename TInputPoint >
double DGtal::functions::Hull2D::computeVProjDistance ( const TInputPoint &  a,
const TInputPoint &  b,
const TInputPoint &  c,
bool &  isInside 
)

Computes the vertical distance a point according to the segment [a, b]. (i.e the vertical projection distance of on [a,b].

Note
if the segment [a, b] is vertical (i.e a [0]== b [0]) then an infinite value (std::numerics<double>::max()) is returned.
Parameters
[in]aone point of the segment.
[in]ba second point of the segment.
[in]cthe point for which the vertical distance is computed.
[out]isInsideindicates if the projected point is inside the segment or not.

◆ getAngle()

template<typename TInputPoint >
double DGtal::functions::Hull2D::getAngle ( const TInputPoint &  a,
const TInputPoint &  b,
const TInputPoint &  c,
const TInputPoint &  d 
)
inline

Computes the angle between the line (a,b) and (c,d)

Parameters
[in]aone of point defining the first line.
[in]ba second point defining the first line.
[in]ca third point defining the second line.
[in]da third point defining the second line.
Returns
the resulting angle is between 0 and 2M_PI.

◆ getThicknessAntipodalPair()

template<typename TInputPoint >
double DGtal::functions::Hull2D::getThicknessAntipodalPair ( const TInputPoint &  p,
const TInputPoint &  q,
const TInputPoint &  r,
const ThicknessDefinition def 
)

Computes the thickness of an anti podal pair (represented by the segment [ p , q ] and vertex r) according to the given distance def definition.

If the distance definition is HorizontalVerticalThickness, it returns the minimal distance between the vertical/horizontal projection of r on ( p , q ).

If the distance definition is EuclideanThickness, it returns the distance between r and its projection on the line ( p ,q ).

Parameters
[in]pthe first point of the edge anti podal pair.
[in]qthe second point of the edge anti podal pair.
[in]rthe vertex of the anti podal pair.
[in]defdefinition of the thickness used in the estimation (i.e HorizontalVerticalThickness or EuclideanThickness).

Referenced by testConvexHullCompThickness().

◆ grahamConvexHullAlgorithm() [1/2]

template<typename ForwardIterator , typename OutputIterator , typename Predicate , typename PolarComparator >
void DGtal::functions::Hull2D::grahamConvexHullAlgorithm ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the convex hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Graham's algorithm [Graham, 1972 : [63]].

  • choose a pole and sort the points in order of increasing angle about the pole (with a counter-clockwise orientation).
  • scan the sorted list of points and remove some points so that the given predicate returns 'true' for all sets of three consecutive points.
    See also
    Hull2D::closedGrahamScanFromVertex
    Postcondition
    The first point of the resulting list of extremal points is guaranteed to be the one with maximal x-coordinate and y-coordinate.
    Warning
    The predicate must be chosen so that is returns 'true' for counter-clockwise oriented 3-point sets. Otherwise, the procedure only returns the last convex hull edge.
    Parameters
    itbbegin iterator
    iteend iterator
    resoutput iterator used to export the retrieved points
    aPredicateany ternary predicate
    Template Parameters
    ForwardIteratora model of forward and readable iterator
    OutputIteratora model of incremental and writable iterator
    Predicatea model of ternary predicate

◆ grahamConvexHullAlgorithm() [2/2]

template<typename ForwardIterator , typename OutputIterator , typename Predicate , typename PolarComparator >
void DGtal::functions::Hull2D::grahamConvexHullAlgorithm ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate,
PolarComparator &  aPolarComparator 
)

Procedure that retrieves the vertices of the convex hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Graham's algorithm [Graham, 1972 : [63]].

  • choose a pole and sort the points in order of increasing angle about the pole with the given comparator.
  • scan the sorted list of points and remove some points so that the given predicate returns 'true' for all sets of three consecutive points.
    See also
    Hull2D::closedGrahamScanFromVertex
    Postcondition
    The first point of the resulting list of extremal points is guaranteed to be the one with maximal x-coordinate and y-coordinate.
    Warning
    The orientation of the predicate and of the polar comparator should be the same. Otherwise, the procedure only returns the last convex hull edge. For instance, you may use a predicate that returns 'true' for three points counter-clockwise oriented together with PolarPointComparatorBy2x2DetComputer.
    Parameters
    itbbegin iterator
    iteend iterator
    resoutput iterator used to export the retrieved points
    aPredicateany ternary predicate
    aPolarComparatorany polar comparator
    Template Parameters
    ForwardIteratora model of forward and readable iterator
    OutputIteratora model of incremental and writable iterator
    Predicatea model of ternary predicate
    PolarComparatora model of CPolarPointComparator2D.

◆ isCoLinearOpp()

template<typename TInputPoint >
bool DGtal::functions::Hull2D::isCoLinearOpp ( const TInputPoint &  a,
const TInputPoint &  b,
const TInputPoint &  c,
const TInputPoint &  d 
)
inline

Determine if the two vectors (a,b) and (c,d) are co linear from opposite directions

Parameters
[in]aone of point defining the first line.
[in]ba second point defining the first line.
[in]ca third point defining the second line.
[in]da third point defining the second line.
Returns
true if the two vectors are co linear from opposite directions, else it returns false.

◆ melkmanConvexHullAlgorithm()

template<typename ForwardIterator , typename OutputIterator , typename Functor >
void DGtal::functions::Hull2D::melkmanConvexHullAlgorithm ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
Functor aFunctor 
)

Procedure that retrieves the vertices of the hull of a set of 2D points given by the range [ itb , ite ). This procedure follows the well-known Melkman algorithm [Melkman, 1979 : [90]].

See also
MelkmanConvexHull
Parameters
itbbegin iterator
iteend iterator
resoutput iterator used to export the retrieved points
aFunctoraFunctor
Template Parameters
ForwardIteratora model of forward and readable iterator
OutputIteratora model of incrementable and writable iterator
Functora model of COrientationFunctor2

◆ openGrahamScan()

template<typename ForwardIterator , typename OutputIterator , typename Predicate >
void DGtal::functions::Hull2D::openGrahamScan ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
OutputIterator  res,
const Predicate &  aPredicate 
)

Procedure that retrieves the vertices of the convex hull of a weakly externally visible polygon (WEVP) in linear-time.

See also
Hull2D::buildHullWithStack
Parameters
itbbegin iterator
iteend iterator
resoutput iterator used to export the retrieved points
aPredicateany ternary predicate
Template Parameters
ForwardIteratora model of forward and readable iterator
OutputIteratora model of incremental and writable iterator
Predicatea model of ternary predicate

◆ updateHullWithAdaptedStack()

template<typename Stack , typename Point , typename Predicate >
void DGtal::functions::Hull2D::updateHullWithAdaptedStack ( Stack  aStack,
const Point aNewPoint,
const Predicate &  aPredicate 
)

Procedure that calls Hull2D::updateHullWithStack on a copy of the stack object used to retrieved the hull vertices. Useful when the first argument is a stack adapter returned by a function: it must be copied before being passed by reference in Hull2D::updateHullWithStack.

Parameters
aStackstack copy
aNewPointnew point
aPredicatepoint predicate
Template Parameters
Stacka model of CStack
Pointa model of point
Predicatea model of ternary predicate

◆ updateHullWithStack()

template<typename Stack , typename Point , typename Predicate >
void DGtal::functions::Hull2D::updateHullWithStack ( Stack &  aStack,
const Point aNewPoint,
const Predicate &  aPredicate 
)

Procedure that updates the hull when an extra point aNewPoint is considered: while the last three consecutive points, ie. aNewPoint and the last two points of the container are not oriented such that the predicate aPredicate returns 'true',
the last point of the container is removed.

See also
Hull2D::buildHullWithStack
Parameters
aStackreference to the stack of retrieved vertices
aNewPointnew point
aPredicatepoint predicate
Template Parameters
Stacka model of CStack
Pointa model of point
Predicatea model of ternary predicate

Variable Documentation

◆ angleTolerance

constexpr double DGtal::functions::Hull2D::angleTolerance = 1e-6
constexpr

Definition at line 76 of file Hull2DHelpers.h.