DGtal  1.4.beta
MelkmanConvexHull.h
1 
17 #pragma once
18 
31 #if defined(MelkmanConvexHull_RECURSES)
32 #error Recursive header files inclusion detected in MelkmanConvexHull.h
33 #else // defined(MelkmanConvexHull_RECURSES)
34 
35 #define MelkmanConvexHull_RECURSES
36 
37 #if !defined MelkmanConvexHull_h
38 
39 #define MelkmanConvexHull_h
40 
42 // Inclusions
43 #include <iostream>
44 #include "DGtal/base/Common.h"
45 #include "DGtal/base/Alias.h"
46 #include "DGtal/base/CountedConstPtrOrConstPtr.h"
47 #include "DGtal/base/IteratorCirculatorTraits.h"
48 #include "DGtal/base/FrontInsertionSequenceToStackAdapter.h"
49 #include "DGtal/base/BackInsertionSequenceToStackAdapter.h"
50 
51 #include "DGtal/geometry/tools/determinant/COrientationFunctor2.h"
52 #include "DGtal/geometry/tools/determinant/PredicateFromOrientationFunctor2.h"
53 
54 #include "DGtal/geometry/tools/Hull2DHelpers.h"
56 
57 namespace DGtal
58 {
59 
61  // template class MelkmanConvexHull
87  template <typename TPoint,
88  typename TOrientationFunctor >
90  {
91  // ----------------------- Types ------------------------------------------
92  public:
93 
98 
102  typedef TPoint Point;
106  typedef TOrientationFunctor Functor;
108  //the two types of points must be the same
109  BOOST_STATIC_ASSERT (( boost::is_same< Point, typename Functor::Point >::value ));
110 
119 
123  typedef typename std::deque<Point>::const_iterator ConstIterator;
124 
125  // ----------------------- Standard services ------------------------------
126  public:
127 
128  MelkmanConvexHull( Alias<Functor> aFunctor);
130 
136  MelkmanConvexHull( const MelkmanConvexHull & mch ) = default;
137 
138  // ----------------------- Interface --------------------------------------
139  public:
140 
147  void add ( const Point& aPoint );
148 
157  ConstIterator begin() const;
158 
163  ConstIterator end() const;
164 
169  void selfDisplay ( std::ostream & out ) const;
170 
175  bool isValid() const;
176 
183  Self & operator= ( const Self & mch );
184 
189  const Point & operator[](unsigned int i) const;
190 
194  unsigned int size() const;
195 
199  void clear();
200 
207  void reverse();
208 
209 
210  // ------------------------- Private Datas --------------------------------
211  private:
216  std::deque<Point> myContainer;
233 
234  // ------------------------- Internals ------------------------------------
235  private:
236 
237  }; // end of class MelkmanConvexHull
238 
245  template <typename TPoint, typename TOrientationFunctor>
246  std::ostream&
247  operator<< ( std::ostream & out, const MelkmanConvexHull<TPoint, TOrientationFunctor> & object );
248 
249  namespace functions
250  {
251  namespace Hull2D
252  {
271  template <typename ForwardIterator,
272  typename OutputIterator,
273  typename Functor >
274  void melkmanConvexHullAlgorithm(const ForwardIterator& itb,
275  const ForwardIterator& ite,
276  OutputIterator res,
277  Functor& aFunctor );
278  } //namespace Hull2D
279  } //namespace functions
280 
281 } // namespace DGtal
282 
283 
285 // Includes inline functions.
286 #include "DGtal/geometry/tools/MelkmanConvexHull.ih"
287 
288 // //
290 
291 #endif // !defined MelkmanConvexHull_h
292 
293 #undef MelkmanConvexHull_RECURSES
294 #endif // else defined(MelkmanConvexHull_RECURSES)
DGtal::MelkmanConvexHull::myContainer
std::deque< Point > myContainer
Definition: MelkmanConvexHull.h:216
DGtal::MelkmanConvexHull::reverse
void reverse()
DGtal::MelkmanConvexHull::Functor
TOrientationFunctor Functor
Definition: MelkmanConvexHull.h:106
DGtal::MelkmanConvexHull::operator[]
const Point & operator[](unsigned int i) const
aPoint
const Point aPoint(3, 4)
DGtal::functions::Hull2D::melkmanConvexHullAlgorithm
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 ,...
Functor
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
Definition: testConvexHull2DReverse.cpp:51
DGtal::MelkmanConvexHull::Point
TPoint Point
Definition: MelkmanConvexHull.h:102
DGtal::MelkmanConvexHull::myBackwardPredicate
BackwardPredicate myBackwardPredicate
Definition: MelkmanConvexHull.h:220
DGtal::MelkmanConvexHull::isValid
bool isValid() const
DGtal::MelkmanConvexHull::size
unsigned int size() const
DGtal::operator<<
std::ostream & operator<<(std::ostream &out, const ATu0v1< TKSpace, TLinearAlgebra > &object)
DGtal::MelkmanConvexHull::BOOST_STATIC_ASSERT
BOOST_STATIC_ASSERT((boost::is_same< Point, typename Functor::Point >::value))
DGtal::MelkmanConvexHull::ForwardPredicate
PredicateFromOrientationFunctor2< Functor, true, false > ForwardPredicate
Definition: MelkmanConvexHull.h:118
DGtal::MelkmanConvexHull::myDefaultFunctor
Functor myDefaultFunctor
Definition: MelkmanConvexHull.h:228
DGtal::MelkmanConvexHull::begin
ConstIterator begin() const
DGtal::MelkmanConvexHull::Self
MelkmanConvexHull< TPoint, TOrientationFunctor > Self
Definition: MelkmanConvexHull.h:97
DGtal::MelkmanConvexHull::selfDisplay
void selfDisplay(std::ostream &out) const
DGtal
DGtal is the top-level namespace which contains all DGtal functions and types.
DGtal::MelkmanConvexHull::end
ConstIterator end() const
DGtal::concepts::COrientationFunctor2
Aim: This concept is a refinement of COrientationFunctor, useful for simple algebraic curves that can...
Definition: COrientationFunctor2.h:84
DGtal::InHalfPlaneBySimple3x3Matrix
Aim: Class that implements an orientation functor, ie. it provides a way to compute the orientation o...
Definition: InHalfPlaneBySimple3x3Matrix.h:91
DGtal::MelkmanConvexHull::operator=
Self & operator=(const Self &mch)
DGtal::MelkmanConvexHull::ConstIterator
std::deque< Point >::const_iterator ConstIterator
Definition: MelkmanConvexHull.h:123
DGtal::MelkmanConvexHull::BOOST_CONCEPT_ASSERT
BOOST_CONCEPT_ASSERT((concepts::COrientationFunctor2< Functor >))
DGtal::MelkmanConvexHull::BackwardPredicate
PredicateFromOrientationFunctor2< Functor, false, false > BackwardPredicate
Definition: MelkmanConvexHull.h:114
DGtal::MelkmanConvexHull
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of...
Definition: MelkmanConvexHull.h:89
DGtal::PredicateFromOrientationFunctor2< Functor, false, false >
DGtal::Alias
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
Definition: Alias.h:182
DGtal::MelkmanConvexHull::add
void add(const Point &aPoint)
DGtal::MelkmanConvexHull::myForwardPredicate
ForwardPredicate myForwardPredicate
Definition: MelkmanConvexHull.h:224
DGtal::MelkmanConvexHull::clear
void clear()
DGtal::MelkmanConvexHull::MelkmanConvexHull
MelkmanConvexHull()
DGtal::MelkmanConvexHull::myFirstPoint
Point myFirstPoint
Definition: MelkmanConvexHull.h:232