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"
25#include "itkListSample.h"
26
27namespace itk::Statistics
28{
37template <typename TImage>
54
84
85template <typename TImage>
86class ITK_TEMPLATE_EXPORT JointDomainImageToListSampleAdaptor
87 : public ListSample<typename ImageJointDomainTraits<TImage>::MeasurementVectorType>
88{
89public:
90 ITK_DISALLOW_COPY_AND_MOVE(JointDomainImageToListSampleAdaptor);
91
94
96
99
101
107
109 itkOverrideGetNameOfClassMacro(JointDomainImageToListSampleAdaptor);
110
112 itkNewMacro(Self);
113
116
118
123 using typename Superclass::InstanceIdentifier;
124
126 using ImageType = TImage;
129
130 using ImagePointer = typename ImageType::Pointer;
131 using ImageConstPointer = typename ImageType::ConstPointer;
132 using PixelType = typename ImageType::PixelType;
133 using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
134 using ImageIndexType = typename ImageType::IndexType;
135 using ImageSizeType = typename ImageType::SizeType;
136 using ImageRegionType = typename ImageType::RegionType;
138
140 void
141 SetImage(const TImage * image);
142
144 const TImage *
145 GetImage() const;
146
149 Size() const override;
150
154
157 GetTotalFrequency() const override;
158
160
162
163 using InstanceIdentifierVectorType = std::vector<InstanceIdentifier>;
165
167 void
169
174
177 itkSetMacro(UsePixelContainer, bool);
178 itkGetConstMacro(UsePixelContainer, bool);
179 itkBooleanMacro(UsePixelContainer);
181 // void PrintSelf(std::ostream& os, Indent indent) const override;
182
188 {
190
191 public:
192 ConstIterator() = delete;
193
194 ConstIterator(const JointDomainImageToListSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
195
200
203 {
205 return *this;
206 }
207
210 {
211 return 1;
212 }
213
216 {
218 return this->m_MeasurementVectorCache;
219 }
220
223 {
225 }
226
229 {
231 return *this;
232 }
233
234 bool
235 operator==(const ConstIterator & it) const
236 {
238 }
239
241
242 protected:
243 // This method should only be available to the ListSample class
248
249 private:
253 };
254
259 class Iterator : public ConstIterator
260 {
262
263 public:
264 Iterator() = delete;
265 Iterator(const Self * adaptor) = delete;
266 Iterator(const ConstIterator & it) = delete;
268 operator=(const ConstIterator & it) = delete;
269
270 Iterator(Self * adaptor)
271 : ConstIterator(adaptor)
272 {}
273
274 Iterator(const Iterator & iter)
275 : ConstIterator(iter)
276 {}
277
278 Iterator &
279 operator=(const Iterator & iter)
280 {
281 this->ConstIterator::operator=(iter);
282 return *this;
283 }
284
285 protected:
287 : ConstIterator(adaptor, iid)
288 {}
289 };
290
292 Iterator
294 {
295 const Iterator iter(this, 0);
296
297 return iter;
298 }
299
301 Iterator
303 {
304 const Iterator iter(this, m_Image->GetPixelContainer()->Size());
305
306 return iter;
307 }
308
310 ConstIterator
311 Begin() const
312 {
313 ConstIterator iter(this, 0);
314
315 return iter;
316 }
317
319 ConstIterator
320 End() const
321 {
322 ConstIterator iter(this, m_Image->GetPixelContainer()->Size());
323
324 return iter;
325 }
326
327protected:
330 void
331 PrintSelf(std::ostream & os, Indent indent) const override;
332
333private:
341
343}; // end of class JointDomainImageToListSampleAdaptor
344} // namespace itk::Statistics
345
346#ifndef ITK_MANUAL_INSTANTIATION
347# include "itkJointDomainImageToListSampleAdaptor.hxx"
348#endif
349
350#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
typename 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