ITK  6.0.0
Insight Toolkit
itkImageToListSampleAdaptor.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 itkImageToListSampleAdaptor_h
19#define itkImageToListSampleAdaptor_h
20
21#include <typeinfo>
22#include <utility>
23
24#include "itkImage.h"
25#include "itkPixelTraits.h"
26#include "itkListSample.h"
27#include "itkSmartPointer.h"
30
31namespace itk
32{
33namespace Statistics
34{
53template <typename TImage>
54class ITK_TEMPLATE_EXPORT ImageToListSampleAdaptor
55 : public ListSample<typename MeasurementVectorPixelTraits<typename TImage::PixelType>::MeasurementVectorType>
56{
57public:
58 ITK_DISALLOW_COPY_AND_MOVE(ImageToListSampleAdaptor);
59
63 using Superclass =
68
70 itkOverrideGetNameOfClassMacro(ImageToListSampleAdaptor);
71
73 itkNewMacro(Self);
74
76 using ImageType = TImage;
80 using PixelType = typename ImageType::PixelType;
81 using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
82
87
88
96
97 using typename Superclass::AbsoluteFrequencyType;
98 using typename Superclass::TotalAbsoluteFrequencyType;
99 using typename Superclass::MeasurementVectorSizeType;
100 using typename Superclass::InstanceIdentifier;
103
105 void
106 SetImage(const TImage * image);
107
109 const TImage *
110 GetImage() const;
111
114 Size() const override;
115
118 GetMeasurementVector(InstanceIdentifier id) const override;
119
121 GetMeasurementVectorSize() const override
122 {
123 // some filter are expected that this method returns something even if the
124 // input is not set. This won't be the right value for a variable length vector
125 // but it's better than an exception.
126 if (m_Image.IsNull())
127 {
128 return Superclass::GetMeasurementVectorSize();
129 }
130
131 return m_Image->GetNumberOfComponentsPerPixel();
132 }
133
136 GetFrequency(InstanceIdentifier id) const override;
137
140 GetTotalFrequency() const override;
141
146 class ConstIterator
148 friend class ImageToListSampleAdaptor;
149
150 public:
151 ConstIterator() = delete;
153 ConstIterator(const ImageToListSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
155 ConstIterator(const ConstIterator & iter)
156 : m_Iter(iter.m_Iter)
157 , m_InstanceIdentifier(iter.m_InstanceIdentifier)
158 {}
159
161 operator=(const ConstIterator & iter)
162 {
163 m_Iter = iter.m_Iter;
164 m_InstanceIdentifier = iter.m_InstanceIdentifier;
165 return *this;
166 }
167
169 GetFrequency() const
170 {
171 return 1;
172 }
173
175 GetMeasurementVector() const
176 {
177 MeasurementVectorTraits::Assign(this->m_MeasurementVectorCache, m_Iter.Get());
178 return this->m_MeasurementVectorCache;
179 }
180
182 GetInstanceIdentifier() const
183 {
184 return m_InstanceIdentifier;
185 }
186
188 operator++()
189 {
190 ++m_Iter;
191 ++m_InstanceIdentifier;
192 return *this;
193 }
194
195 bool
196 operator==(const ConstIterator & it) const
197 {
198 return (m_Iter == it.m_Iter);
199 }
201 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION(ConstIterator);
202
203 protected:
204 // This method should only be available to the ListSample class
206 : m_Iter(std::move(iter))
207 , m_InstanceIdentifier(iid)
208 {}
209
210 private:
212 mutable MeasurementVectorType m_MeasurementVectorCache;
213 InstanceIdentifier m_InstanceIdentifier;
214 };
215
220 class Iterator : public ConstIterator
222 friend class ImageToListSampleAdaptor;
223
224 public:
225 Iterator() = delete;
226 Iterator(const Self * adaptor) = delete;
228 Iterator(const ConstIterator & it) = delete;
230 operator=(const ConstIterator & it) = delete;
232 Iterator(Self * adaptor)
233 : ConstIterator(adaptor)
234 {}
236 Iterator(const Iterator & iter)
237 : ConstIterator(iter)
238 {}
239
241 operator=(const Iterator & iter)
242 {
243 this->ConstIterator::operator=(iter);
244 return *this;
245 }
246
247 protected:
249 : ConstIterator(iter, iid)
250 {}
251 };
252
255 Begin()
256 {
257 const ImagePointer nonConstImage = const_cast<ImageType *>(m_Image.GetPointer());
258 ImageIteratorType imageIterator(nonConstImage, nonConstImage->GetLargestPossibleRegion());
259 imageIterator.GoToBegin();
260 const Iterator iter(imageIterator, 0);
261 return iter;
262 }
267 End()
268 {
269 const ImagePointer nonConstImage = const_cast<ImageType *>(m_Image.GetPointer());
270 const typename ImageType::RegionType & largestRegion = nonConstImage->GetLargestPossibleRegion();
271 ImageIteratorType imageIterator(nonConstImage, largestRegion);
272 imageIterator.GoToEnd();
273 const Iterator iter(imageIterator, largestRegion.GetNumberOfPixels());
276 return iter;
277 }
278
281 Begin() const
282 {
283 ImageConstIteratorType imageConstIterator(m_Image, m_Image->GetLargestPossibleRegion());
284 imageConstIterator.GoToBegin();
285 const ConstIterator iter(imageConstIterator, 0);
288 return iter;
289 }
290
293 End() const
294 {
295 const typename ImageType::RegionType & largestRegion = m_Image->GetLargestPossibleRegion();
296 ImageConstIteratorType imageConstIterator(m_Image, largestRegion);
297 imageConstIterator.GoToEnd();
298 const ConstIterator iter(imageConstIterator, largestRegion.GetNumberOfPixels());
301 return iter;
302 }
303
304protected:
306 ~ImageToListSampleAdaptor() override = default;
307 void
308 PrintSelf(std::ostream & os, Indent indent) const override;
309
310private:
312 mutable MeasurementVectorType m_MeasurementVectorInternal{};
313
314}; // end of class ImageToListSampleAdaptor
315} // end of namespace Statistics
316} // end of namespace itk
317
318#ifndef ITK_MANUAL_INSTANTIATION
319# include "itkImageToListSampleAdaptor.hxx"
320#endif
321
322#endif
Base class for all data objects in ITK.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
Traits for a pixel that define the dimension and component type.
Iterator(const ImageConstIteratorType &iter, InstanceIdentifier iid)=delete
This class provides ListSample interface to ITK Image.
typename MeasurementVectorTraitsType::ValueType MeasurementType
typename MeasurementPixelTraitsType::MeasurementVectorType MeasurementVectorType
typename ImageType::ConstPointer ImageConstPointer
typename ImageType::PixelContainerConstPointer PixelContainerConstPointer
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:52
typename TMeasurementVector::ValueType ValueType
static void Assign(TArrayType &m, const TArrayType &v)
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
bool operator==(const Index< VDimension > &one, const Index< VDimension > &two)
Definition: itkIndex.h:552
STL namespace.
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:70