ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkImageRegionReverseConstIterator.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkImageRegionReverseConstIterator_h
19 #define itkImageRegionReverseConstIterator_h
20 
22 #include "itkImageRegionIterator.h"
23 
24 namespace itk
25 {
102 template< typename TImage >
104 {
105 public:
106 
110 
114 
117  typedef typename Superclass::SizeType SizeType;
118 
122 
125 
129 
134  typedef typename PixelContainer::Pointer PixelContainerPointer;
135 
138 
141 
145 
148 
151  {
152  m_SpanBeginOffset = 0;
153  m_SpanEndOffset = 0;
154  }
155 
159  Superclass(ptr, region)
160  {
162  m_SpanEndOffset = this->m_BeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
163  }
165 
174  {
175  IndexType ind = this->GetIndex();
176 
177  m_SpanBeginOffset = this->m_Offset + static_cast< OffsetValueType >( this->m_Region.GetSize()[0] )
178  - ( ind[0] - this->m_Region.GetIndex()[0] );
179  m_SpanEndOffset = m_SpanBeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
180  }
181 
185  {
186  IndexType ind = this->GetIndex();
187 
188  m_SpanBeginOffset = this->m_Offset + static_cast< OffsetValueType >( this->m_Region.GetSize()[0] )
189  - ( ind[0] - this->m_Region.GetIndex()[0] );
190  m_SpanEndOffset = m_SpanBeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
191  }
192 
196  {
197  IndexType ind = this->GetIndex();
198 
199  m_SpanBeginOffset = this->m_Offset + static_cast< OffsetValueType >( this->m_Region.GetSize()[0] )
200  - ( ind[0] - this->m_Region.GetIndex()[0] );
201  m_SpanEndOffset = m_SpanBeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
202  }
203 
206  void GoToBegin()
207  {
209 
210  // reset the span offsets
212  m_SpanEndOffset = this->m_BeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
213  }
214 
217  void GoToEnd()
218  {
220 
221  // reset the span offsets
223  m_SpanBeginOffset = m_SpanEndOffset + static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
224  }
225 
229  itkLegacyMacro(Self Begin(void) const);
230 
235  itkLegacyMacro(Self End(void) const);
236 
240  void SetIndex(const IndexType & ind) ITK_OVERRIDE
241  {
243  m_SpanBeginOffset = this->m_Offset + static_cast< OffsetValueType >( this->m_Region.GetSize()[0] )
244  - ( ind[0] - this->m_Region.GetIndex()[0] );
245  m_SpanEndOffset = m_SpanBeginOffset - static_cast< OffsetValueType >( this->m_Region.GetSize()[0] );
246  }
248 
257  Self &
259  {
260  if ( --this->m_Offset <= m_SpanEndOffset )
261  {
262  // We have past the beginning of the span (row), need to wrap around.
263 
264  // First move forward one pixel, because we are going to use a different
265  // algorithm to compute the next pixel
266  this->m_Offset++;
267 
268  // Get the index of the first pixel on the span (row)
270  ind = this->m_Image->ComputeIndex( static_cast< OffsetValueType >( this->m_Offset ) );
271 
273  startIndex = this->m_Region.GetIndex();
275  size = this->m_Region.GetSize();
276 
277  // Deccrement along a row, then wrap at the beginning of the region row.
278  bool done;
279  unsigned int dim;
280 
281  // Check to see if we are past the first pixel in the region
282  // Note that --ind[0] moves to the previous pixel along the row.
283  done = ( --ind[0] == startIndex[0] - 1 );
284  for ( unsigned int i = 1; done && i < this->ImageIteratorDimension; i++ )
285  {
286  done = ( ind[i] == startIndex[i] );
287  }
288 
289  // if the iterator is outside the region (but not past region begin) then
290  // we need to wrap around the region
291  dim = 0;
292  if ( !done )
293  {
294  while ( ( dim < this->ImageIteratorDimension - 1 )
295  && ( ind[dim] < startIndex[dim] ) )
296  {
297  ind[dim] = startIndex[dim] + static_cast< OffsetValueType >( size[dim] ) - 1;
298  ind[++dim]--;
299  }
300  }
301  this->m_Offset = this->m_Image->ComputeOffset(ind);
302  m_SpanBeginOffset = this->m_Offset;
303  m_SpanEndOffset = m_SpanBeginOffset - static_cast< OffsetValueType >( size[0] );
304  }
305  return *this;
306  }
307 
317  {
318  if ( ++this->m_Offset >= m_SpanBeginOffset )
319  {
320  // We have reached the end of the span (row), need to wrap around.
321 
322  // First back up one pixel, because we are going to use a different
323  // algorithm to compute the next pixel
324  --this->m_Offset;
325 
326  // Get the index of the last pixel on the span (row)
328  ind = this->m_Image->ComputeIndex( static_cast< OffsetValueType >( this->m_Offset ) );
329 
330  const typename ImageIterator< TImage >::IndexType &
331  startIndex = this->m_Region.GetIndex();
332  const typename ImageIterator< TImage >::SizeType &
333  size = this->m_Region.GetSize();
334 
335  // Increment along a row, then wrap at the end of the region row.
336  bool done;
337  unsigned int dim;
338 
339  // Check to see if we are past the last pixel in the region
340  // Note that ++ind[0] moves to the next pixel along the row.
341  done = ( ++ind[0] == startIndex[0] + static_cast< OffsetValueType >( size[0] ) );
342  for ( unsigned int i = 1; done && i < this->ImageIteratorDimension; i++ )
343  {
344  done = ( ind[i] == startIndex[i] + static_cast< OffsetValueType >( size[i] ) - 1 );
345  }
346 
347  // if the iterator is outside the region (but not past region end) then
348  // we need to wrap around the region
349  dim = 0;
350  if ( !done )
351  {
352  while ( ( dim < this->ImageIteratorDimension - 1 )
353  && ( ind[dim] > startIndex[dim] + static_cast< OffsetValueType >( size[dim] ) - 1 ) )
354  {
355  ind[dim] = startIndex[dim];
356  ind[++dim]++;
357  }
358  }
359  this->m_Offset = this->m_Image->ComputeOffset(ind);
360  m_SpanBeginOffset = this->m_Offset;
361  m_SpanEndOffset = this->m_Offset - static_cast< OffsetValueType >( size[0] );
362  }
363  return *this;
364  }
365 
366 protected:
367  SizeValueType m_SpanBeginOffset; // offset to last pixel in the row
368  SizeValueType m_SpanEndOffset; // offset to one pixel before the row
369 };
370 } // end namespace itk
371 
372 #ifndef ITK_MANUAL_INSTANTIATION
373 #include "itkImageRegionReverseConstIterator.hxx"
374 #endif
375 
376 #endif
Superclass::SizeType SizeType
signed long OffsetValueType
Definition: itkIntTypes.h:154
ImageRegionReverseConstIterator(const ImageReverseConstIterator< TImage > &it)
static const unsigned int ImageIteratorDimension
A multi-dimensional image iterator designed to walk a specified image region in reverse.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
ImageRegionReverseConstIterator(const ImageConstIterator< TImage > &it)
A multi-dimensional image iterator templated over image type.
ImageRegionReverseConstIterator(const ImageRegionIterator< TImage > &it)
Multi-dimensional image iterator.
virtual void SetIndex(const IndexType &ind)
ImageRegionReverseConstIterator(const ImageType *ptr, const RegionType &region)
Superclass::IndexType IndexType
itkLegacyMacro(Self Begin(void) const)
A multi-dimensional iterator templated over image type that walks a region of pixels.