ITK  6.0.0
Insight Toolkit
itkSLICImageFilter.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 itkSLICImageFilter_h
19#define itkSLICImageFilter_h
20
22#include <mutex>
23
24namespace itk
25{
26
60template <typename TInputImage, typename TOutputImage, typename TDistancePixel = float>
61class ITK_TEMPLATE_EXPORT SLICImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
62{
63public:
64 ITK_DISALLOW_COPY_AND_MOVE(SLICImageFilter);
65
71
73 itkNewMacro(Self);
74
76 itkOverrideGetNameOfClassMacro(SLICImageFilter);
77
79 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
80
81
83 using InputImageType = TInputImage;
84 using InputPixelType = typename InputImageType::PixelType;
85 using OutputImageType = TOutputImage;
86 using OutputPixelType = typename OutputImageType::PixelType;
87 using DistanceType = TDistancePixel;
89
93
94 using ClusterComponentType = double;
95 using ClusterType = vnl_vector_ref<ClusterComponentType>;
96
98
100
107 itkSetMacro(SpatialProximityWeight, double);
108 itkGetConstMacro(SpatialProximityWeight, double);
109
114 itkSetMacro(MaximumNumberOfIterations, unsigned int);
115 itkGetConstMacro(MaximumNumberOfIterations, unsigned int);
116
124 itkSetMacro(SuperGridSize, SuperGridSizeType);
125 itkGetConstMacro(SuperGridSize, SuperGridSizeType);
126 void
127 SetSuperGridSize(unsigned int factor);
128 void
129 SetSuperGridSize(unsigned int i, unsigned int factor);
130
138 itkSetMacro(InitializationPerturbation, bool);
139 itkGetMacro(InitializationPerturbation, bool);
140 itkBooleanMacro(InitializationPerturbation);
141
142
150 itkSetMacro(EnforceConnectivity, bool);
151 itkGetMacro(EnforceConnectivity, bool);
152 itkBooleanMacro(EnforceConnectivity);
153
154
161 itkGetConstMacro(AverageResidual, double);
162
163protected:
165 ~SLICImageFilter() override = default;
166
167 void
168 PrintSelf(std::ostream & os, Indent indent) const override;
169
171 void
173
174 void
176
177 void
179
180 void
181 ThreadedUpdateClusters(const OutputImageRegionType & updateRegionForThread);
182
183 void
185
186 void
188
189 void
191
192 void
193 GenerateData() override;
194
195
196 void
198
200 Distance(const ClusterType & cluster1, const ClusterType & cluster2);
201
203 Distance(const ClusterType & cluster, const InputPixelType & _v, const PointType & pt);
204
205private:
206 SuperGridSizeType m_SuperGridSize{};
207 unsigned int m_MaximumNumberOfIterations{};
208 double m_SpatialProximityWeight{ 10.0 };
209
211 std::vector<ClusterComponentType> m_Clusters{};
212 std::vector<ClusterComponentType> m_OldClusters{};
213
214
215 void
217 OutputPixelType requiredLabel,
218 OutputPixelType outputLabel,
219 std::vector<IndexType> & indexStack);
220
222 {
223 size_t count;
224 vnl_vector<ClusterComponentType> cluster;
225 };
226
227 using UpdateClusterMap = std::map<size_t, UpdateCluster>;
228
230
231 std::vector<UpdateClusterMap> m_UpdateClusterPerThread{};
232
233 typename DistanceImageType::Pointer m_DistanceImage{};
234 typename MarkerImageType::Pointer m_MarkerImage{};
235
236 bool m_EnforceConnectivity{ true };
237
238 bool m_InitializationPerturbation{ true };
239
240 double m_AverageResidual{};
241 std::mutex m_Mutex{};
242};
243} // end namespace itk
244
245#ifndef ITK_MANUAL_INSTANTIATION
246# include "itkSLICImageFilter.hxx"
247#endif
248
249#endif // itkSLICImageFilter_h
A templated class holding a point in n-Dimensional image space.
Base class for all data objects in ITK.
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Simple Linear Iterative Clustering (SLIC) super-pixel segmentation.
typename InputImageType::PixelType InputPixelType
void EnlargeOutputRequestedRegion(DataObject *output) override
typename OutputImageType::PixelType OutputPixelType
void SingleThreadedConnectivity()
void ThreadedUpdateClusters(const OutputImageRegionType &updateRegionForThread)
void ThreadedPerturbClusters(SizeValueType clusterIndex)
void BeforeThreadedGenerateData() override
void PrintSelf(std::ostream &os, Indent indent) const override
void AfterThreadedGenerateData() override
void RelabelConnectedRegion(const IndexType &seed, OutputPixelType requiredLabel, OutputPixelType outputLabel, std::vector< IndexType > &indexStack)
std::map< vcl_size_t, UpdateCluster > UpdateClusterMap
void GenerateData() override
TDistancePixel DistanceType
DistanceType Distance(const ClusterType &cluster1, const ClusterType &cluster2)
DistanceType Distance(const ClusterType &cluster, const InputPixelType &_v, const PointType &pt)
typename InputImageType::PointType PointType
void SetSuperGridSize(unsigned int i, unsigned int factor)
~SLICImageFilter() override=default
typename InputImageType::IndexType IndexType
vnl_vector_ref< ClusterComponentType > ClusterType
void SetSuperGridSize(unsigned int factor)
void ThreadedConnectivity(SizeValueType clusterIndex)
void ThreadedUpdateDistanceAndLabel(const OutputImageRegionType &outputRegionForThread)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86
vnl_vector< ClusterComponentType > cluster