DGtal  1.5.beta
exampleConvexHull2D.cpp File Reference
#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"
Include dependency graph for exampleConvexHull2D.cpp:

Go to the source code of this file.

Functions

template<typename ForwardIterator , typename Board >
void drawPolygon (const ForwardIterator &itb, const ForwardIterator &ite, Board &aBoard, bool isClosed=true)
 
void convexHull ()
 
int main (int argc, char **argv)
 

Detailed Description

This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.

Author
Tristan Roussillon (trist.nosp@m.an.r.nosp@m.oussi.nosp@m.llon.nosp@m.@liri.nosp@m.s.cn.nosp@m.rs.fr ) Laboratoire d'InfoRmatique en Image et Systèmes d'information - LIRIS (CNRS, UMR 5205), CNRS, France
Date
2013/12/17

An example file named exampleConvexHull2D.

This file is part of the DGtal library.

Definition in file exampleConvexHull2D.cpp.

Function Documentation

◆ convexHull()

void convexHull ( )

Algorithms that computes the convex hull of a point set

[Hull2D-Namespace]

[Hull2D-Namespace]

[Hull2D-Functor]

[Hull2D-Functor]

[Hull2D-StrictPredicateCCW]

[Hull2D-StrictPredicateCCW]

[Hull2D-AndrewAlgo]

[Hull2D-AndrewAlgo] [Hull2D-Caliper-computeBasic]

[Hull2D-Caliper-computeBasic]

[Hull2D-Caliper-computeAnti]

[Hull2D-Caliper-computeAnti]

[Hull2D-Caliper-display]

[Hull2D-Caliper-display]

[Hull2D-LargePredicateCCW]

[Hull2D-LargePredicateCCW]

[Hull2D-StrictPredicateCW]

[Hull2D-StrictPredicateCW]

[Hull2D-GrahamAlgo]

[Hull2D-GrahamAlgo]

[Hull2D-OnLineMelkmanAlgo]

[Hull2D-OnLineMelkmanAlgo]

[Hull2D-OffLineMelkmanAlgo]

[Hull2D-OffLineMelkmanAlgo]

Examples
geometry/tools/exampleConvexHull2D.cpp.

Definition at line 122 of file exampleConvexHull2D.cpp.

