ITK  5.4.0
Insight Toolkit
itkImageConstIterator.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 itkImageConstIterator_h
19 #define itkImageConstIterator_h
20 
21 #include "itkImage.h"
22 #include "itkIndex.h"
23 #include "itkNumericTraits.h"
24 
25 namespace itk
26 {
83 template <typename TImage>
84 class ITK_TEMPLATE_EXPORT ImageConstIterator
85 {
86 public:
89 
94  static constexpr unsigned int ImageIteratorDimension = TImage::ImageDimension;
95 
97  itkVirtualGetNameOfClassMacro(ImageConstIterator);
98 
100  using IndexType = typename TImage::IndexType;
101 
103  using SizeType = typename TImage::SizeType;
104 
106  using OffsetType = typename TImage::OffsetType;
107 
109  using RegionType = typename TImage::RegionType;
110 
112  using ImageType = TImage;
113 
117  using PixelContainer = typename TImage::PixelContainer;
119 
121  using InternalPixelType = typename TImage::InternalPixelType;
122 
124  using PixelType = typename TImage::PixelType;
125 
128  using AccessorType = typename TImage::AccessorType;
129  using AccessorFunctorType = typename TImage::AccessorFunctorType;
130 
134  : m_Region()
135  , m_PixelAccessor()
136  , m_PixelAccessorFunctor()
137  {
138  m_Image = nullptr;
139  m_Buffer = nullptr;
140  m_Offset = 0;
141  m_BeginOffset = 0;
142  m_EndOffset = 0;
143  m_PixelAccessorFunctor.SetBegin(m_Buffer);
144  }
148  virtual ~ImageConstIterator() = default;
149 
153  {
154  m_Image = it.m_Image; // copy the smart pointer
155 
156  m_Region = it.m_Region;
157 
158  m_Buffer = it.m_Buffer;
159  m_Offset = it.m_Offset;
160  m_BeginOffset = it.m_BeginOffset;
161  m_EndOffset = it.m_EndOffset;
162  m_PixelAccessor = it.m_PixelAccessor;
163  m_PixelAccessorFunctor = it.m_PixelAccessorFunctor;
164  m_PixelAccessorFunctor.SetBegin(m_Buffer);
165  }
166 
169  ImageConstIterator(const ImageType * ptr, const RegionType & region)
170  {
171  m_Image = ptr;
172  m_Buffer = m_Image->GetBufferPointer();
175  SetRegion(region);
176 
177  m_PixelAccessor = ptr->GetPixelAccessor();
178  m_PixelAccessorFunctor.SetPixelAccessor(m_PixelAccessor);
179  m_PixelAccessorFunctor.SetBegin(m_Buffer);
180  }
181 
184  Self &
185  operator=(const Self & it)
186  {
187  if (this != &it)
188  {
189  m_Image = it.m_Image; // copy the smart pointer
190  m_Region = it.m_Region;
193  m_Buffer = it.m_Buffer;
194  m_Offset = it.m_Offset;
195  m_BeginOffset = it.m_BeginOffset;
196  m_EndOffset = it.m_EndOffset;
197  m_PixelAccessor = it.m_PixelAccessor;
198  m_PixelAccessorFunctor = it.m_PixelAccessorFunctor;
199  m_PixelAccessorFunctor.SetBegin(m_Buffer);
200  }
201  return *this;
202  }
203 
205  virtual void
206  SetRegion(const RegionType & region)
207  {
208  m_Region = region;
209 
210  if (region.GetNumberOfPixels() > 0) // If region is non-empty
211  {
212  const RegionType & bufferedRegion = m_Image->GetBufferedRegion();
213  itkAssertOrThrowMacro((bufferedRegion.IsInside(m_Region)),
214  "Region " << m_Region << " is outside of buffered region " << bufferedRegion);
215  }
216 
217  // Compute the start offset
218  m_Offset = m_Image->ComputeOffset(m_Region.GetIndex());
219  m_BeginOffset = m_Offset;
220 
221  // Compute the end offset. If any component of m_Region.GetSize()
222  // is zero, the region is not valid and we set the EndOffset
223  // to be same as BeginOffset so that iterator end condition is met
224  // immediately.
225  IndexType ind(m_Region.GetIndex());
226  SizeType size(m_Region.GetSize());
227  if (m_Region.GetNumberOfPixels() == 0)
228  {
229  // region is empty, probably has a size of 0 along one dimension
230  m_EndOffset = m_BeginOffset;
231  }
232  else
233  {
234  for (unsigned int i = 0; i < TImage::ImageDimension; ++i)
235  {
236  ind[i] += (static_cast<IndexValueType>(size[i]) - 1);
237  }
238  m_EndOffset = m_Image->ComputeOffset(ind);
239  ++m_EndOffset;
240  }
241  }
242 
244  static unsigned int
246  {
247  return TImage::ImageDimension;
248  }
249 
252  bool
253  operator==(const Self & it) const
254  {
255  // two iterators are the same if they "point to" the same memory location
256  return (m_Buffer + m_Offset) == (it.m_Buffer + it.m_Offset);
257  }
258 
259  ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION(Self);
260 
263  bool
264  operator<=(const Self & it) const
265  {
266  // an iterator is "less than" another if it "points to" a lower
267  // memory location
268  return (m_Buffer + m_Offset) <= (it.m_Buffer + it.m_Offset);
269  }
270 
273  bool
274  operator<(const Self & it) const
275  {
276  // an iterator is "less than" another if it "points to" a lower
277  // memory location
278  return (m_Buffer + m_Offset) < (it.m_Buffer + it.m_Offset);
279  }
280 
283  bool
284  operator>=(const Self & it) const
285  {
286  // an iterator is "greater than" another if it "points to" a higher
287  // memory location
288  return (m_Buffer + m_Offset) >= (it.m_Buffer + it.m_Offset);
289  }
290 
293  bool
294  operator>(const Self & it) const
295  {
296  // an iterator is "greater than" another if it "points to" a higher
297  // memory location
298  return (m_Buffer + m_Offset) > (it.m_Buffer + it.m_Offset);
299  }
300 
305  const IndexType
306  GetIndex() const
307  {
308  return m_Image->ComputeIndex(static_cast<OffsetValueType>(m_Offset));
309  }
310 
313  virtual void
314  SetIndex(const IndexType & ind)
315  {
316  m_Offset = m_Image->ComputeOffset(ind);
317  }
318 
321  const RegionType &
322  GetRegion() const
323  {
324  return m_Region;
325  }
326 
328  const ImageType *
329  GetImage() const
330  {
331  return m_Image.GetPointer();
332  }
333 
335  PixelType
336  Get() const
337  {
338  return m_PixelAccessorFunctor.Get(*(m_Buffer + m_Offset));
339  }
340 
344  const PixelType &
345  Value() const
346  {
347  return *(m_Buffer + m_Offset);
348  }
349 
352  void
354  {
355  m_Offset = m_BeginOffset;
356  }
357 
360  void
362  {
363  m_Offset = m_EndOffset;
364  }
365 
368  bool
369  IsAtBegin() const
370  {
371  return (m_Offset == m_BeginOffset);
372  }
373 
376  bool
377  IsAtEnd() const
378  {
379  return (m_Offset == m_EndOffset);
380  }
381 
382 protected: // made protected so other iterators can access
383  typename TImage::ConstWeakPointer m_Image{};
384 
385  RegionType m_Region{}; // region to iterate over
386 
387  OffsetValueType m_Offset{};
388  OffsetValueType m_BeginOffset{}; // offset to first pixel in region
389  OffsetValueType m_EndOffset{}; // offset to one pixel past last pixel in region
390 
391  const InternalPixelType * m_Buffer{};
392 
393  AccessorType m_PixelAccessor{};
394  AccessorFunctorType m_PixelAccessorFunctor{};
395 };
396 } // end namespace itk
397 
398 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ImageConstIterator::ImageConstIterator
ImageConstIterator()
Definition: itkImageConstIterator.h:133
itk::ImageConstIterator< ImageType >::SizeType
typename ImageType ::SizeType SizeType
Definition: itkImageConstIterator.h:103
itk::ImageConstIterator< ImageType >::PixelType
typename ImageType ::PixelType PixelType
Definition: itkImageConstIterator.h:124
itk::ImageConstIterator::m_PixelAccessorFunctor
AccessorFunctorType m_PixelAccessorFunctor
Definition: itkImageConstIterator.h:394
itk::operator<
bool operator<(const Index< VDimension > &one, const Index< VDimension > &two)
Definition: itkIndex.h:559
itk::operator<=
bool operator<=(const Index< VDimension > &one, const Index< VDimension > &two)
Definition: itkIndex.h:573
itk::ImageConstIterator::m_BeginOffset
OffsetValueType m_BeginOffset
Definition: itkImageConstIterator.h:388
itk::ImageConstIterator::m_Buffer
const InternalPixelType * m_Buffer
Definition: itkImageConstIterator.h:391
itk::ImageConstIterator< ImageType >::OffsetType
typename ImageType ::OffsetType OffsetType
Definition: itkImageConstIterator.h:106
itk::ImageConstIterator< ImageType >::AccessorType
typename ImageType ::AccessorType AccessorType
Definition: itkImageConstIterator.h:128
itk::ImageConstIterator::GetImage
const ImageType * GetImage() const
Definition: itkImageConstIterator.h:329
itk::ImageConstIterator< ImageType >::RegionType
typename ImageType ::RegionType RegionType
Definition: itkImageConstIterator.h:109
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::ImageConstIterator::SetRegion
virtual void SetRegion(const RegionType &region)
Definition: itkImageConstIterator.h:206
itk::ImageConstIterator::operator=
Self & operator=(const Self &it)
Definition: itkImageConstIterator.h:185
itk::ImageConstIterator::ImageConstIterator
ImageConstIterator(const Self &it)
Definition: itkImageConstIterator.h:152
itk::ImageConstIterator::GoToBegin
void GoToBegin()
Definition: itkImageConstIterator.h:353
itk::ImageConstIterator::IsAtEnd
bool IsAtEnd() const
Definition: itkImageConstIterator.h:377
itk::ImageConstIterator::m_Region
RegionType m_Region
Definition: itkImageConstIterator.h:385
itk::ImageConstIterator< ImageType >::PixelContainerPointer
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkImageConstIterator.h:118
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageConstIterator::operator>=
bool operator>=(const Self &it) const
Definition: itkImageConstIterator.h:284
itk::ImageConstIterator::GetIndex
const IndexType GetIndex() const
Definition: itkImageConstIterator.h:306
itk::ImageConstIterator::GetImageIteratorDimension
static unsigned int GetImageIteratorDimension()
Definition: itkImageConstIterator.h:245
itk::ImageConstIterator::m_Image
TImage::ConstWeakPointer m_Image
Definition: itkImageConstIterator.h:383
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageConstIterator::GoToEnd
void GoToEnd()
Definition: itkImageConstIterator.h:361
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:94
itk::ImageConstIterator::m_Offset
OffsetValueType m_Offset
Definition: itkImageConstIterator.h:387
itk::ImageConstIterator::m_PixelAccessor
AccessorType m_PixelAccessor
Definition: itkImageConstIterator.h:393
itk::ImageConstIterator< ImageType >::InternalPixelType
typename ImageType ::InternalPixelType InternalPixelType
Definition: itkImageConstIterator.h:121
itkIndex.h
itk::ImageConstIterator::IsAtBegin
bool IsAtBegin() const
Definition: itkImageConstIterator.h:369
itk::ImageConstIterator::Value
const PixelType & Value() const
Definition: itkImageConstIterator.h:345
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageConstIterator::Get
PixelType Get() const
Definition: itkImageConstIterator.h:336
itk::ImageConstIterator
A multi-dimensional image iterator templated over image type.
Definition: itkImageConstIterator.h:84
itk::ImageConstIterator< ImageType >::ImageType
ImageType ImageType
Definition: itkImageConstIterator.h:112
itk::ImageConstIterator< ImageType >::PixelContainer
typename ImageType ::PixelContainer PixelContainer
Definition: itkImageConstIterator.h:117
itk::ImageConstIterator::GetRegion
const RegionType & GetRegion() const
Definition: itkImageConstIterator.h:322
itk::ImageConstIterator::SetIndex
virtual void SetIndex(const IndexType &ind)
Definition: itkImageConstIterator.h:314
itkNumericTraits.h
itk::ImageConstIterator< ImageType >::IndexType
typename ImageType ::IndexType IndexType
Definition: itkImageConstIterator.h:100
itk::ImageConstIterator::ImageConstIterator
ImageConstIterator(const ImageType *ptr, const RegionType &region)
Definition: itkImageConstIterator.h:169
itk::ImageConstIterator::operator>
bool operator>(const Self &it) const
Definition: itkImageConstIterator.h:294
itk::ImageConstIterator::m_EndOffset
OffsetValueType m_EndOffset
Definition: itkImageConstIterator.h:389
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::ImageConstIterator::operator==
bool operator==(const Self &it) const
Definition: itkImageConstIterator.h:253