ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkImageToNeighborhoodSampleAdaptor.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 itkImageToNeighborhoodSampleAdaptor_h
19#define itkImageToNeighborhoodSampleAdaptor_h
20
21#include <typeinfo>
22#include <vector>
23#include <iostream>
24
25#include "itkImage.h"
26#include "itkListSample.h"
27#include "itkSmartPointer.h"
29#include "itkMacro.h"
32
33namespace itk
34{
35namespace Statistics
36{
37
53
54template <typename TImage, typename TBoundaryCondition>
55class ITK_TEMPLATE_EXPORT ImageToNeighborhoodSampleAdaptor
56 : public ListSample<std::vector<ConstNeighborhoodIterator<TImage, TBoundaryCondition>>>
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(ImageToNeighborhoodSampleAdaptor);
60
63
65
68
70 itkOverrideGetNameOfClassMacro(ImageToNeighborhoodSampleAdaptor);
71
73 itkNewMacro(Self);
74
76 using ImageType = TImage;
77 using ImagePointer = typename ImageType::Pointer;
78 using ImageConstPointer = typename ImageType::ConstPointer;
79 using IndexType = typename ImageType::IndexType;
80 using OffsetType = typename ImageType::OffsetType;
81 using OffsetValueType = typename ImageType::OffsetValueType;
82 using PixelType = typename ImageType::PixelType;
83 using PixelContainerConstPointer = typename ImageType::PixelContainerConstPointer;
84 using RegionType = typename ImageType::RegionType;
85 using OffsetTableType = typename RegionType::OffsetTableType;
86 using SizeType = typename ImageType::SizeType;
88
96
99 using MeasurementVectorType = typename std::vector<ConstNeighborhoodIterator<TImage, TBoundaryCondition>>;
100 using ValueType = typename MeasurementVectorType::value_type;
102
106 using typename Superclass::InstanceIdentifier;
107
109 void
110 SetImage(const TImage * image);
111
113 const TImage *
114 GetImage() const;
115
117 void
119
122 GetRadius() const;
123
125 void
126 SetRegion(const RegionType & region);
127
130 GetRegion() const;
131
132 void
133 SetUseImageRegion(const bool flag);
134
136 itkGetConstMacro(UseImageRegion, bool);
137
139 itkBooleanMacro(UseImageRegion);
140
141
144 Size() const override;
145
149
153
156 GetTotalFrequency() const override;
157
164 {
166
167 public:
168 ConstIterator() = delete;
169
170 ConstIterator(const ImageToNeighborhoodSampleAdaptor * adaptor) { *this = adaptor->Begin(); }
171
176
178 operator=(const ConstIterator & iter) = default;
179
182 {
183 return 1;
184 }
185
188 {
189 return this->m_MeasurementVectorCache;
190 }
191
194 {
196 }
197
200 {
203 return *this;
204 }
205
206 bool
207 operator==(const ConstIterator & it) const
208 {
210 }
211
213
214 protected:
215 // This method should only be available to the ListSample class
218 {
219 this->m_MeasurementVectorCache.clear();
220 this->m_MeasurementVectorCache.push_back(iter);
221 }
222
223 private:
226 };
227
233 class Iterator : public ConstIterator
234 {
235
237
238 public:
239 Iterator() = delete;
240 Iterator(const Self * adaptor) = delete;
241 Iterator(const ConstIterator & it) = delete;
243 operator=(const ConstIterator & it) = delete;
244
245 Iterator(Self * adaptor)
246 : ConstIterator(adaptor)
247 {}
248
249 Iterator(const Iterator & iter)
250 : ConstIterator(iter)
251 {}
252
253 Iterator &
254 operator=(const Iterator & iter)
255 {
256 this->ConstIterator::operator=(iter);
257 return *this;
258 }
259
260 protected:
261 // This copy constructor is actually used in Iterator Begin()!
265 };
266
268 Iterator
270 {
272 nIterator.GoToBegin();
273 const Iterator iter(nIterator, 0);
274 return iter;
275 }
276
278 Iterator
280 {
282 nIterator.GoToEnd();
283 const Iterator iter(nIterator, m_Region.GetNumberOfPixels());
284 return iter;
285 }
286
288 ConstIterator
289 Begin() const
290 {
292 nIterator.GoToBegin();
293 ConstIterator iter(nIterator, 0);
294 return iter;
295 }
296
298 ConstIterator
299 End() const
300 {
302 nIterator.GoToEnd();
303 ConstIterator iter(nIterator, m_Region.GetNumberOfPixels());
304 return iter;
305 }
306
307protected:
310 void
311 PrintSelf(std::ostream & os, Indent indent) const override;
312
313private:
320 bool m_UseImageRegion{ true };
322
323}; // end of class ImageToNeighborhoodSampleAdaptor
324
325} // end of namespace Statistics
326
327template <typename TImage, typename TBoundaryCondition>
328std::ostream &
329operator<<(std::ostream & os, const std::vector<itk::ConstNeighborhoodIterator<TImage, TBoundaryCondition>> & mv);
330
331} // end of namespace itk
332
333#ifndef ITK_MANUAL_INSTANTIATION
334# include "itkImageToNeighborhoodSampleAdaptor.hxx"
335#endif
336
337#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Control indentation during Print() invocation.
Definition itkIndent.h:50
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Implements transparent reference counting.
ConstIterator & operator=(const ConstIterator &iter)=default
ConstIterator & operator=(const ConstIterator &it)=delete
Iterator(NeighborhoodIteratorType iter, InstanceIdentifier iid)
InstanceIdentifier Size() const override
typename std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > MeasurementVectorType
typename NeighborhoodIteratorType::SizeType NeighborhoodSizeType
ConstNeighborhoodIterator< TImage, TBoundaryCondition > NeighborhoodIteratorType
const MeasurementVectorType & GetMeasurementVector(InstanceIdentifier id) const override
typename ImageType::PixelContainerConstPointer PixelContainerConstPointer
void PrintSelf(std::ostream &os, Indent indent) const override
NeighborhoodIterator< TImage, TBoundaryCondition > NonConstNeighborhoodIteratorType
typename NeighborhoodIteratorType::NeighborhoodType NeighborhoodType
NeighborhoodRadiusType GetRadius() const
void SetRadius(const NeighborhoodRadiusType &radius)
AbsoluteFrequencyType GetFrequency(InstanceIdentifier id) const override
typename NeighborhoodIteratorType::RadiusType NeighborhoodRadiusType
ListSample< std::vector< ConstNeighborhoodIterator< TImage, TBoundaryCondition > > > Superclass
TotalAbsoluteFrequencyType GetTotalFrequency() const override
void SetRegion(const RegionType &region)
typename NeighborhoodIteratorType::IndexType NeighborhoodIndexType
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)