123 {
124  //Digitization of a disk of radius 6
125  Ball2D<Z2i::Space> ball(Z2i::Point(0,0), 6);
126  Z2i::Domain domain(ball.getLowerBound(), ball.getUpperBound());
127  Z2i::DigitalSet pointSet(domain);
128  Shapes<Z2i::Domain>::euclideanShaper(pointSet, ball);
129 
131  using namespace functions::Hull2D;
133 
136  Functor functor;
138 
139  { //convex hull in counter-clockwise order
140  vector<Z2i::Point> res;
141 
144  StrictPredicate predicate( functor );
146  //according to the last two template arguments, neither strictly negative values, nor zeros are accepted:
147  //the predicate returns 'true' only for strictly positive values returned by the underlying functor.
148 
150  andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
155 
157  Z2i::Point antipodalP, antipodalQ, antipodalS;
158  th = DGtal::functions::Hull2D::computeHullThickness(res.begin(), res.end(), DGtal::functions::Hull2D::HorizontalVerticalThickness, antipodalP, antipodalQ, antipodalS);
160 
161 
162  trace.info() <<" ConvexHull HV thickness: " << th << std::endl;
163  //display
164  Board2D board;
165  drawPolygon( res.begin(), res.end(), board );
168  board.drawCircle( antipodalS[0], antipodalS[1], 0.2) ;
170  board.drawCircle(antipodalP[0], antipodalP[1], 0.2);
171  board.drawCircle(antipodalQ[0], antipodalQ[1], 0.2);
172  board.drawLine(antipodalP[0], antipodalP[1], antipodalQ[0], antipodalQ[1]);
174 
175  board.saveSVG( "ConvexHullCCW.svg" );
176 #ifdef WITH_CAIRO
177  board.saveCairo("ConvexHullCCW.png", Board2D::CairoPNG);
178 #endif
179  }
180 
181  { //convex hull in counter-clockwise order with all the points lying on the edges
182  vector<Z2i::Point> res;
183 
186  LargePredicate predicate( functor );
188  //according to the last template argument, zeros are accepted so that
189  //the predicate returns 'true' for all the positive values returned by the underlying functor.
190 
191  //andrew algorithm
192  andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
193 
194  //display
195  Board2D board;
196  drawPolygon( res.begin(), res.end(), board );
197  board.saveSVG( "ConvexHullCCWWithPointsOnEdges.svg" );
198 #ifdef WITH_CAIRO
199  board.saveCairo("ConvexHullCCWWithPointsOnEdges.png", Board2D::CairoPNG);
200 #endif
201 
202  }
203 
204  { //convex hull in clockwise order
205  vector<Z2i::Point> res;
206 
209  StrictPredicate predicate( functor );
211  //according to the last two argument template,
212  //the predicate returns 'true' only for strictly negative values returned by the underlying functor.
213 
214  //andrew algorithm
215  andrewConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
216 
217  //display
218  Board2D board;
219  drawPolygon( res.begin(), res.end(), board );
220  board.saveSVG( "ConvexHullCW.svg" );
221 #ifdef WITH_CAIRO
222  board.saveCairo("ConvexHullCW.png", Board2D::CairoPNG);
223 #endif
224  }
225 
226  { //convex hull in counter-clockwise order
227  vector<Z2i::Point> res;
228 
229  //geometric predicate
231  StrictPredicate predicate( functor );
232 
234  grahamConvexHullAlgorithm( pointSet.begin(), pointSet.end(), back_inserter( res ), predicate );
236 
237  //display
238  Board2D board;
239  drawPolygon( res.begin(), res.end(), board );
240  board.saveSVG( "ConvexHullCCWbis.svg" );
241 #ifdef WITH_CAIRO
242  board.saveCairo("ConvexHullCCWbis.png", Board2D::CairoPNG);
243 #endif
244  }
245 
246  { //convex hull of a simple polygonal line that is not weakly externally visible
247  vector<Z2i::Point> polygonalLine;
248  polygonalLine.push_back(Z2i::Point(0,0));
249  polygonalLine.push_back(Z2i::Point(0,4));
250  polygonalLine.push_back(Z2i::Point(1,4));
251  polygonalLine.push_back(Z2i::Point(1,1));
252  polygonalLine.push_back(Z2i::Point(3,1));
253  polygonalLine.push_back(Z2i::Point(2,2));
254  polygonalLine.push_back(Z2i::Point(3,4));
255  polygonalLine.push_back(Z2i::Point(4,4));
256  polygonalLine.push_back(Z2i::Point(4,0));
257 
258  vector<Z2i::Point> resGraham, res;
259 
261  StrictPredicate predicate( functor );
262  closedGrahamScanFromAnyPoint( polygonalLine.begin(), polygonalLine.end(), back_inserter( resGraham ), predicate );
263 
266  for (std::vector<Z2i::Point>::const_iterator
267  it = polygonalLine.begin(),
268  itEnd = polygonalLine.end();
269  it != itEnd; ++it)
270  ch.add( *it );
272 
274  melkmanConvexHullAlgorithm( polygonalLine.begin(), polygonalLine.end(), back_inserter( res ), functor );
276 
277  //display
278  Board2D board;
279  drawPolygon( polygonalLine.begin(), polygonalLine.end(), board, true );
280  board.saveSVG( "SimplePolygonalLine.svg" );
281 #ifdef WITH_CAIRO
282  board.saveCairo("SimplePolygonalLine.png", Board2D::CairoPNG);
283 #endif
284  board.clear();
285  drawPolygon( resGraham.begin(), resGraham.end(), board );
286  board.saveSVG( "SimplePolygonalLineGraham.svg" );
287 #ifdef WITH_CAIRO
288  board.saveCairo("SimplePolygonalLineGraham.png", Board2D::CairoPNG);
289 #endif
290  board.clear();
291  drawPolygon( res.begin(), res.end(), board );
292  board.saveSVG( "SimplePolygonalLineMelkman.svg" );
293 #ifdef WITH_CAIRO
294  board.saveCairo("SimplePolygonalLineMelkman.png", Board2D::CairoPNG);
295 #endif
296  board.clear();
297  drawPolygon( ch.begin(), ch.end(), board );
298  board.saveSVG( "SimplePolygonalLineMelkman2.svg" );
299 #ifdef WITH_CAIRO
300  board.saveCairo("SimplePolygonalLineMelkman2.png", Board2D::CairoPNG);
301 #endif
302  }
303 
304  { //order of the points for andrew algorithm
305  vector<Z2i::Point> res;
306  std::copy( pointSet.begin(), pointSet.end(), back_inserter( res ) );
307 
308  std::sort( res.begin(), res.end() );
309 
310  //display
311  Board2D board;
312  drawPolygon( res.begin(), res.end(), board, false );
313  board.saveSVG( "AndrewWEVP.svg" );
314 #ifdef WITH_CAIRO
315  board.saveCairo("AndrewWEVP.png", Board2D::CairoPNG);
316 #endif
317  }
318 
319  { //order of the points for graham algorithm
320  vector<Z2i::Point> res;
321  std::copy( pointSet.begin(), pointSet.end(), back_inserter( res ) );
322 
323  //find an extremal point
324  //NB: we choose the point of greatest x-coordinate
325  //so that the sort step (by a polar comparator)
326  //returns a weakly externally visible polygon
327  std::vector<Z2i::Point>::iterator itMax
328  = std::max_element( res.begin(), res.end() );
329 
330  //sort around this point with a polar comparator
332  comparator.setPole( *itMax );
333  std::sort( res.begin(), res.end(), comparator );
334 
335  //display
336  Board2D board;
337  drawPolygon( res.begin(), res.end(), board, false );
338  board.saveSVG( "GrahamWEVP.svg" );
339 #ifdef WITH_CAIRO
340  board.saveCairo("GrahamWEVP.png", Board2D::CairoPNG);
341 #endif
342  }
343 
344 
345 }
Aim: Model of the concept StarShaped represents any circle in the plane.
Definition: Ball2D.h:61
Aim: This class specializes a 'Board' class so as to display DGtal objects more naturally (with <<)....
Definition: Board2D.h:71
static const Color Red
Definition: Color.h:416
static const Color Blue
Definition: Color.h:419
Aim: A wrapper class around a STL associative container for storing sets of digital points within som...
Aim: Class that implements an orientation functor, ie. it provides a way to compute the orientation o...
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of...
Aim: A utility class for constructing different shapes (balls, diamonds, and others).
std::ostream & info()
Aim: Class that implements a binary point predicate, which is able to compare the position of two giv...
Board & setPenColor(const DGtal::Color &color)
Definition: Board.cpp:297
void drawLine(double x1, double y1, double x2, double y2, int depthValue=-1)
Definition: Board.cpp:367
void clear(const DGtal::Color &color=DGtal::Color::None)
Definition: Board.cpp:151
void drawCircle(double x, double y, double radius, int depthValue=-1)
Definition: Board.cpp:450
void saveSVG(const char *filename, PageSize size=Board::BoundingBox, double margin=10.0) const
Definition: Board.cpp:1011
void saveCairo(const char *filename, CairoType type=CairoPNG, PageSize size=Board::BoundingBox, double margin=10.0) const
Definition: Board.cpp:1138
void drawPolygon(const ForwardIterator &itb, const ForwardIterator &ite, Board &aBoard, bool isClosed=true)
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...
Trace trace
Definition: Common.h:153
DGtal::MelkmanConvexHull< Point, Functor > ch
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
Domain domain

