ITK  6.0.0
Insight Toolkit
itkImageToNeighborhoodSampleAdaptor.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 itkImageToNeighborhoodSampleAdaptor_h
19 #define itkImageToNeighborhoodSampleAdaptor_h
20 
21 #include <typeinfo>
22 #include <vector>
23 #include <iostream>
24 
25 #include "itkImage.h"
26 #include "itkListSample.h"
27 #include "itkSmartPointer.h"
29 #include "itkMacro.h"
32 
33 namespace itk
34 {
35 namespace Statistics
36 {
37 
54 template <typename TImage, typename TBoundaryCondition>
55 class ITK_TEMPLATE_EXPORT ImageToNeighborhoodSampleAdaptor
56  : public ListSample<std::vector<ConstNeighborhoodIterator<TImage, TBoundaryCondition>>>
57 {
58 public:
59  ITK_DISALLOW_COPY_AND_MOVE(ImageToNeighborhoodSampleAdaptor);
60 
63 
65 
68 
70  itkOverrideGetNameOfClassMacro(ImageToNeighborhoodSampleAdaptor);
71 
73  itkNewMacro(Self);
74 
76  using ImageType = TImage;
77  using ImagePointer = typename ImageType::Pointer;
79  using IndexType = typename ImageType::IndexType;
80  using OffsetType = typename ImageType::OffsetType;
82  using PixelType = typename ImageType::PixelType;
83  using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
85  using OffsetTableType = typename RegionType::OffsetTableType;
86  using SizeType = typename ImageType::SizeType;
88 
96 
99  using MeasurementVectorType = typename std::vector<ConstNeighborhoodIterator<TImage, TBoundaryCondition>>;
100  using ValueType = typename MeasurementVectorType::value_type;
102 
103  using typename Superclass::AbsoluteFrequencyType;
104  using typename Superclass::TotalAbsoluteFrequencyType;
105  using typename Superclass::MeasurementVectorSizeType;
106  using typename Superclass::InstanceIdentifier;
107 
109  void
110  SetImage(const TImage * image);
111 
113  const TImage *
114  GetImage() const;
115 
117  void
118  SetRadius(const NeighborhoodRadiusType & radius);
119 
122  GetRadius() const;
123 
125  void
126  SetRegion(const RegionType & region);
127 
129  RegionType
130  GetRegion() const;
131 
132  void
133  SetUseImageRegion(const bool flag);
134 
136  itkGetConstMacro(UseImageRegion, bool);
137 
139  itkBooleanMacro(UseImageRegion);
140 
141 
144  Size() const override;
145 
147  const MeasurementVectorType &
148  GetMeasurementVector(InstanceIdentifier id) const override;
149 
152  GetFrequency(InstanceIdentifier id) const override;
153 
156  GetTotalFrequency() const override;
157 
164  {
166 
167  public:
168  ConstIterator() = delete;
169 
170  ConstIterator(const ImageToNeighborhoodSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
171 
173  {
174  m_MeasurementVectorCache = iter.m_MeasurementVectorCache;
175  m_InstanceIdentifier = iter.m_InstanceIdentifier;
176  }
177 
178  ConstIterator &
179  operator=(const ConstIterator & iter)
180  {
181  m_MeasurementVectorCache = iter.m_MeasurementVectorCache;
182  m_InstanceIdentifier = iter.m_InstanceIdentifier;
183  return *this;
184  }
185 
187  GetFrequency() const
188  {
189  return 1;
190  }
191 
192  const MeasurementVectorType &
194  {
195  return this->m_MeasurementVectorCache;
196  }
197 
200  {
201  return m_InstanceIdentifier;
202  }
203 
204  ConstIterator &
206  {
207  ++(m_MeasurementVectorCache[0]);
208  ++m_InstanceIdentifier;
209  return *this;
210  }
211 
212  bool
213  operator==(const ConstIterator & it) const
214  {
215  return (m_MeasurementVectorCache[0] == it.m_MeasurementVectorCache[0]);
216  }
217 
218  ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION(ConstIterator);
219 
220  protected:
221  // This method should only be available to the ListSample class
223  {
224  this->m_MeasurementVectorCache.clear();
225  this->m_MeasurementVectorCache.push_back(iter);
226  m_InstanceIdentifier = iid;
227  }
228 
229  private:
232  };
233 
239  class Iterator : public ConstIterator
240  {
241 
243 
244  public:
245  Iterator() = delete;
246  Iterator(const Self * adaptor) = delete;
247  Iterator(const ConstIterator & it) = delete;
248  ConstIterator &
249  operator=(const ConstIterator & it) = delete;
250 
251  Iterator(Self * adaptor)
252  : ConstIterator(adaptor)
253  {}
254 
255  Iterator(const Iterator & iter)
256  : ConstIterator(iter)
257  {}
258 
259  Iterator &
260  operator=(const Iterator & iter)
261  {
262  this->ConstIterator::operator=(iter);
263  return *this;
264  }
265 
266  protected:
267  // This copy constructor is actually used in Iterator Begin()!
269  : ConstIterator(iter, iid)
270  {}
271  };
272 
274  Iterator
276  {
277  NeighborhoodIteratorType nIterator(m_Radius, m_Image, m_Region);
278  nIterator.GoToBegin();
279  Iterator iter(nIterator, 0);
280  return iter;
281  }
285  Iterator
286  End()
287  {
288  NeighborhoodIteratorType nIterator(m_Radius, m_Image, m_Region);
289  nIterator.GoToEnd();
290  Iterator iter(nIterator, m_Region.GetNumberOfPixels());
291  return iter;
292  }
297  ConstIterator
298  Begin() const
299  {
300  NeighborhoodIteratorType nIterator(m_Radius, m_Image, m_Region);
301  nIterator.GoToBegin();
302  ConstIterator iter(nIterator, 0);
303  return iter;
304  }
308  ConstIterator
309  End() const
310  {
311  NeighborhoodIteratorType nIterator(m_Radius, m_Image, m_Region);
312  nIterator.GoToEnd();
313  ConstIterator iter(nIterator, m_Region.GetNumberOfPixels());
314  return iter;
315  }
318 protected:
320  ~ImageToNeighborhoodSampleAdaptor() override = default;
321  void
322  PrintSelf(std::ostream & os, Indent indent) const override;
323 
324 private:
325  ImageConstPointer m_Image{};
326  mutable MeasurementVectorType m_MeasurementVectorInternal{};
327  mutable InstanceIdentifier m_InstanceIdentifierInternal{};
328  mutable IndexType m_NeighborIndexInternal{};
330  RegionType m_Region{};
331  bool m_UseImageRegion{ true };
332  OffsetTableType m_OffsetTable{};
333 
334 }; // end of class ImageToNeighborhoodSampleAdaptor
335 
336 } // end of namespace Statistics
337 
338 template <typename TImage, typename TBoundaryCondition>
339 std::ostream &
340 operator<<(std::ostream & os, const std::vector<itk::ConstNeighborhoodIterator<TImage, TBoundaryCondition>> & mv);
341 
342 } // end of namespace itk
343 
344 #ifndef ITK_MANUAL_INSTANTIATION
345 # include "itkImageToNeighborhoodSampleAdaptor.hxx"
346 #endif
347 
348 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Statistics::ImageToNeighborhoodSampleAdaptor::OffsetType
typename ImageType::OffsetType OffsetType
Definition: itkImageToNeighborhoodSampleAdaptor.h:80
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::operator=
ConstIterator & operator=(const ConstIterator &iter)
Definition: itkImageToNeighborhoodSampleAdaptor.h:179
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ImageType
TImage ImageType
Definition: itkImageToNeighborhoodSampleAdaptor.h:76
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Iterator::operator=
Iterator & operator=(const Iterator &iter)
Definition: itkImageToNeighborhoodSampleAdaptor.h:260
itkConstNeighborhoodIterator.h
itk::Statistics::ImageToNeighborhoodSampleAdaptor::PixelType
typename ImageType::PixelType PixelType
Definition: itkImageToNeighborhoodSampleAdaptor.h:82
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::Statistics::ImageToNeighborhoodSampleAdaptor::End
ConstIterator End() const
Definition: itkImageToNeighborhoodSampleAdaptor.h:309
itkNeighborhoodIterator.h
itk::operator<<
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:216
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::operator==
bool operator==(const ConstIterator &it) const
Definition: itkImageToNeighborhoodSampleAdaptor.h:213
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::ConstIterator
ConstIterator(const ConstIterator &iter)
Definition: itkImageToNeighborhoodSampleAdaptor.h:172
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ValueType
typename MeasurementVectorType::value_type ValueType
Definition: itkImageToNeighborhoodSampleAdaptor.h:100
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Iterator::Iterator
Iterator(const Iterator &iter)
Definition: itkImageToNeighborhoodSampleAdaptor.h:255
itk::Statistics::ImageToNeighborhoodSampleAdaptor::NeighborhoodType
typename NeighborhoodIteratorType::NeighborhoodType NeighborhoodType
Definition: itkImageToNeighborhoodSampleAdaptor.h:92
itk::Statistics::ImageToNeighborhoodSampleAdaptor::SizeType
typename ImageType::SizeType SizeType
Definition: itkImageToNeighborhoodSampleAdaptor.h:86
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Iterator::Iterator
Iterator(Self *adaptor)
Definition: itkImageToNeighborhoodSampleAdaptor.h:251
itk::Statistics::ImageToNeighborhoodSampleAdaptor::PixelContainerConstPointer
typename ImageType::PixelContainerConstPointer PixelContainerConstPointer
Definition: itkImageToNeighborhoodSampleAdaptor.h:83
itk::Statistics::Sample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > >::TotalAbsoluteFrequencyType
NumericTraits< AbsoluteFrequencyType >::AccumulateType TotalAbsoluteFrequencyType
Definition: itkSample.h:87
itk::Neighborhood
A light-weight container object for storing an N-dimensional neighborhood of values.
Definition: itkNeighborhood.h:54
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ImagePointer
typename ImageType::Pointer ImagePointer
Definition: itkImageToNeighborhoodSampleAdaptor.h:77
itk::Statistics::ListSample
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:51
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::GetFrequency
AbsoluteFrequencyType GetFrequency() const
Definition: itkImageToNeighborhoodSampleAdaptor.h:187
itkImage.h
itk::SmartPointer< Self >
itk::Statistics::ImageToNeighborhoodSampleAdaptor::NeighborhoodRadiusType
typename NeighborhoodIteratorType::RadiusType NeighborhoodRadiusType
Definition: itkImageToNeighborhoodSampleAdaptor.h:93
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Statistics::Sample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > >::InstanceIdentifier
typename MeasurementVectorTraits::InstanceIdentifier InstanceIdentifier
Definition: itkSample.h:91
itkImageRegionIteratorWithIndex.h
itk::Statistics::ImageToNeighborhoodSampleAdaptor::NeighborhoodIndexType
typename NeighborhoodIteratorType::IndexType NeighborhoodIndexType
Definition: itkImageToNeighborhoodSampleAdaptor.h:94
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::Statistics::ImageToNeighborhoodSampleAdaptor
This class provides ListSample interface to ITK Image.
Definition: itkImageToNeighborhoodSampleAdaptor.h:55
itk::ConstNeighborhoodIterator::GoToBegin
void GoToBegin()
itk::Statistics::ImageToNeighborhoodSampleAdaptor::OffsetTableType
typename RegionType::OffsetTableType OffsetTableType
Definition: itkImageToNeighborhoodSampleAdaptor.h:85
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::ConstIterator
ConstIterator(const ImageToNeighborhoodSampleAdaptor *adaptor)
Definition: itkImageToNeighborhoodSampleAdaptor.h:170
itkMacro.h
itk::Statistics::ImageToNeighborhoodSampleAdaptor::MeasurementType
ValueType MeasurementType
Definition: itkImageToNeighborhoodSampleAdaptor.h:101
itk::Statistics::ImageToNeighborhoodSampleAdaptor::MeasurementVectorType
typename std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > MeasurementVectorType
Definition: itkImageToNeighborhoodSampleAdaptor.h:99
itk::Statistics::Sample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > >::AbsoluteFrequencyType
MeasurementVectorTraits::AbsoluteFrequencyType AbsoluteFrequencyType
Definition: itkSample.h:84
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Iterator
Iterator.
Definition: itkImageToNeighborhoodSampleAdaptor.h:239
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::m_InstanceIdentifier
InstanceIdentifier m_InstanceIdentifier
Definition: itkImageToNeighborhoodSampleAdaptor.h:231
itk::ConstNeighborhoodIterator::GoToEnd
void GoToEnd()
itk::Statistics::ImageToNeighborhoodSampleAdaptor::End
Iterator End()
Definition: itkImageToNeighborhoodSampleAdaptor.h:286
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:97
itkListSample.h
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::GetInstanceIdentifier
InstanceIdentifier GetInstanceIdentifier() const
Definition: itkImageToNeighborhoodSampleAdaptor.h:199
itk::ImageRegionIteratorWithIndex
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Definition: itkImageRegionIteratorWithIndex.h:73
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::m_MeasurementVectorCache
MeasurementVectorType m_MeasurementVectorCache
Definition: itkImageToNeighborhoodSampleAdaptor.h:230
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::operator++
ConstIterator & operator++()
Definition: itkImageToNeighborhoodSampleAdaptor.h:205
itk::Statistics::ImageToNeighborhoodSampleAdaptor::IndexType
typename ImageType::IndexType IndexType
Definition: itkImageToNeighborhoodSampleAdaptor.h:79
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::GetMeasurementVector
const MeasurementVectorType & GetMeasurementVector() const
Definition: itkImageToNeighborhoodSampleAdaptor.h:193
itk::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
itk::Statistics::ImageToNeighborhoodSampleAdaptor::NeighborhoodSizeType
typename NeighborhoodIteratorType::SizeType NeighborhoodSizeType
Definition: itkImageToNeighborhoodSampleAdaptor.h:95
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator
Const Iterator.
Definition: itkImageToNeighborhoodSampleAdaptor.h:163
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::Statistics::ImageToNeighborhoodSampleAdaptor::OffsetValueType
typename ImageType::OffsetValueType OffsetValueType
Definition: itkImageToNeighborhoodSampleAdaptor.h:81
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Begin
Iterator Begin()
Definition: itkImageToNeighborhoodSampleAdaptor.h:275
itkSmartPointer.h
itk::Statistics::ImageToNeighborhoodSampleAdaptor::RegionType
typename ImageType::RegionType RegionType
Definition: itkImageToNeighborhoodSampleAdaptor.h:84
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ConstIterator::ConstIterator
ConstIterator(NeighborhoodIteratorType iter, InstanceIdentifier iid)
Definition: itkImageToNeighborhoodSampleAdaptor.h:222
itk::Statistics::ImageToNeighborhoodSampleAdaptor::ImageConstPointer
typename ImageType::ConstPointer ImageConstPointer
Definition: itkImageToNeighborhoodSampleAdaptor.h:78
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Begin
ConstIterator Begin() const
Definition: itkImageToNeighborhoodSampleAdaptor.h:298
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::Statistics::ImageToNeighborhoodSampleAdaptor::Iterator::Iterator
Iterator(NeighborhoodIteratorType iter, InstanceIdentifier iid)
Definition: itkImageToNeighborhoodSampleAdaptor.h:268