DGtal 1.3.0
|
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of a simple polygonal line (without self-intersection) [Melkman, 1987: [85]]. More...
#include <DGtal/geometry/tools/MelkmanConvexHull.h>
Public Types | |
typedef MelkmanConvexHull< TPoint, TOrientationFunctor > | Self |
typedef TPoint | Point |
typedef TOrientationFunctor | Functor |
typedef PredicateFromOrientationFunctor2< Functor, false, false > | BackwardPredicate |
typedef PredicateFromOrientationFunctor2< Functor, true, false > | ForwardPredicate |
typedef std::deque< Point >::const_iterator | ConstIterator |
Public Member Functions | |
BOOST_CONCEPT_ASSERT ((concepts::COrientationFunctor2< Functor >)) | |
BOOST_STATIC_ASSERT ((boost::is_same< Point, typename Functor::Point >::value)) | |
MelkmanConvexHull (Alias< Functor > aFunctor) | |
MelkmanConvexHull () | |
MelkmanConvexHull (const MelkmanConvexHull &mch)=default | |
void | add (const Point &aPoint) |
ConstIterator | begin () const |
ConstIterator | end () const |
void | selfDisplay (std::ostream &out) const |
bool | isValid () const |
Self & | operator= (const Self &mch) |
const Point & | operator[] (unsigned int i) const |
unsigned int | size () const |
void | clear () |
void | reverse () |
Private Attributes | |
std::deque< Point > | myContainer |
BackwardPredicate | myBackwardPredicate |
ForwardPredicate | myForwardPredicate |
Functor | myDefaultFunctor |
Point | myFirstPoint |
Aim: This class implements the on-line algorithm of Melkman for the computation of the convex hull of a simple polygonal line (without self-intersection) [Melkman, 1987: [85]].
This algorithm is based on a deque, which stores the vertices of the convex hull (for convenience, the first and last vertex contained in the deque are the same point). Since we assume that the input points form a simple polygonal line, a new point cannot be located in the cone formed by the first and last edges of the current convex hull. As a consequence, it is enough to update the convex hull by a Graham scan from the front and/or from the back of the deque; it is never required to remove/insert points in the middle of the container. See Convex hull and alpha-shape in the plane for more details.
Note that if the input points do not form a simple polygonal line, the behavior is not defined.
TPoint | a model of point |
TOrientationFunctor | a model of COrientationFunctor2 (whose inner type 'Point' match to 'TPoint') |
Definition at line 89 of file MelkmanConvexHull.h.
typedef PredicateFromOrientationFunctor2<Functor,false,false> DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::BackwardPredicate |
Type of predicate devoted to the backward scan
Definition at line 114 of file MelkmanConvexHull.h.
typedef std::deque<Point>::const_iterator DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::ConstIterator |
Type of iterator on the convex hull vertices
Definition at line 123 of file MelkmanConvexHull.h.
typedef PredicateFromOrientationFunctor2<Functor,true,false> DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::ForwardPredicate |
Type of predicate devoted to the forward scan
Definition at line 118 of file MelkmanConvexHull.h.
typedef TOrientationFunctor DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::Functor |
Type of orientation functor
Definition at line 106 of file MelkmanConvexHull.h.
typedef TPoint DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::Point |
Type of point
Definition at line 102 of file MelkmanConvexHull.h.
typedef MelkmanConvexHull<TPoint, TOrientationFunctor> DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::Self |
Self type
Definition at line 97 of file MelkmanConvexHull.h.
DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::MelkmanConvexHull | ( | Alias< Functor > | aFunctor | ) |
DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::MelkmanConvexHull | ( | ) |
|
default |
Default copy constructor.
mch | the object to copy. |
void DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::add | ( | const Point & | aPoint | ) |
Consider a new point and possibly update the convex hull.
aPoint | an extra point |
ConstIterator DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::begin | ( | ) | const |
Begin iterator
DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::BOOST_CONCEPT_ASSERT | ( | (concepts::COrientationFunctor2< Functor >) | ) |
DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::BOOST_STATIC_ASSERT | ( | (boost::is_same< Point, typename Functor::Point >::value) | ) |
void DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::clear | ( | ) |
clear the current content of the convex hull.
ConstIterator DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::end | ( | ) | const |
End iterator
bool DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::isValid | ( | ) | const |
Checks the validity/consistency of the object.
Self & DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::operator= | ( | const Self & | mch | ) |
Assignement Operator
mch | the object to copy. |
const Point & DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::operator[] | ( | unsigned int | i | ) | const |
i | the index of the considered point. |
void DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::reverse | ( | ) |
Reverse the convex hull container allowing to change the order of adding points from the front or from the back in reference to an input contour. Such a reverse is important to avoid wrong convexhull and to allow convex hull extension from two directions.
void DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::selfDisplay | ( | std::ostream & | out | ) | const |
Writes/Displays the object on an output stream.
out | the output stream where the object is written. |
unsigned int DGtal::MelkmanConvexHull< TPoint, TOrientationFunctor >::size | ( | ) | const |
|
private |
Predicate devoted to the backward scan
Definition at line 220 of file MelkmanConvexHull.h.
|
private |
Deque container, which stores the vertices of the convex hull NB: the first and last point is the same.
Definition at line 216 of file MelkmanConvexHull.h.
|
private |
Used to define a default functor to allow default constructor
Definition at line 228 of file MelkmanConvexHull.h.
|
private |
first point used to reverse the convexhull container.
Definition at line 232 of file MelkmanConvexHull.h.
|
private |
Predicate devoted to the forward scan
Definition at line 224 of file MelkmanConvexHull.h.