DGtal 1.4.0
Loading...
Searching...
No Matches
MelkmanConvexHull.ih
1/**
2 * This program is free software: you can redistribute it and/or modify
3 * it under the terms of the GNU Lesser General Public License as
4 * published by the Free Software Foundation, either version 3 of the
5 * License, or (at your option) any later version.
6 *
7 * This program is distributed in the hope that it will be useful,
8 * but WITHOUT ANY WARRANTY; without even the implied warranty of
9 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
10 * GNU General Public License for more details.
11 *
12 * You should have received a copy of the GNU General Public License
13 * along with this program. If not, see <http://www.gnu.org/licenses/>.
14 *
15 **/
16
17/**
18 * @file MelkmanConvexHull.ih
19 * @author Tristan Roussillon (\c tristan.roussillon@liris.cnrs.fr )
20 * Laboratoire d'InfoRmatique en Image et Systèmes d'information - LIRIS (CNRS, UMR 5205), CNRS, France
21 *
22 * @date 2013/12/20
23 *
24 * Implementation of inline methods defined in MelkmanConvexHull.h
25 *
26 * This file is part of the DGtal library.
27 */
28
29
30//////////////////////////////////////////////////////////////////////////////
31#include <cstdlib>
32#include <boost/utility.hpp>
33#include <boost/next_prior.hpp>
34//////////////////////////////////////////////////////////////////////////////
35
36///////////////////////////////////////////////////////////////////////////////
37// IMPLEMENTATION of inline methods.
38///////////////////////////////////////////////////////////////////////////////
39
40///////////////////////////////////////////////////////////////////////////////
41// ----------------------------------------------------------------------------
42template <typename TPoint, typename TOrientationFunctor>
43inline
44DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::MelkmanConvexHull( Alias<Functor> aFunctor )
45 : myContainer(),
46 myBackwardPredicate( aFunctor ),
47 myForwardPredicate( aFunctor )
48{
49}
50
51// ----------------------------------------------------------------------------
52template <typename TPoint, typename TOrientationFunctor>
53inline
54DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::MelkmanConvexHull()
55 : myContainer(),
56 myBackwardPredicate(myDefaultFunctor),
57 myForwardPredicate(myDefaultFunctor)
58{
59}
60
61// ----------------------------------------------------------------------------
62template <typename TPoint, typename TOrientationFunctor>
63inline
64void
65DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::add(const Point& aPoint)
66{
67
68 using namespace DGtal::functions::Hull2D;
69
70 if (myContainer.size() < 4)
71 {
72 if (myContainer.size() == 0)
73 {
74 myFirstPoint = aPoint;
75
76 myContainer.push_back( aPoint );
77 myContainer.push_front( aPoint );
78 }
79 else if (myContainer.size() == 2)
80 {
81 myContainer.pop_back();
82 myContainer.push_back( aPoint );
83 myContainer.push_front( aPoint );
84 }
85 else if (myContainer.size() == 3)
86 {
87 //required to deal with the case where the k first input points are aligned.
88 if ( !myBackwardPredicate( *boost::next(myContainer.rbegin()), myContainer.back(), aPoint ) )
89 myContainer.pop_back();
90 myContainer.push_back( aPoint );
91 if ( !myForwardPredicate( *boost::next(myContainer.begin()), myContainer.front(), aPoint ) )
92 myContainer.pop_front();
93 myContainer.push_front( aPoint );
94 }
95 }
96 else
97 {
98 if ( ( !myBackwardPredicate( *boost::next(myContainer.rbegin()), myContainer.back(), aPoint ) ||
99 !myForwardPredicate( *boost::next(myContainer.begin()), myContainer.front(), aPoint ) ) )
100 {
101 //backward scan
102 updateHullWithAdaptedStack( backStack(myContainer), aPoint, myBackwardPredicate );
103 myContainer.push_back( aPoint );
104
105 //forward scan
106 updateHullWithAdaptedStack( frontStack(myContainer), aPoint, myForwardPredicate );
107 myContainer.push_front( aPoint );
108 }
109 }
110}
111
112// ----------------------------------------------------------------------------
113template <typename TPoint, typename TOrientationFunctor>
114inline
115typename DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::ConstIterator
116DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::begin() const
117{
118 if (myContainer.size() == 0)
119 return myContainer.end();
120 else
121 return boost::next(myContainer.begin());
122}
123
124// ----------------------------------------------------------------------------
125template <typename TPoint, typename TOrientationFunctor>
126inline
127typename DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::ConstIterator
128DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::end() const
129{
130 return myContainer.end();
131}
132
133// ----------------------------------------------------------------------------
134template <typename TPoint, typename TOrientationFunctor>
135inline
136void
137DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::selfDisplay ( std::ostream & out ) const
138{
139 out << "[MelkmanConvexHull]" << " #";
140 if ( myContainer.size() == 0 )
141 out << " 0 " << std::endl;
142 else
143 out << myContainer.size() - 1 << std::endl;
144 std::copy( myContainer.begin(), myContainer.end(),
145 std::ostream_iterator<Point>( out, "," ) );
146 out << std::endl;
147}
148
149// ----------------------------------------------------------------------------
150template <typename TPoint, typename TOrientationFunctor>
151inline
152bool
153DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::isValid() const
154{
155 return true;
156}
157
158// ----------------------------------------------------------------------------
159template <typename TPoint, typename TOrientationFunctor>
160inline
161DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>&
162DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::operator= (const Self & mch)
163{
164 myContainer = mch.myContainer;
165 return *this;
166}
167
168// ----------------------------------------------------------------------------
169template <typename TPoint, typename TOrientationFunctor>
170inline
171const TPoint &
172DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::operator[](unsigned int i) const
173{
174 ASSERT( i < myContainer.size() );
175 return myContainer[i];
176}
177
178// ----------------------------------------------------------------------------
179template <typename TPoint, typename TOrientationFunctor>
180inline
181size_t
182DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::size() const
183{
184 // by definition the first and last points of the deque are the same.
185 return ( myContainer.size()-1u );
186
187}
188// ----------------------------------------------------------------------------
189template <typename TPoint, typename TOrientationFunctor>
190inline
191void
192DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::clear()
193{
194 myContainer.clear();
195}
196
197// ----------------------------------------------------------------------------
198template <typename TPoint, typename TOrientationFunctor>
199inline
200void
201DGtal::MelkmanConvexHull<TPoint, TOrientationFunctor>::reverse()
202{
203 // if convexhull is reduced into a single segment no need to reverse anything.
204 if(myContainer.size()<=3)
205 return;
206 TPoint theLast = myContainer.back();
207 myContainer.pop_back();
208 bool foundFirst = myContainer.back() == myFirstPoint;
209 while(!foundFirst){
210 myContainer.push_front(myContainer.back());
211 myContainer.pop_back();
212 foundFirst = myContainer.back() == myFirstPoint;
213 }
214 myContainer.push_front(myContainer.back());
215 myFirstPoint = theLast;
216}
217
218
219///////////////////////////////////////////////////////////////////////////////
220// Implementation of inline functions //
221
222template <typename TPoint, typename TOrientationFunctor>
223inline
224std::ostream&
225DGtal::operator<< ( std::ostream & out,
226 const MelkmanConvexHull<TPoint, TOrientationFunctor> & object )
227{
228 object.selfDisplay( out );
229 return out;
230}
231
232///////////////////////////////////////////////////////////////////////////////
233template <typename ForwardIterator,
234 typename OutputIterator,
235 typename Functor >
236inline
237void DGtal::functions::Hull2D::
238melkmanConvexHullAlgorithm(const ForwardIterator& itb, const ForwardIterator& ite,
239 OutputIterator res,
240 Functor& aFunctor )
241{
242 BOOST_CONCEPT_ASSERT(( boost_concepts::ForwardTraversalConcept<ForwardIterator> ));
243 BOOST_CONCEPT_ASSERT(( boost_concepts::ReadableIteratorConcept<ForwardIterator> ));
244 typedef typename DGtal::IteratorCirculatorTraits<ForwardIterator>::Value Point;
245 BOOST_CONCEPT_ASSERT(( boost_concepts::IncrementableIteratorConcept<OutputIterator> ));
246 BOOST_CONCEPT_ASSERT(( boost_concepts::WritableIteratorConcept<OutputIterator,Point> ));
247
248 if ( itb != ite )
249 {
250
251 //convex hull computation with Melkman's algorithm
252 DGtal::MelkmanConvexHull<Point, Functor> ch( aFunctor );
253 for (ForwardIterator it = itb; it != ite; ++it)
254 ch.add( *it );
255
256 std::copy( ch.begin(), ch.end(), res );
257 }
258}
259// //
260///////////////////////////////////////////////////////////////////////////////
261
262