ITK  6.0.0
Insight Toolkit
itkJointDomainImageToListSampleAdaptor.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 itkJointDomainImageToListSampleAdaptor_h
19#define itkJointDomainImageToListSampleAdaptor_h
20
21#include "itkPoint.h"
22#include "itkPixelTraits.h"
25#include "itkListSample.h"
26
27namespace itk
28{
29namespace Statistics
30{
39template <typename TImage>
41{
45
46 static constexpr unsigned int ImageDimension = TImage::ImageDimension;
47 static constexpr unsigned int Dimension = TImage::ImageDimension + PixelTraitsType::Dimension;
48
49 using CoordinateRepType = float;
53
55}; // end of ImageJointDomainTraits
56
87template <typename TImage>
88class ITK_TEMPLATE_EXPORT JointDomainImageToListSampleAdaptor
89 : public ListSample<typename ImageJointDomainTraits<TImage>::MeasurementVectorType>
90{
91public:
92 ITK_DISALLOW_COPY_AND_MOVE(JointDomainImageToListSampleAdaptor);
93
96
98
101
103
109
111 itkOverrideGetNameOfClassMacro(JointDomainImageToListSampleAdaptor);
112
114 itkNewMacro(Self);
115
117 static constexpr unsigned int MeasurementVectorSize = ImageJointDomainTraitsType::Dimension;
118
119 using typename Superclass::MeasurementVectorSizeType;
120
123 using typename Superclass::AbsoluteFrequencyType;
124 using typename Superclass::TotalAbsoluteFrequencyType;
125 using typename Superclass::InstanceIdentifier;
126
128 using ImageType = TImage;
131
134 using PixelType = typename ImageType::PixelType;
135 using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
140
142 void
143 SetImage(const TImage * image);
144
146 const TImage *
147 GetImage() const;
148
151 Size() const override;
152
156
159 GetTotalFrequency() const override;
160
161 static constexpr unsigned int RangeDomainDimension = itk::PixelTraits<typename TImage::PixelType>::Dimension;
162
164
165 using InstanceIdentifierVectorType = std::vector<InstanceIdentifier>;
167
169 void
171
176
178 itkSetMacro(UsePixelContainer, bool);
179 itkGetConstMacro(UsePixelContainer, bool);
180 itkBooleanMacro(UsePixelContainer);
183 // void PrintSelf(std::ostream& os, Indent indent) const override;
184
190 {
192
193 public:
194 ConstIterator() = delete;
195
196 ConstIterator(const JointDomainImageToListSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
197
199 {
200 m_InstanceIdentifier = iter.m_InstanceIdentifier;
201 m_Adaptor = iter.m_Adaptor;
202 }
203
206 {
207 m_InstanceIdentifier = iter.m_InstanceIdentifier;
208 return *this;
209 }
210
213 {
214 return 1;
215 }
216
219 {
220 m_MeasurementVectorCache = m_Adaptor->GetMeasurementVector(m_InstanceIdentifier);
221 return this->m_MeasurementVectorCache;
222 }
223
226 {
227 return m_InstanceIdentifier;
228 }
229
232 {
233 ++m_InstanceIdentifier;
234 return *this;
235 }
236
237 bool
238 operator==(const ConstIterator & it) const
239 {
240 return (m_InstanceIdentifier == it.m_InstanceIdentifier);
241 }
242
244
245 protected:
246 // This method should only be available to the ListSample class
248 {
249 m_Adaptor = adaptor;
250 m_InstanceIdentifier = iid;
251 }
252
253 private:
257 };
258
263 class Iterator : public ConstIterator
264 {
266
267 public:
268 Iterator() = delete;
269 Iterator(const Self * adaptor) = delete;
270 Iterator(const ConstIterator & it) = delete;
272 operator=(const ConstIterator & it) = delete;
273
274 Iterator(Self * adaptor)
275 : ConstIterator(adaptor)
276 {}
277
278 Iterator(const Iterator & iter)
279 : ConstIterator(iter)
280 {}
281
282 Iterator &
283 operator=(const Iterator & iter)
284 {
285 this->ConstIterator::operator=(iter);
286 return *this;
287 }
288
289 protected:
291 : ConstIterator(adaptor, iid)
292 {}
293 };
294
296 Iterator
298 {
299 Iterator iter(this, 0);
300
301 return iter;
302 }
303
305 Iterator
307 {
308 Iterator iter(this, m_Image->GetPixelContainer()->Size());
309
310 return iter;
311 }
312
314 ConstIterator
315 Begin() const
316 {
317 ConstIterator iter(this, 0);
318
319 return iter;
320 }
321
323 ConstIterator
324 End() const
325 {
326 ConstIterator iter(this, m_Image->GetPixelContainer()->Size());
327
328 return iter;
329 }
330
331protected:
334 void
335 PrintSelf(std::ostream & os, Indent indent) const override;
336
337private:
338 NormalizationFactorsType m_NormalizationFactors{};
339 mutable MeasurementVectorType m_TempVector{};
340 mutable PointType m_TempPoint{};
341 mutable ImageIndexType m_TempIndex{};
342 mutable RangeDomainMeasurementVectorType m_TempRangeVector{};
344 bool m_UsePixelContainer{};
345
346 PixelContainerConstPointer m_PixelContainer{};
347}; // end of class JointDomainImageToListSampleAdaptor
348} // end of namespace Statistics
349} // end of namespace itk
350
351#ifndef ITK_MANUAL_INSTANTIATION
352# include "itkJointDomainImageToListSampleAdaptor.hxx"
353#endif
354
355#endif
Base class for all data objects in ITK.
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
A multi-dimensional iterator templated over image type that walks a region of pixels.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Trait to determine what datatype is needed if the specified pixel types are "joined" into a single ve...
Base class for most ITK classes.
Definition: itkObject.h:62
Traits for a pixel that define the dimension and component type.
typename TPixelType::ValueType ValueType
static constexpr unsigned int Dimension
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:54
ConstIterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
Iterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
ConstIterator & operator=(const ConstIterator &it)=delete
This adaptor returns measurement vectors composed of an image pixel's range domain value (pixel value...
typename ImageType::PixelContainerConstPointer PixelContainerConstPointer
typename ImageJointDomainTraitsType::RangeDomainMeasurementType RangeDomainMeasurementType
TotalAbsoluteFrequencyType GetTotalFrequency() const override
InstanceIdentifier Size() const override
void PrintSelf(std::ostream &os, Indent indent) const override
const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const override
typename ImageJointDomainTraitsType::MeasurementType MeasurementType
typename ImageJointDomainTraitsType::CoordinateRepType CoordinateRepType
AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const override
void SetNormalizationFactors(NormalizationFactorsType &factors)
typename ImageJointDomainTraitsType::MeasurementVectorType MeasurementVectorType
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:52
typename MeasurementVectorTraits::InstanceIdentifier InstanceIdentifier
Definition: itkSample.h:91
NumericTraits< AbsoluteFrequencyType >::AccumulateType TotalAbsoluteFrequencyType
Definition: itkSample.h:87
MeasurementVectorTraits::AbsoluteFrequencyType AbsoluteFrequencyType
Definition: itkSample.h:84
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
This class provides the type definition for the measurement vector in the joint domain (range domain ...
typename PixelTraitsType::ValueType RangeDomainMeasurementType