ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkConstNeighborhoodIterator.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 itkConstNeighborhoodIterator_h
19#define itkConstNeighborhoodIterator_h
20
21#include <vector>
22#include <cstring>
23#include <iostream>
24#include "itkImage.h"
25#include "itkNeighborhood.h"
26#include "itkMacro.h"
28
29namespace itk
30{
50template <typename TImage, typename TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TImage>>
51class ITK_TEMPLATE_EXPORT ConstNeighborhoodIterator
52 : public Neighborhood<typename TImage::InternalPixelType *, TImage::ImageDimension>
53{
54public:
56 using InternalPixelType = typename TImage::InternalPixelType;
57 using PixelType = typename TImage::PixelType;
58
60 using DimensionValueType = unsigned int;
61
63 static constexpr DimensionValueType Dimension = TImage::ImageDimension;
64
68
70 using typename Superclass::OffsetType;
71 using typename Superclass::RadiusType;
72 using typename Superclass::SizeType;
73 using typename Superclass::Iterator;
74 using typename Superclass::ConstIterator;
75
77 using ImageType = TImage;
78 using RegionType = typename TImage::RegionType;
81
82
84 using BoundaryConditionType = TBoundaryCondition;
85
86 using OutputImageType = typename BoundaryConditionType::OutputImageType;
87
90
94 using NeighborhoodAccessorFunctorType = typename ImageType::NeighborhoodAccessorFunctorType;
95
99
102
104 ~ConstNeighborhoodIterator() override = default;
105
108
111 ConstNeighborhoodIterator(const SizeType & radius, const ImageType * ptr, const RegionType & region)
112 {
113 this->Initialize(radius, ptr, region);
114 for (DimensionValueType i = 0; i < Dimension; ++i)
115 {
116 m_InBounds[i] = false;
117 }
119 m_NeighborhoodAccessorFunctor = ptr->GetNeighborhoodAccessor();
120 m_NeighborhoodAccessorFunctor.SetBegin(ptr->GetBufferPointer());
121 }
122
124 Self &
125 operator=(const Self & orig);
126
128 void
129 PrintSelf(std::ostream &, Indent) const override;
130
135
138 GetBound() const
139 {
140 return m_Bound;
141 }
142
147 {
148 return m_Bound[n];
149 }
150
152 const InternalPixelType *
154 {
155 return (this->operator[]((this->Size()) >> 1));
156 }
157
160 PixelType
162 {
164 }
165
167 const ImageType *
169 {
170 return m_ConstImage;
171 }
172
175 IndexType
176 GetIndex() const
177 {
178 return m_Loop;
179 }
180
181 inline IndexType
183 {
184 return m_Loop + o;
185 }
186
189 NeighborhoodType
191
195 {
197 {
198 return (m_NeighborhoodAccessorFunctor.Get(this->operator[](i)));
199 }
200
201 OffsetType internalIndex;
202 OffsetType offset;
203
204 return this->IndexInBounds(i, internalIndex, offset)
205 ? m_NeighborhoodAccessorFunctor.Get(this->operator[](i))
206 : m_NeighborhoodAccessorFunctor.BoundaryCondition(internalIndex, offset, this, m_BoundaryCondition);
207 }
208
214 PixelType
215 GetPixel(NeighborIndexType n, bool & IsInBounds) const;
216
220 GetPixel(const OffsetType & o) const
221 {
222 bool inbounds;
223
224 return (this->GetPixel(this->GetNeighborhoodIndex(o), inbounds));
225 }
226
232 PixelType
233 GetPixel(const OffsetType & o, bool & IsInBounds) const
234 {
235 return (this->GetPixel(this->GetNeighborhoodIndex(o), IsInBounds));
236 }
237
241 PixelType
242 GetNext(const unsigned int axis, NeighborIndexType i) const
243 {
244 return (this->GetPixel(this->GetCenterNeighborhoodIndex() + (i * this->GetStride(axis))));
245 }
246
250 PixelType
251 GetNext(const unsigned int axis) const
252 {
253 return (this->GetPixel(this->GetCenterNeighborhoodIndex() + this->GetStride(axis)));
254 }
255
259 PixelType
260 GetPrevious(const unsigned int axis, NeighborIndexType i) const
261 {
262 return (this->GetPixel(this->GetCenterNeighborhoodIndex() - (i * this->GetStride(axis))));
263 }
264
268 PixelType
269 GetPrevious(const unsigned int axis) const
270 {
271 return (this->GetPixel(this->GetCenterNeighborhoodIndex() - this->GetStride(axis)));
272 }
273
276 IndexType
277 GetIndex(const OffsetType & o) const
278 {
279 return (this->GetIndex() + o);
280 }
281
284 IndexType
286 {
287 return (this->GetIndex() + this->GetOffset(i));
288 }
289
291 RegionType
292 GetRegion() const
293 {
294 return m_Region;
295 }
296
299 IndexType
301 {
302 return m_BeginIndex;
303 }
304
307 RegionType
309
313 {
314 return m_WrapOffset;
315 }
316
324 {
325 return m_WrapOffset[n];
326 }
327
329 void
331
334 void
336
339 void
340 Initialize(const SizeType & radius, const ImageType * ptr, const RegionType & region);
341
344 bool
345 IsAtBegin() const
346 {
347 return (this->GetCenterPointer() == m_Begin);
348 }
349
352 bool
353 IsAtEnd() const
354 {
355 if (this->GetCenterPointer() > m_End)
356 {
357 ExceptionObject e(__FILE__, __LINE__);
358 std::ostringstream msg;
359 msg << "In method IsAtEnd, CenterPointer = " << this->GetCenterPointer() << " is greater than End = " << m_End
360 << std::endl
361 << " " << *this;
362 e.SetDescription(msg.str().c_str());
363 throw e;
364 }
365 return (this->GetCenterPointer() == m_End);
366 }
367
372 Self &
374
379 Self &
381
385 bool
386 operator==(const Self & it) const
387 {
388 return it.GetCenterPointer() == this->GetCenterPointer();
389 }
390
392
396 bool
397 operator<(const Self & it) const
398 {
399 return this->GetCenterPointer() < it.GetCenterPointer();
400 }
401
405 bool
406 operator<=(const Self & it) const
407 {
408 return this->GetCenterPointer() <= it.GetCenterPointer();
409 }
410
414 bool
415 operator>(const Self & it) const
416 {
417 return this->GetCenterPointer() > it.GetCenterPointer();
418 }
419
423 bool
424 operator>=(const Self & it) const
425 {
426 return this->GetCenterPointer() >= it.GetCenterPointer();
427 }
428
433 void
434 SetLocation(const IndexType & position)
435 {
436 this->SetLoop(position);
437 this->SetPixelPointers(position);
438 }
439
443 Self &
445
449 Self &
451
454 operator-(const Self & b) const
455 {
456 return m_Loop - b.m_Loop;
457 }
458
462 bool
463 InBounds() const;
464
476 bool
477 IndexInBounds(const NeighborIndexType n, OffsetType & internalIndex, OffsetType & offset) const;
478
481 bool
483
489 void
494
497 void
502
504 void
505 SetBoundaryCondition(const TBoundaryCondition & c)
506 {
508 }
509
511 ImageBoundaryConditionPointerType
513 {
514 return m_BoundaryCondition;
515 }
516
518 void
523
524 void
529
530 void
535
536 bool
541
543 void
544 SetRegion(const RegionType & region);
545
546protected:
549 void
551 {
552 m_Loop = p;
553 m_IsInBoundsValid = false;
554 }
555
559 void
561
566 void
568
571 void
573 {
574 m_BeginIndex = start;
575 }
576
579 void
581
585
588
590 const InternalPixelType * m_Begin{ nullptr };
591
593 typename ImageType::ConstWeakPointer m_ConstImage{};
594
596 const InternalPixelType * m_End{ nullptr };
597
601
604
607
614
616 TBoundaryCondition m_InternalBoundaryCondition{};
617
623
626 mutable bool m_InBounds[Dimension]{ false };
627
629 mutable bool m_IsInBounds{ false };
630
634 mutable bool m_IsInBoundsValid{ false };
635
638
641
644
647};
648
649template <typename TImage>
650inline ConstNeighborhoodIterator<TImage>
653{
655 ret += ind;
656 return ret;
657}
658
659template <typename TImage>
660inline ConstNeighborhoodIterator<TImage>
663{
664 return (it + ind);
665}
666
667template <typename TImage>
668inline ConstNeighborhoodIterator<TImage>
671{
673 ret -= ind;
674 return ret;
675}
676} // namespace itk
677
678#ifndef ITK_MANUAL_INSTANTIATION
679# include "itkConstNeighborhoodIterator.hxx"
680#endif
681
682#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
IndexValueType GetBound(NeighborIndexType n) const
NeighborhoodType GetNeighborhood() const
void SetRegion(const RegionType &region)
PixelType GetNext(const unsigned int axis) const
~ConstNeighborhoodIterator() override=default
const InternalPixelType * GetCenterPointer() const
typename ImageType::NeighborhoodAccessorFunctorType NeighborhoodAccessorFunctorType
Neighborhood< PixelType, Self::Dimension > NeighborhoodType
Self & operator=(const Self &orig)
bool IndexInBounds(const NeighborIndexType n) const
void OverrideBoundaryCondition(const ImageBoundaryConditionPointerType i)
PixelType GetPixel(const NeighborIndexType i) const
void SetBeginIndex(const IndexType &start)
typename TImage::InternalPixelType InternalPixelType
const ImageBoundaryCondition< ImageType, OutputImageType > * ImageBoundaryConditionConstPointerType
void SetBoundaryCondition(const TBoundaryCondition &c)
Self & operator+=(const OffsetType &)
Neighborhood< InternalPixelType *, Self::Dimension > Superclass
OffsetType ComputeInternalIndex(const NeighborIndexType n) const
RegionType GetBoundingBoxAsImageRegion() const
PixelType GetPixel(const OffsetType &o, bool &IsInBounds) const
ConstNeighborhoodIterator(const ConstNeighborhoodIterator &)
void Initialize(const SizeType &radius, const ImageType *ptr, const RegionType &region)
OffsetValueType GetWrapOffset(NeighborIndexType n) const
void SetBound(const SizeType &)
void PrintSelf(std::ostream &, Indent) const override
IndexType GetIndex(NeighborIndexType i) const
PixelType GetPrevious(const unsigned int axis, NeighborIndexType i) const
IndexType GetIndex(const OffsetType &o) const
ImageBoundaryCondition< ImageType, OutputImageType > * ImageBoundaryConditionPointerType
Self & operator-=(const OffsetType &)
PixelType GetPrevious(const unsigned int axis) const
typename NeighborhoodType::NeighborIndexType NeighborIndexType
ConstNeighborhoodIterator(const SizeType &radius, const ImageType *ptr, const RegionType &region)
PixelType GetPixel(NeighborIndexType n, bool &IsInBounds) const
bool IndexInBounds(const NeighborIndexType n, OffsetType &internalIndex, OffsetType &offset) const
OffsetType operator-(const Self &b) const
PixelType GetPixel(const OffsetType &o) const
ImageBoundaryConditionPointerType GetBoundaryCondition() const
typename BoundaryConditionType::OutputImageType OutputImageType
IndexType GetFastIndexPlusOffset(const OffsetType &o) const
void SetLocation(const IndexType &position)
PixelType GetNext(const unsigned int axis, NeighborIndexType i) const
Standard exception handling object.
virtual void SetDescription(const std::string &s)
A virtual base object that defines an interface to a class of boundary condition objects for use by n...
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual NeighborIndexType GetNeighborhoodIndex(const OffsetType &) const
typename AllocatorType::const_iterator ConstIterator
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long OffsetValueType
Definition itkIntTypes.h:97
ConstNeighborhoodIterator< TImage > operator-(const ConstNeighborhoodIterator< TImage > &it, const typename ConstNeighborhoodIterator< TImage >::OffsetType &ind)
ConstNeighborhoodIterator< TImage > operator+(const ConstNeighborhoodIterator< TImage > &it, const typename ConstNeighborhoodIterator< TImage >::OffsetType &ind)
long IndexValueType
Definition itkIntTypes.h:93
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69