ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkImageToListSampleAdaptor.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 itkImageToListSampleAdaptor_h
19 #define itkImageToListSampleAdaptor_h
20 
21 #include <typeinfo>
22 
23 #include "itkImage.h"
24 #include "itkPixelTraits.h"
25 #include "itkListSample.h"
26 #include "itkSmartPointer.h"
27 #include "itkImageRegionIterator.h"
29 
30 namespace itk
31 {
32 namespace Statistics
33 {
52 template< typename TImage >
53 class ITK_TEMPLATE_EXPORT ImageToListSampleAdaptor:
54  public ListSample< typename MeasurementVectorPixelTraits< typename TImage::PixelType >::MeasurementVectorType >
55 {
56 public:
59 
60  typedef ListSample< typename MeasurementVectorPixelTraits<
61  typename TImage::PixelType >::MeasurementVectorType >
63 
66 
68  itkTypeMacro(ImageToListSampleAdaptor, ListSample);
69 
71  itkNewMacro(Self);
72 
74  typedef TImage ImageType;
75  typedef typename ImageType::Pointer ImagePointer;
76  typedef typename ImageType::ConstPointer ImageConstPointer;
77  typedef typename ImageType::IndexType IndexType;
78  typedef typename ImageType::PixelType PixelType;
79  typedef typename ImageType::PixelContainerConstPointer PixelContainerConstPointer;
80 
85 
86 
91 
94 
95  typedef typename Superclass::AbsoluteFrequencyType AbsoluteFrequencyType;
96  typedef typename Superclass::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType;
97  typedef typename Superclass::MeasurementVectorSizeType MeasurementVectorSizeType;
98  typedef typename Superclass::InstanceIdentifier InstanceIdentifier;
99 
101 
103  void SetImage(const TImage *image);
104 
106  const TImage * GetImage() const;
107 
109  InstanceIdentifier Size() const ITK_OVERRIDE;
110 
112  virtual const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const ITK_OVERRIDE;
113 
114  virtual MeasurementVectorSizeType GetMeasurementVectorSize() const ITK_OVERRIDE
115  {
116  // some filter are expected that this method returns something even if the
117  // input is not set. This won't be the right value for a variable length vector
118  // but it's better than an exception.
119  if( m_Image.IsNull() )
120  {
121  return Superclass::GetMeasurementVectorSize();
122  }
123  else
124  {
125  return m_Image->GetNumberOfComponentsPerPixel();
126  }
127  }
128 
130  AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const ITK_OVERRIDE;
131 
133  TotalAbsoluteFrequencyType GetTotalFrequency() const ITK_OVERRIDE;
134 
140  {
142 
143  public:
144 
146  {
147  *this = adaptor->Begin();
148  }
149 
150  ConstIterator(const ConstIterator & iter) :
151  m_Iter(iter.m_Iter),
152  m_InstanceIdentifier(iter.m_InstanceIdentifier)
153  {}
154 
155  ConstIterator & operator=(const ConstIterator & iter)
156  {
157  m_Iter = iter.m_Iter;
158  m_InstanceIdentifier = iter.m_InstanceIdentifier;
159  return *this;
160  }
161 
163  {
164  return 1;
165  }
166 
168  {
169  MeasurementVectorTraits::Assign( this->m_MeasurementVectorCache, m_Iter.Get() );
170  return this->m_MeasurementVectorCache;
171  }
172 
174  {
175  return m_InstanceIdentifier;
176  }
177 
178  ConstIterator & operator++()
179  {
180  ++m_Iter;
181  ++m_InstanceIdentifier;
182  return *this;
183  }
184 
185  bool operator!=(const ConstIterator & it)
186  {
187  return ( m_Iter != it.m_Iter );
188  }
189 
190  bool operator==(const ConstIterator & it)
191  {
192  return ( m_Iter == it.m_Iter );
193  }
194 
195  protected:
196  // This method should only be available to the ListSample class
198  m_Iter(iter),
199  m_InstanceIdentifier(iid)
200  {}
201 
202  private:
203  ConstIterator() ITK_DELETED_FUNCTION;
205  mutable MeasurementVectorType m_MeasurementVectorCache;
206  InstanceIdentifier m_InstanceIdentifier;
207  };
208 
213  class Iterator:
214  public ConstIterator
215  {
217 
218  public:
219 
220  Iterator(Self *adaptor) :
221  ConstIterator(adaptor)
222  {}
223 
224  Iterator(const Iterator & iter):
225  ConstIterator(iter)
226  {}
227 
228  Iterator & operator=(const Iterator & iter)
229  {
230  this->ConstIterator::operator=(iter);
231  return *this;
232  }
233 
234  protected:
236  ConstIterator(iter, iid)
237  {}
238 
239  private:
240  // To ensure const-correctness these method must not be in the public API.
241  // The are purposly not implemented, since they should never be called.
242  Iterator() ITK_DELETED_FUNCTION;
243  Iterator(const Self *adaptor) ITK_DELETED_FUNCTION;
244  Iterator(const ImageConstIteratorType & iter, InstanceIdentifier iid) ITK_DELETED_FUNCTION;
245  Iterator(const ConstIterator & it) ITK_DELETED_FUNCTION;
246  ConstIterator & operator=(const ConstIterator & it) ITK_DELETED_FUNCTION;
247  };
248 
250  Iterator Begin()
251  {
252  ImagePointer nonConstImage = const_cast< ImageType * >( m_Image.GetPointer() );
253  ImageIteratorType imageIterator( nonConstImage, nonConstImage->GetLargestPossibleRegion() );
254  imageIterator.GoToBegin();
255  Iterator iter(imageIterator, 0);
256  return iter;
257  }
259 
262  {
263  ImagePointer nonConstImage = const_cast< ImageType * >( m_Image.GetPointer() );
264  const typename ImageType::RegionType & largestRegion = nonConstImage->GetLargestPossibleRegion();
265  ImageIteratorType imageIterator( nonConstImage, largestRegion );
266  imageIterator.GoToEnd();
267  Iterator iter( imageIterator, largestRegion.GetNumberOfPixels() );
269 
270  return iter;
271  }
272 
274  ConstIterator Begin() const
275  {
276  ImageConstIteratorType imageConstIterator( m_Image, m_Image->GetLargestPossibleRegion() );
277  imageConstIterator.GoToBegin();
278  ConstIterator iter(imageConstIterator, 0);
280 
281  return iter;
282  }
283 
285  ConstIterator End() const
286  {
287  const typename ImageType::RegionType & largestRegion = m_Image->GetLargestPossibleRegion();
288  ImageConstIteratorType imageConstIterator( m_Image, largestRegion );
289  imageConstIterator.GoToEnd();
290  ConstIterator iter( imageConstIterator, largestRegion.GetNumberOfPixels() );
292 
293  return iter;
294  }
295 
296 protected:
298  virtual ~ImageToListSampleAdaptor() ITK_OVERRIDE {}
299  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
300 
301 private:
302  ITK_DISALLOW_COPY_AND_ASSIGN(ImageToListSampleAdaptor);
303 
306 
307 }; // end of class ImageToListSampleAdaptor
308 } // end of namespace Statistics
309 } // end of namespace itk
310 
311 #ifndef ITK_MANUAL_INSTANTIATION
312 #include "itkImageToListSampleAdaptor.hxx"
313 #endif
314 
315 #endif
MeasurementVectorTraitsType::ValueType MeasurementType
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
ImageRegionConstIterator< ImageType > ImageConstIteratorType
ImageRegionIterator< ImageType > ImageIteratorType
static void Assign(TArrayType &m, const TArrayType &v)
Traits for a pixel that define the dimension and component type.
MeasurementVectorTraitsTypes< MeasurementVectorType > MeasurementVectorTraitsType
Superclass::AbsoluteFrequencyType AbsoluteFrequencyType
MeasurementPixelTraitsType::MeasurementVectorType MeasurementVectorType
Iterator(const ImageIteratorType &iter, InstanceIdentifier iid)
Superclass::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType
Superclass::MeasurementVectorSizeType MeasurementVectorSizeType
ListSample< typename MeasurementVectorPixelTraits< typename TImage::PixelType >::MeasurementVectorType > 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
This class provides ListSample interface to ITK Image.
PixelTraits< typename TImage::PixelType > PixelTraitsType
Base class for all data objects in ITK.
MeasurementVectorPixelTraits< PixelType > MeasurementPixelTraitsType
A multi-dimensional iterator templated over image type that walks a region of pixels.
ConstIterator(const ImageConstIteratorType &iter, InstanceIdentifier iid)
ImageType::PixelContainerConstPointer PixelContainerConstPointer