ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkImageToNeighborhoodSampleAdaptor.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 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 namespace Statistics {
35 
51  template < typename TImage, typename TBoundaryCondition >
52 class ITK_TEMPLATE_EXPORT ImageToNeighborhoodSampleAdaptor :
53  public ListSample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > >
54 {
55 public:
58 
61 
64 
67 
69  itkNewMacro(Self);
70 
72  typedef TImage ImageType;
73  typedef typename ImageType::Pointer ImagePointer;
74  typedef typename ImageType::ConstPointer ImageConstPointer;
75  typedef typename ImageType::IndexType IndexType;
76  typedef typename ImageType::OffsetType OffsetType;
78  typedef typename ImageType::PixelType PixelType;
79  typedef typename ImageType::PixelContainerConstPointer PixelContainerConstPointer;
80  typedef typename ImageType::RegionType RegionType;
81  typedef typename RegionType::OffsetTableType OffsetTableType;
82  typedef typename ImageType::SizeType SizeType;
84 
92 
95  typedef typename std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > >
97  typedef typename MeasurementVectorType::value_type ValueType;
99 
100  typedef typename Superclass::AbsoluteFrequencyType AbsoluteFrequencyType;
101  typedef typename Superclass::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType;
102  typedef typename Superclass::MeasurementVectorSizeType MeasurementVectorSizeType;
103  typedef typename Superclass::InstanceIdentifier InstanceIdentifier;
104 
106  void SetImage(const TImage* image);
107 
109  const TImage* GetImage() const;
110 
112  void SetRadius(const NeighborhoodRadiusType& radius);
113 
115  NeighborhoodRadiusType GetRadius() const;
116 
118  void SetRegion(const RegionType& region);
119 
121  RegionType GetRegion() const;
122 
123  void SetUseImageRegion(const bool& flag);
124 
126  itkGetConstMacro( UseImageRegion, bool );
127 
129  itkBooleanMacro( UseImageRegion );
130 
131 
133  InstanceIdentifier Size() const ITK_OVERRIDE;
134 
136  virtual const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const ITK_OVERRIDE;
137 
139  AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const ITK_OVERRIDE;
140 
142  TotalAbsoluteFrequencyType GetTotalFrequency() const ITK_OVERRIDE;
143 
149  {
151  public:
152 
154  {
155  *this = adaptor->Begin();
156  }
157 
158  ConstIterator(const ConstIterator &iter)
159  {
160  m_MeasurementVectorCache = iter.m_MeasurementVectorCache;
161  m_InstanceIdentifier = iter.m_InstanceIdentifier;
162  }
163 
164  ConstIterator& operator=( const ConstIterator & iter )
165  {
166  m_MeasurementVectorCache = iter.m_MeasurementVectorCache;
167  m_InstanceIdentifier = iter.m_InstanceIdentifier;
168  return *this;
169  }
170 
172  {
173  return 1;
174  }
175 
177  {
178  return this->m_MeasurementVectorCache;
179  }
180 
182  {
183  return m_InstanceIdentifier;
184  }
185 
186  ConstIterator& operator++()
187  {
188  ++(m_MeasurementVectorCache[0]);
189  ++m_InstanceIdentifier;
190  return *this;
191  }
192 
193  bool operator!=(const ConstIterator &it)
194  {
195  return (m_MeasurementVectorCache[0] != it.m_MeasurementVectorCache[0]);
196  }
197 
198  bool operator==(const ConstIterator &it)
199  {
200  return (m_MeasurementVectorCache[0] == it.m_MeasurementVectorCache[0]);
201  }
202 
203  protected:
204  // This method should only be available to the ListSample class
207  InstanceIdentifier iid)
208  {
209  this->m_MeasurementVectorCache.clear();
210  this->m_MeasurementVectorCache.push_back(iter);
211  m_InstanceIdentifier = iid;
212  }
213 
214  private:
215  ConstIterator() ITK_DELETED_FUNCTION;
216  mutable MeasurementVectorType m_MeasurementVectorCache;
217  InstanceIdentifier m_InstanceIdentifier;
218  };
219 
224  class Iterator : public ConstIterator
225  {
226 
228 
229  public:
230 
231  Iterator(Self * adaptor):ConstIterator(adaptor)
232  {
233  }
234 
235  Iterator(const Iterator &iter):ConstIterator( iter )
236  {
237  }
238 
239  Iterator& operator =(const Iterator & iter)
240  {
241  this->ConstIterator::operator=( iter );
242  return *this;
243  }
244 
245 #if !(defined(_MSC_VER) && (_MSC_VER <= 1200))
246  protected:
247 #endif
248  //This copy constructor is actually used in Iterator Begin()!
249  Iterator(NeighborhoodIteratorType iter, InstanceIdentifier iid):ConstIterator( iter, iid )
250  {
251  }
252 
253  private:
254  // To ensure const-correctness these method must not be in the public API.
255  // The are not implemented, since they should never be called.
256  Iterator() ITK_DELETED_FUNCTION;
257  Iterator(const Self * adaptor) ITK_DELETED_FUNCTION;
258  Iterator(const ConstIterator & it) ITK_DELETED_FUNCTION;
259  ConstIterator& operator=(const ConstIterator& it) ITK_DELETED_FUNCTION;
260  };
261 
263  Iterator Begin()
264  {
265  NeighborhoodIteratorType nIterator( m_Radius, m_Image, m_Region);
266  nIterator.GoToBegin();
267  Iterator iter(nIterator, 0);
268  return iter;
269  }
271 
274  {
275  NeighborhoodIteratorType nIterator( m_Radius, m_Image, m_Region);
276  nIterator.GoToEnd();
277  Iterator iter(nIterator, m_Region.GetNumberOfPixels());
278  return iter;
279  }
281 
282 
284  ConstIterator Begin() const
285  {
286  NeighborhoodIteratorType nIterator( m_Radius, m_Image, m_Region);
287  nIterator.GoToBegin();
288  ConstIterator iter(nIterator, 0);
289  return iter;
290  }
292 
294  ConstIterator End() const
295  {
296  NeighborhoodIteratorType nIterator( m_Radius, m_Image, m_Region);
297  nIterator.GoToEnd();
298  ConstIterator iter(nIterator, m_Region.GetNumberOfPixels());
299  return iter;
300  }
302 
303 protected:
305  virtual ~ImageToNeighborhoodSampleAdaptor() ITK_OVERRIDE {}
306  void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
307 
308 private:
309  ITK_DISALLOW_COPY_AND_ASSIGN(ImageToNeighborhoodSampleAdaptor);
310 
319 
320 }; // end of class ImageToNeighborhoodSampleAdaptor
321 
322 } // end of namespace Statistics
323 
324 template <typename TImage, typename TBoundaryCondition>
325  std::ostream & operator<<(std::ostream &os,
327 
328 } // end of namespace itk
329 
330 #ifndef ITK_MANUAL_INSTANTIATION
331 #include "itkImageToNeighborhoodSampleAdaptor.hxx"
332 #endif
333 
334 #endif
ConstIterator(NeighborhoodIteratorType iter, InstanceIdentifier iid)
ConstNeighborhoodIterator< TImage, TBoundaryCondition > NeighborhoodIteratorType
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
signed long OffsetValueType
Definition: itkIntTypes.h:154
A light-weight container object for storing an N-dimensional neighborhood of values.
std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > MeasurementVectorType
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:192
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
ListSample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > > Superclass
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:51
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Iterator(NeighborhoodIteratorType iter, InstanceIdentifier iid)
This class provides ListSample interface to ITK Image.
NeighborhoodIterator< TImage, TBoundaryCondition > NonConstNeighborhoodIteratorType
Base class for all data objects in ITK.
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.