References DGtal::DigitalSetByAssociativeContainer< TDomain, TContainer >::begin(), DGtal::Color::Blue, ch, LibBoard::Board::clear(), DGtal::functions::Hull2D::computeHullThickness(), domain, LibBoard::Board::drawCircle(), LibBoard::Board::drawLine(), drawPolygon(), DGtal::DigitalSetByAssociativeContainer< TDomain, TContainer >::end(), DGtal::Ball2D< TSpace >::getLowerBound(), DGtal::Ball2D< TSpace >::getUpperBound(), DGtal::functions::Hull2D::HorizontalVerticalThickness, DGtal::Trace::info(), DGtal::Color::Red, LibBoard::Board::saveCairo(), LibBoard::Board::saveSVG(), LibBoard::Board::setPenColor(), DGtal::functors::PolarPointComparatorBy2x2DetComputer< TPoint, TDetComputer >::setPole(), and DGtal::trace.

Referenced by main().

◆ drawPolygon()

template<typename ForwardIterator , typename Board >
void drawPolygon ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
Board &  aBoard,
bool  isClosed = true 
)

Definition at line 79 of file exampleConvexHull2D.cpp.

81 {
85 
86  ForwardIterator it = itb;
87  if (it != ite)
88  {//for the first point
89  Point p1 = *it;
90  Point p = p1;
91  //display
92  aBoard << SetMode( p.className(), "Grid" )
93  << CustomStyle( p.className()+"/Grid", new CustomPenColor( Color::Red ) )
94  << p
95  << CustomStyle( p.className()+"/Grid", new CustomPenColor( Color::Black ) );
96 
97  Point prev = p;
98  for (++it; it != ite; ++it, prev = p)
99  {//for all other points
100  //conversion
101  p = *it;
102  //display
103  aBoard << p;
104  if (prev != p)
105  aBoard.drawArrow(prev[0], prev[1], p[0], p[1]);
106  }
107 
108  //display the last edge too
109  if (isClosed)
110  {
111  if (prev != p1)
112  aBoard.drawArrow(prev[0], prev[1], p1[0], p1[1]);
113  }
114  }
115 }
Custom style class redefining the pen color. You may use Board2D::Color::None for transparent color.
Definition: Board2D.h:313
Modifier class in a Board2D stream. Useful to choose your own mode for a given class....
Definition: Board2D.h:247
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.
MyPointD Point
Definition: testClone2.cpp:383

Referenced by convexHull().

◆ main()

int main ( int  argc,
char **  argv 
)

Main procedure. Draw the convex hull of a point set coming from the digitization of a disk.

Parameters
argcnumber of arguments
argvarray of arguments

Definition at line 354 of file exampleConvexHull2D.cpp.

355 {
356  trace.beginBlock ( "Example exampleConvexHull2D" );
357  trace.info() << "Args:";
358  for ( int i = 0; i < argc; ++i )
359  trace.info() << " " << argv[ i ];
360  trace.info() << endl;
361 
362  convexHull();
363 
364  trace.endBlock();
365 
366  return 0;
367 }
void beginBlock(const std::string &keyword="")
double endBlock()
void convexHull()

References DGtal::Trace::beginBlock(), convexHull(), DGtal::Trace::endBlock(), DGtal::Trace::info(), and DGtal::trace.