DGtal  1.5.beta
geometry/tools/exampleConvexHull2D.cpp

Computation of the convex hull of a planar point set by different algorithms:

See also
Convex hull and alpha-shape in the plane
#include <iostream>
#include "DGtal/base/Common.h"
#include "DGtal/base/IteratorCirculatorTraits.h"
#include "DGtal/helpers/StdDefs.h"
#include "DGtal/geometry/tools/Hull2DHelpers.h"
#include "DGtal/geometry/tools/PolarPointComparatorBy2x2DetComputer.h"
#include "DGtal/geometry/tools/determinant/AvnaimEtAl2x2DetSignComputer.h"
#include "DGtal/geometry/tools/determinant/InHalfPlaneBySimple3x3Matrix.h"
#include "DGtal/shapes/ShapeFactory.h"
#include "DGtal/shapes/Shapes.h"
#include "DGtal/topology/DigitalSetBoundary.h"
#include "DGtal/topology/DigitalSurface.h"
#include "DGtal/graph/DepthFirstVisitor.h"
#include "DGtal/io/boards/Board2D.h"
using namespace std;
using namespace DGtal;
/*
* @brief Procedure that draws on a board a polygon given by a range of vertices.
* @param itb begin iterator
* @param ite end iterator
* @param aBoard any board
* @param isClosed bool equal to 'true' for closed polygon, 'false' otherwise.
* @tparam ForwardIterator at least a readable and forward iterator
* @tparam Board equivalent to Board2D
*/
template <typename ForwardIterator, typename Board>
void drawPolygon(const ForwardIterator& itb, const ForwardIterator& ite,
Board& aBoard, bool isClosed = true)
{
ForwardIterator it = itb;
if (it != ite)
{//for the first point
Point p1 = *it;
Point p = p1;
//display
aBoard << SetMode( p.className(), "Grid" )
<< CustomStyle( p.className()+"/Grid", new CustomPenColor( Color::Red ) )
<< p
<< CustomStyle( p.className()+"/Grid", new CustomPenColor( Color::Black ) );
Point prev = p;
for (++it; it != ite; ++it, prev = p)
{//for all other points
//conversion
p = *it;
//display
aBoard << p;
if (prev != p)
aBoard.drawArrow(prev[0], prev[1], p[0], p[1]);
}
//display the last edge too
if (isClosed)
{
if (prev != p1)
aBoard.drawArrow(prev[0], prev[1], p1[0], p1[1]);
}
}
}
void convexHull()
{
//Digitization of a disk of radius 6
Ball2D<Z2i::Space> ball(Z2i::Point(0,0), 6);
Z2i::Domain domain(ball.getLowerBound(), ball.getUpperBound());
Z2i::DigitalSet pointSet(domain);
Shapes<Z2i::Domain>::euclideanShaper(pointSet, ball);
using namespace functions::Hull2D;
typedef InHalfPlaneBySimple3x3Matrix<Z2i::Point, DGtal::int64_t> Functor;
Functor functor;
{ //convex hull in counter-clockwise order
vector<Z2i::Point> res;
typedef PredicateFromOrientationFunctor2<Functor, false, false> StrictPredicate;
StrictPredicate predicate( functor );
//according to the last two template arguments, neither strictly negative values, nor zeros are accepted:
//the predicate returns 'true' only for strictly positive values returned by the underlying functor.
andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
Z2i::Point antipodalP, antipodalQ, antipodalS;
th = DGtal::functions::Hull2D::computeHullThickness(res.begin(), res.end(), DGtal::functions::Hull2D::HorizontalVerticalThickness, antipodalP, antipodalQ, antipodalS);
trace.info() <<" ConvexHull HV thickness: " << th << std::endl;
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board );
board.setPenColor(DGtal::Color::Red);
board.drawCircle( antipodalS[0], antipodalS[1], 0.2) ;
board.setPenColor(DGtal::Color::Blue);
board.drawCircle(antipodalP[0], antipodalP[1], 0.2);
board.drawCircle(antipodalQ[0], antipodalQ[1], 0.2);
board.drawLine(antipodalP[0], antipodalP[1], antipodalQ[0], antipodalQ[1]);
board.saveSVG( "ConvexHullCCW.svg" );
#ifdef WITH_CAIRO
board.saveCairo("ConvexHullCCW.png", Board2D::CairoPNG);
#endif
}
{ //convex hull in counter-clockwise order with all the points lying on the edges
vector<Z2i::Point> res;
typedef PredicateFromOrientationFunctor2<Functor, false, true> LargePredicate;
LargePredicate predicate( functor );
//according to the last template argument, zeros are accepted so that
//the predicate returns 'true' for all the positive values returned by the underlying functor.
//andrew algorithm
andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board );
board.saveSVG( "ConvexHullCCWWithPointsOnEdges.svg" );
#ifdef WITH_CAIRO
board.saveCairo("ConvexHullCCWWithPointsOnEdges.png", Board2D::CairoPNG);
#endif
}
{ //convex hull in clockwise order
vector<Z2i::Point> res;
typedef PredicateFromOrientationFunctor2<Functor, true, false> StrictPredicate;
StrictPredicate predicate( functor );
//according to the last two argument template,
//the predicate returns 'true' only for strictly negative values returned by the underlying functor.
//andrew algorithm
andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board );
board.saveSVG( "ConvexHullCW.svg" );
#ifdef WITH_CAIRO
board.saveCairo("ConvexHullCW.png", Board2D::CairoPNG);
#endif
}
{ //convex hull in counter-clockwise order
vector<Z2i::Point> res;
//geometric predicate
typedef PredicateFromOrientationFunctor2<Functor, false, false> StrictPredicate;
StrictPredicate predicate( functor );
grahamConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board );
board.saveSVG( "ConvexHullCCWbis.svg" );
#ifdef WITH_CAIRO
board.saveCairo("ConvexHullCCWbis.png", Board2D::CairoPNG);
#endif
}
{ //convex hull of a simple polygonal line that is not weakly externally visible
vector<Z2i::Point> polygonalLine;
polygonalLine.push_back(Z2i::Point(0,0));
polygonalLine.push_back(Z2i::Point(0,4));
polygonalLine.push_back(Z2i::Point(1,4));
polygonalLine.push_back(Z2i::Point(1,1));
polygonalLine.push_back(Z2i::Point(3,1));
polygonalLine.push_back(Z2i::Point(2,2));
polygonalLine.push_back(Z2i::Point(3,4));
polygonalLine.push_back(Z2i::Point(4,4));
polygonalLine.push_back(Z2i::Point(4,0));
vector<Z2i::Point> resGraham, res;
typedef PredicateFromOrientationFunctor2<Functor, false, false> StrictPredicate;
StrictPredicate predicate( functor );
closedGrahamScanFromAnyPoint( polygonalLine.begin(), polygonalLine.end(), back_inserter( resGraham ), predicate );
for (std::vector<Z2i::Point>::const_iterator
it = polygonalLine.begin(),
itEnd = polygonalLine.end();
it != itEnd; ++it)
ch.add( *it );
melkmanConvexHullAlgorithm( polygonalLine.begin(), polygonalLine.end(), back_inserter( res ), functor );
//display
Board2D board;
drawPolygon( polygonalLine.begin(), polygonalLine.end(), board, true );
board.saveSVG( "SimplePolygonalLine.svg" );
#ifdef WITH_CAIRO
board.saveCairo("SimplePolygonalLine.png", Board2D::CairoPNG);
#endif
board.clear();
drawPolygon( resGraham.begin(), resGraham.end(), board );
board.saveSVG( "SimplePolygonalLineGraham.svg" );
#ifdef WITH_CAIRO
board.saveCairo("SimplePolygonalLineGraham.png", Board2D::CairoPNG);
#endif
board.clear();
drawPolygon( res.begin(), res.end(), board );
board.saveSVG( "SimplePolygonalLineMelkman.svg" );
#ifdef WITH_CAIRO
board.saveCairo("SimplePolygonalLineMelkman.png", Board2D::CairoPNG);
#endif
board.clear();
drawPolygon( ch.begin(), ch.end(), board );
board.saveSVG( "SimplePolygonalLineMelkman2.svg" );
#ifdef WITH_CAIRO
board.saveCairo("SimplePolygonalLineMelkman2.png", Board2D::CairoPNG);
#endif
}
{ //order of the points for andrew algorithm
vector<Z2i::Point> res;
std::copy( pointSet.begin(), pointSet.end(), back_inserter( res ) );
std::sort( res.begin(), res.end() );
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board, false );
board.saveSVG( "AndrewWEVP.svg" );
#ifdef WITH_CAIRO
board.saveCairo("AndrewWEVP.png", Board2D::CairoPNG);
#endif
}
{ //order of the points for graham algorithm
vector<Z2i::Point> res;
std::copy( pointSet.begin(), pointSet.end(), back_inserter( res ) );
//find an extremal point
//NB: we choose the point of greatest x-coordinate
//so that the sort step (by a polar comparator)
//returns a weakly externally visible polygon
std::vector<Z2i::Point>::iterator itMax
= std::max_element( res.begin(), res.end() );
//sort around this point with a polar comparator
functors::PolarPointComparatorBy2x2DetComputer<Z2i::Point> comparator;
comparator.setPole( *itMax );
std::sort( res.begin(), res.end(), comparator );
//display
Board2D board;
drawPolygon( res.begin(), res.end(), board, false );
board.saveSVG( "GrahamWEVP.svg" );
#ifdef WITH_CAIRO
board.saveCairo("GrahamWEVP.png", Board2D::CairoPNG);
#endif
}
}
int main( int argc, char** argv )
{
trace.beginBlock ( "Example exampleConvexHull2D" );
trace.info() << "Args:";
for ( int i = 0; i < argc; ++i )
trace.info() << " " << argv[ i ];
trace.info() << endl;
return 0;
}
// //
static const Color Red
Definition: Color.h:416
static const Color Blue
Definition: Color.h:419
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of...
void beginBlock(const std::string &keyword="")
std::ostream & info()
double endBlock()
void drawPolygon(const ForwardIterator &itb, const ForwardIterator &ite, Board &aBoard, bool isClosed=true)
void convexHull()
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 ,...
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...
double computeHullThickness(const ForwardIterator &itb, const ForwardIterator &ite, const ThicknessDefinition &def)
Procedure to compute the convex hull thickness given from different definitions (Horizontal/vertical ...
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 ,...
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 [ i...
DGtal is the top-level namespace which contains all DGtal functions and types.
Trace trace
Definition: Common.h:153
Go to http://www.boost.org/doc/libs/1_52_0/libs/iterator/doc/ForwardTraversal.html.
Go to http://www.boost.org/doc/libs/1_52_0/libs/iterator/doc/ReadableIterator.html.
int main(int argc, char **argv)
MyPointD Point
Definition: testClone2.cpp:383
DGtal::MelkmanConvexHull< Point, Functor > ch
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
Domain domain
HyperRectDomain< Space > Domain
Z2i::DigitalSet DigitalSet