ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkJointDomainImageToListSampleAdaptor.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 itkJointDomainImageToListSampleAdaptor_h
19 #define itkJointDomainImageToListSampleAdaptor_h
20 
21 #include "itkPoint.h"
22 #include "itkPixelTraits.h"
24 #include "itkImageRegionIterator.h"
25 #include "itkListSample.h"
26 
27 namespace itk
28 {
29 namespace Statistics
30 {
39 template< typename TImage >
44 
45  itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
46  itkStaticConstMacro(Dimension,
47  unsigned int,
48  TImage::ImageDimension
50 
51  typedef float CoordinateRepType;
55 
58 }; // end of ImageJointDomainTraits
59 
90 template< typename TImage >
92  public ListSample< typename ImageJointDomainTraits< TImage >::MeasurementVectorType >
93 {
94 public:
95 
98 
99  typedef ListSample<
101 
104 
106 
112 
115 
117  itkNewMacro(Self);
118 
120  itkStaticConstMacro(MeasurementVectorSize, unsigned int,
122 
124 
130 
132  typedef TImage ImageType;
135 
136  typedef typename ImageType::Pointer ImagePointer;
137  typedef typename ImageType::ConstPointer ImageConstPointer;
138  typedef typename ImageType::PixelType PixelType;
139  typedef typename ImageType::PixelContainerConstPointer PixelContainerConstPointer;
140  typedef typename ImageType::IndexType ImageIndexType;
141  typedef typename ImageType::SizeType ImageSizeType;
142  typedef typename ImageType::RegionType ImageRegionType;
144 
146  void SetImage(const TImage *image);
147 
149  const TImage * GetImage() const;
150 
152  InstanceIdentifier Size() const ITK_OVERRIDE;
153 
156 
158  TotalAbsoluteFrequencyType GetTotalFrequency() const ITK_OVERRIDE;
159 
160  itkStaticConstMacro(RangeDomainDimension, unsigned int,
161  itk::PixelTraits< typename TImage::PixelType >::Dimension);
162 
164  itkGetStaticConstMacro(RangeDomainDimension) > RangeDomainMeasurementVectorType;
165 
167  typedef FixedArray< float, itkGetStaticConstMacro(MeasurementVectorSize) > NormalizationFactorsType;
168 
170  void SetNormalizationFactors(NormalizationFactorsType & factors);
171 
174  const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const ITK_OVERRIDE;
175 
177  itkSetMacro(UsePixelContainer, bool);
178  itkGetConstMacro(UsePixelContainer, bool);
179  itkBooleanMacro(UsePixelContainer);
181 
182  // void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
183 
189  {
191 
192 public:
193 
195  {
196  *this = adaptor->Begin();
197  }
198 
199  ConstIterator(const ConstIterator & iter)
200  {
201  m_InstanceIdentifier = iter.m_InstanceIdentifier;
202  m_Adaptor = iter.m_Adaptor;
203  }
204 
205  ConstIterator & operator=(const ConstIterator & iter)
206  {
207  m_InstanceIdentifier = iter.m_InstanceIdentifier;
208  return *this;
209  }
210 
212  {
213  return 1;
214  }
215 
217  {
218  m_MeasurementVectorCache = m_Adaptor->GetMeasurementVector(m_InstanceIdentifier);
219  return this->m_MeasurementVectorCache;
220  }
221 
223  {
224  return m_InstanceIdentifier;
225  }
226 
227  ConstIterator & operator++()
228  {
229  ++m_InstanceIdentifier;
230  return *this;
231  }
232 
233  bool operator!=(const ConstIterator & it)
234  {
235  return ( m_InstanceIdentifier != it.m_InstanceIdentifier );
236  }
237 
238  bool operator==(const ConstIterator & it)
239  {
240  return ( m_InstanceIdentifier == it.m_InstanceIdentifier );
241  }
242 
243 protected:
244  // This method should only be available to the ListSample class
247  InstanceIdentifier iid)
248  {
249  m_Adaptor = adaptor;
250  m_InstanceIdentifier = iid;
251  }
252 
253  // This method is purposely not implemented
254  ConstIterator();
255 
256 private:
260  };
261 
266  class Iterator:public ConstIterator
267  {
269 
270 public:
271 
272  Iterator(Self *adaptor):ConstIterator(adaptor)
273  {}
274 
275  Iterator(const Iterator & iter):ConstIterator(iter)
276  {}
277 
278  Iterator & operator=(const Iterator & iter)
279  {
280  this->ConstIterator::operator=(iter);
281  return *this;
282  }
283 
284 protected:
285  // To ensure const-correctness these method must not be in the public API.
286  // The are purposly not implemented, since they should never be called.
287  Iterator();
288  Iterator(const Self *adaptor);
289  Iterator(const ConstIterator & it);
290  ConstIterator & operator=(const ConstIterator & it);
291 
294  InstanceIdentifier iid):ConstIterator(adaptor, iid)
295  {}
296 
297 private:
298  };
299 
302  {
303  Iterator iter(this, 0);
304 
305  return iter;
306  }
307 
310  {
311  Iterator iter( this, m_Image->GetPixelContainer()->Size() );
312 
313  return iter;
314  }
315 
318  {
319  ConstIterator iter(this, 0);
320 
321  return iter;
322  }
323 
326  {
327  ConstIterator iter( this, m_Image->GetPixelContainer()->Size() );
328 
329  return iter;
330  }
331 
332 protected:
335  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
336 
337 private:
338  JointDomainImageToListSampleAdaptor(const Self &); //purposely not implemented
339  void operator=(const Self &); //purposely not implemented
340 
348 
350 }; // end of class JointDomainImageToListSampleAdaptor
351 } // end of namespace Statistics
352 } // end of namespace itk
353 
354 #ifndef ITK_MANUAL_INSTANTIATION
355 #include "itkJointDomainImageToListSampleAdaptor.hxx"
356 #endif
357 
358 #endif
AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const override
Point< CoordinateRepType, itkGetStaticConstMacro(ImageDimension) > PointType
TPixelType::ValueType ValueType
JoinTraits< RangeDomainMeasurementType, CoordinateRepType > JoinTraitsType
Superclass::TotalAbsoluteFrequencyType TotalAbsoluteFrequencyType
Definition: itkListSample.h:71
ImageJointDomainTraitsType::MeasurementVectorType MeasurementVectorType
InstanceIdentifier Size() const override
ImageJointDomainTraitsType::RangeDomainMeasurementType RangeDomainMeasurementType
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
Traits for a pixel that define the dimension and component type.
void SetNormalizationFactors(NormalizationFactorsType &factors)
Superclass::AbsoluteFrequencyType AbsoluteFrequencyType
Definition: itkListSample.h:70
FixedArray< MeasurementType, itkGetStaticConstMacro(Dimension) > MeasurementVectorType
Trait to determine what datatype is needed if the specified pixel types are &quot;joined&quot; into a single ve...
Superclass::InstanceIdentifier InstanceIdentifier
Definition: itkListSample.h:72
Iterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
Superclass::MeasurementVectorSizeType MeasurementVectorSizeType
Definition: itkListSample.h:68
ListSample< typename ImageJointDomainTraits< TImage >::MeasurementVectorType > Superclass
void PrintSelf(std::ostream &os, Indent indent) const override
This adaptor returns measurement vectors composed of an image pixel&#39;s range domain value (pixel value...
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
PixelTraits< typename TImage::PixelType > PixelTraitsType
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:51
static const unsigned int Dimension
This class provides the type definition for the measurement vector in the joint domain (range domain ...
Base class for all data objects in ITK.
A multi-dimensional iterator templated over image type that walks a region of pixels.
ConstIterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const override
TotalAbsoluteFrequencyType GetTotalFrequency() const override