DGtal  1.2.0
Functions
exampleAlphaShape.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/geometry/tools/determinant/InGeneralizedDiskOfGivenRadius.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 exampleAlphaShape.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 alphaShape ()
 
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 exampleAlphaShape.

This file is part of the DGtal library.

Definition in file exampleAlphaShape.cpp.

Function Documentation

◆ alphaShape()

void alphaShape ( )

Algorithms that computes the alpha-shape of a point set

[Hull2D-RadiusPredicateInf]

[Hull2D-RadiusPredicateInf]

[Hull2D-ClosedGrahamScan]

[Hull2D-ClosedGrahamScan]

[Hull2D-RadiusPredicateM1]

[Hull2D-RadiusPredicateM1]

[Hull2D-RadiusPredicateMsqrt5]

[Hull2D-RadiusPredicateMsqrt5]

[Hull2D-RadiusPredicateM5]

[Hull2D-RadiusPredicateM5]

[Hull2D-RadiusPredicateP8]

[Hull2D-RadiusPredicateP8]

[Hull2D-RadiusPredicateP9]

[Hull2D-RadiusPredicateP9]

Examples
geometry/tools/exampleAlphaShape.cpp.

Definition at line 123 of file exampleAlphaShape.cpp.

124 {
125 
126  //Digitization of a disk of radius 8
127  Ball2D<Z2i::Space> ball(Z2i::Point(0,0), 8);
128  Z2i::Domain domain(ball.getLowerBound(), ball.getUpperBound());
129  Z2i::DigitalSet digitalSet(domain);
130  Shapes<Z2i::Domain>::euclideanShaper(digitalSet, ball);
131 
132  //Contour of the digital set
133  Z2i::KSpace kspace;
134  kspace.init(domain.lowerBound()-Z2i::Point(1,1), domain.upperBound()+Z2i::Point(1,1), true);
135  typedef DigitalSetBoundary<Z2i::KSpace, Z2i::DigitalSet> DigitalSurfaceContainer;
136  typedef DigitalSurface<DigitalSurfaceContainer> CustomDigitalSurface;
137  DigitalSurfaceContainer digitalSurfaceContainer( kspace, digitalSet );
138  CustomDigitalSurface digitalSurface( digitalSurfaceContainer );
139 
140  //Grid curve
141  Z2i::Curve gridCurve;
143  CustomVisitor visitor( digitalSurface, *digitalSurface.begin() );
144  while ( ! visitor.finished() )
145  {
146  gridCurve.pushBack( visitor.current().first );
147  visitor.expand();
148  }
149 
150  //Point set defined as the set of (not duplicated) inner points
151  typedef Z2i::Curve::InnerPointsRange PointRange;
152  PointRange pointsRange = gridCurve.getInnerPointsRange();
153  vector<Z2i::Point> border;
154  unique_copy( pointsRange.begin(), pointsRange.end(), back_inserter( border ) );
155 
156  //namespace for hull functions
157  using namespace functions::Hull2D;
158 
159  { //alpha = 0
160  trace.info() << " alpha == 0 " << endl;
161  vector<Z2i::Point> res;
162 
166  Functor functor(true, 1, 0); //alpha = 0; 1/alpha -> +inf
168  Predicate predicate( functor );
170 
172  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
174 
175  //display
176  Board2D board;
177  drawPolygon( res.begin(), res.end(), board );
178  board.saveSVG( "AlphaShape0.svg" );
179 #ifdef WITH_CAIRO
180  board.saveCairo("AlphaShape0.png", Board2D::CairoPNG);
181 #endif
182 
183  }
184 
185  { //alpha = 0
186  trace.info() << " alpha == 0 " << endl;
187  vector<Z2i::Point> res;
188 
189  //comparator and functor
191  Functor functor;
193  Predicate predicate( functor );
194 
195  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
196 
197  //display
198  Board2D board;
199  drawPolygon( res.begin(), res.end(), board );
200  board.saveSVG( "AlphaShape0bis.svg" );
201 #ifdef WITH_CAIRO
202  board.saveCairo("AlphaShape0bis.png", Board2D::CairoPNG);
203 #endif
204 
205  }
206 
207  //negative alpha shape
208  { //alpha = -1
209  trace.info() << " alpha == -1 " << endl;
210  vector<Z2i::Point> res;
211 
215  Functor functor(false, 1, 1); //1/alpha = -sqrt(1/1) = -1
218  Predicate predicate( functor );
219 
220  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
221 
222  //display
223  Board2D board;
224  drawPolygon( res.begin(), res.end(), board );
225  board.saveSVG( "AlphaShapeM1.svg" );
226 #ifdef WITH_CAIRO
227  board.saveCairo("AlphaShapeM1.png", Board2D::CairoPNG);
228 #endif
229 
230  }
231 
232  { //alpha = -sqrt(5)
233  trace.info() << " alpha == -sqrt(5) " << endl;
234  vector<Z2i::Point> res;
235 
239  Functor functor(false, 5, 1); //1/alpha = -sqrt(5)
242  Predicate predicate( functor );
243 
244  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
245 
246  //display
247  Board2D board;
248  drawPolygon( res.begin(), res.end(), board );
249  board.saveSVG( "AlphaShapeMSqrt5.svg" );
250 #ifdef WITH_CAIRO
251  board.saveCairo("AlphaShapeMSqrt5.png", Board2D::CairoPNG);
252 #endif
253 
254  }
255 
256  { //alpha = -5
257  trace.info() << " alpha == -5 " << endl;
258  vector<Z2i::Point> res;
259 
263  Functor functor(false, 25, 1); //1/alpha = -sqrt(25/1) = -5
266  Predicate predicate( functor );
267 
268  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
269 
270  //display
271  Board2D board;
272  drawPolygon( res.begin(), res.end(), board );
273  board.saveSVG( "AlphaShapeM5.svg" );
274 #ifdef WITH_CAIRO
275  board.saveCairo("AlphaShapeM5.png", Board2D::CairoPNG);
276 #endif
277 
278  }
279 
280  //positive alpha shape
281  {
282  trace.info() << " alpha == 8 " << endl;
283  vector<Z2i::Point> res;
284 
288  Functor functor(true, 64, 1); //1/alpha = sqrt(64/1) = 8
290  Predicate predicate( functor );
292 
293  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
294 
295  //display
296  Board2D board;
297  drawPolygon( res.begin(), res.end(), board );
298  board.saveSVG( "AlphaShapeP8.svg" );
299 #ifdef WITH_CAIRO
300  board.saveCairo("AlphaShapeP8.png", Board2D::CairoPNG);
301 #endif
302 
303  }
304 
305  //positive alpha shape
306  {
307  trace.info() << " alpha == 9 " << endl;
308  vector<Z2i::Point> res;
309 
313  Functor functor(true, 81, 1); //1/alpha = sqrt(81/1) = 9
316  Predicate predicate( functor );
317 
318  closedGrahamScanFromAnyPoint( border.begin(), border.end(), back_inserter( res ), predicate );
319 
320  //display
321  Board2D board;
322  drawPolygon( res.begin(), res.end(), board );
323  board.saveSVG( "AlphaShapeP9.svg" );
324 #ifdef WITH_CAIRO
325  board.saveCairo("AlphaShapeP9.png", Board2D::CairoPNG);
326 #endif
327 
328  }
329 
330 }
Aim: Class that provides a way of computing the sign of the determinant of a 2x2 matrix from its four...
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
Aim: model of CConstBidirectionalRange that adapts any range of elements bounded by two iterators [it...
Aim: This class is useful to perform a depth-first exploration of a graph given a starting point or s...
Aim: A model of CDigitalSurfaceContainer which defines the digital surface as the boundary of a given...
Aim: A wrapper class around a STL associative container for storing sets of digital points within som...
Aim: Represents a set of n-1-cells in a nD space, together with adjacency relation between these cell...
Aim: describes, in a cellular space of dimension n, a closed or open sequence of signed d-cells (or d...
Definition: GridCurve.h:173
void pushBack(const SCell &aSCell)
InnerPointsRange getInnerPointsRange() const
Definition: GridCurve.h:462
const Point & lowerBound() const
const Point & upperBound() const
Aim: This class implements an orientation functor that provides a way to determine the position of ...
Aim: Class that implements an orientation functor, ie. it provides a way to compute the orientation o...
Aim: This class is a model of CCellularGridSpaceND. It represents the cubical grid as a cell complex,...
bool init(const Point &lower, const Point &upper, bool isClosed)
Specifies the upper and lower bounds for the maximal cells in this space.
Aim: Small adapter to models of COrientationFunctor2. It is a model of concepts::CPointPredicate....
Aim: A utility class for constructing different shapes (balls, diamonds, and others).
std::ostream & info()
void saveSVG(const char *filename, PageSize size=Board::BoundingBox, double margin=10.0) const
Definition: Board.cpp:1012
void saveCairo(const char *filename, CairoType type=CairoPNG, PageSize size=Board::BoundingBox, double margin=10.0) const
Definition: Board.cpp:1139
void drawPolygon(const ForwardIterator &itb, const ForwardIterator &ite, Board &aBoard, bool isClosed=true)
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...
Trace trace
Definition: Common.h:154
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
Domain domain

References domain, drawPolygon(), DGtal::GridCurve< TKSpace >::getInnerPointsRange(), DGtal::Ball2D< TSpace >::getLowerBound(), DGtal::Ball2D< TSpace >::getUpperBound(), DGtal::Trace::info(), DGtal::KhalimskySpaceND< dim, TInteger >::init(), DGtal::HyperRectDomain< TSpace >::lowerBound(), DGtal::GridCurve< TKSpace >::pushBack(), LibBoard::Board::saveCairo(), LibBoard::Board::saveSVG(), DGtal::trace, and DGtal::HyperRectDomain< TSpace >::upperBound().

Referenced by main().

◆ drawPolygon()

template<typename ForwardIterator , typename Board >
void drawPolygon ( const ForwardIterator &  itb,
const ForwardIterator &  ite,
Board &  aBoard,
bool  isClosed = true 
)
Examples
geometry/tools/exampleAlphaShape.cpp, and geometry/tools/exampleConvexHull2D.cpp.

Definition at line 79 of file exampleAlphaShape.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 alphaShape().

◆ main()

int main ( int  argc,
char **  argv 
)

Main procedure. Draw the convex hull and alpha-shape (for several values of alpha) of a point set coming from the digitization of a disk.

Parameters
argcnumber of arguments
argvarray of arguments

Definition at line 339 of file exampleAlphaShape.cpp.

340 {
341  trace.beginBlock ( "Example exampleAlphaShape" );
342  trace.info() << "Args:";
343  for ( int i = 0; i < argc; ++i )
344  trace.info() << " " << argv[ i ];
345  trace.info() << endl;
346 
347  alphaShape();
348 
349  trace.endBlock();
350 
351  return 0;
352 }
void beginBlock(const std::string &keyword="")
double endBlock()
void alphaShape()

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