ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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"
24#include "itkListSample.h"
25
26namespace itk::Statistics
27{
36template <typename TImage>
53
83
84template <typename TImage>
85class ITK_TEMPLATE_EXPORT JointDomainImageToListSampleAdaptor
86 : public ListSample<typename ImageJointDomainTraits<TImage>::MeasurementVectorType>
87{
88public:
89 ITK_DISALLOW_COPY_AND_MOVE(JointDomainImageToListSampleAdaptor);
90
93
95
98
100
106
108 itkOverrideGetNameOfClassMacro(JointDomainImageToListSampleAdaptor);
109
111 itkNewMacro(Self);
112
115
117
122 using typename Superclass::InstanceIdentifier;
123
125 using ImageType = TImage;
128
129 using ImagePointer = typename ImageType::Pointer;
130 using ImageConstPointer = typename ImageType::ConstPointer;
131 using PixelType = typename ImageType::PixelType;
132 using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
133 using ImageIndexType = typename ImageType::IndexType;
134 using ImageSizeType = typename ImageType::SizeType;
135 using ImageRegionType = typename ImageType::RegionType;
137
139 void
140 SetImage(const TImage * image);
141
143 const TImage *
144 GetImage() const;
145
148 Size() const override;
149
153
156 GetTotalFrequency() const override;
157
159
161
162 using InstanceIdentifierVectorType = std::vector<InstanceIdentifier>;
164
166 void
168
173
176 itkSetMacro(UsePixelContainer, bool);
177 itkGetConstMacro(UsePixelContainer, bool);
178 itkBooleanMacro(UsePixelContainer);
180 // void PrintSelf(std::ostream& os, Indent indent) const override;
181
187 {
189
190 public:
191 ConstIterator() = delete;
192
193 ConstIterator(const JointDomainImageToListSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
194
199
202 {
204 return *this;
205 }
206
209 {
210 return 1;
211 }
212
215 {
217 return this->m_MeasurementVectorCache;
218 }
219
222 {
224 }
225
228 {
230 return *this;
231 }
232
233 bool
234 operator==(const ConstIterator & it) const
235 {
237 }
238
240
241 protected:
242 // This method should only be available to the ListSample class
247
248 private:
252 };
253
258 class Iterator : public ConstIterator
259 {
261
262 public:
263 Iterator() = delete;
264 Iterator(const Self * adaptor) = delete;
265 Iterator(const ConstIterator & it) = delete;
267 operator=(const ConstIterator & it) = delete;
268
269 Iterator(Self * adaptor)
270 : ConstIterator(adaptor)
271 {}
272
273 Iterator(const Iterator & iter)
274 : ConstIterator(iter)
275 {}
276
277 Iterator &
278 operator=(const Iterator & iter)
279 {
280 this->ConstIterator::operator=(iter);
281 return *this;
282 }
283
284 protected:
286 : ConstIterator(adaptor, iid)
287 {}
288 };
289
291 Iterator
293 {
294 const Iterator iter(this, 0);
295
296 return iter;
297 }
298
300 Iterator
302 {
303 const Iterator iter(this, m_Image->GetPixelContainer()->Size());
304
305 return iter;
306 }
307
309 ConstIterator
310 Begin() const
311 {
312 ConstIterator iter(this, 0);
313
314 return iter;
315 }
316
318 ConstIterator
319 End() const
320 {
321 ConstIterator iter(this, m_Image->GetPixelContainer()->Size());
322
323 return iter;
324 }
325
326protected:
329 void
330 PrintSelf(std::ostream & os, Indent indent) const override;
331
332private:
340
342}; // end of class JointDomainImageToListSampleAdaptor
343} // namespace itk::Statistics
344
345#ifndef ITK_MANUAL_INSTANTIATION
346# include "itkJointDomainImageToListSampleAdaptor.hxx"
347#endif
348
349#endif
Simulate a standard C array with copy semantics.
A multi-dimensional iterator templated over image type that walks a region of pixels.
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...
Traits for a pixel that define the dimension and component type.
typename typename TImage::PixelType::ValueType ValueType
A templated class holding a geometric point in n-Dimensional space.
Definition itkPoint.h:54
Implements transparent reference counting.
ConstIterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
Iterator(const JointDomainImageToListSampleAdaptor *adaptor, InstanceIdentifier iid)
ConstIterator & operator=(const ConstIterator &it)=delete
typename ImageType::PixelContainerConstPointer PixelContainerConstPointer
typename ImageJointDomainTraitsType::RangeDomainMeasurementType RangeDomainMeasurementType
TotalAbsoluteFrequencyType GetTotalFrequency() const override
FixedArray< float, Self::MeasurementVectorSize > NormalizationFactorsType
InstanceIdentifier Size() const override
void PrintSelf(std::ostream &os, Indent indent) const override
ListSample< typename ImageJointDomainTraits< TImage >::MeasurementVectorType > Superclass
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
FixedArray< RangeDomainMeasurementType, Self::RangeDomainDimension > RangeDomainMeasurementVectorType
unsigned int MeasurementVectorSizeType
Definition itkSample.h:92
MeasurementVectorTraits::InstanceIdentifier InstanceIdentifier
Definition itkSample.h:89
NumericTraits< AbsoluteFrequencyType >::AccumulateType TotalAbsoluteFrequencyType
Definition itkSample.h:85
MeasurementVectorTraits::AbsoluteFrequencyType AbsoluteFrequencyType
Definition itkSample.h:82
This class provides the type definition for the measurement vector in the joint domain (range domain ...
JoinTraits< RangeDomainMeasurementType, CoordinateRepType > JoinTraitsType
FixedArray< MeasurementType, Self::Dimension > MeasurementVectorType
Point< CoordinateRepType, Self::ImageDimension > PointType
PixelTraits< typename TImage::PixelType > PixelTraitsType