ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkSLICImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
21 #include "itkImageToImageFilter.h"
22 #include <mutex>
23 
24 namespace itk
25 {
26 
59 template < typename TInputImage, typename TOutputImage, typename TDistancePixel = float>
61  public ImageToImageFilter< TInputImage, TOutputImage >
62 {
63 public:
64  ITK_DISALLOW_COPY_AND_ASSIGN(SLICImageFilter);
65 
71 
73  itkNewMacro(Self);
74 
76  itkTypeMacro(SLICImageFilter, ImageToImageFilter);
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 SetSuperGridSize(unsigned int factor);
127  void SetSuperGridSize(unsigned int i, unsigned int factor);
128 
136  itkSetMacro(InitializationPerturbation, bool);
137  itkGetMacro(InitializationPerturbation, bool);
138  itkBooleanMacro(InitializationPerturbation);
139 
140 
148  itkSetMacro(EnforceConnectivity, bool);
149  itkGetMacro(EnforceConnectivity, bool);
150  itkBooleanMacro(EnforceConnectivity);
151 
152 
159  itkGetConstMacro( AverageResidual, double );
160 
161 protected:
162  SLICImageFilter();
163  ~SLICImageFilter() override = default;
164 
165  void PrintSelf(std::ostream & os, Indent indent) const override;
166 
168  void EnlargeOutputRequestedRegion(DataObject *output) override;
169 
170  void BeforeThreadedGenerateData() override;
171 
172  void ThreadedUpdateDistanceAndLabel(const OutputImageRegionType & outputRegionForThread);
173 
174  void ThreadedUpdateClusters(const OutputImageRegionType & outputRegionForThread);
175 
177 
179 
181 
182  void GenerateData() override;
183 
184 
185  void AfterThreadedGenerateData() override;
186 
187  DistanceType Distance(const ClusterType &cluster1, const ClusterType &cluster2);
188 
189  DistanceType Distance(const ClusterType &cluster, const InputPixelType &v, const PointType &pt);
190 
191 private:
192 
195  double m_SpatialProximityWeight{ 10.0 };
196 
198  std::vector<ClusterComponentType> m_Clusters;
199  std::vector<ClusterComponentType> m_OldClusters;
200 
201 
202  void RelabelConnectedRegion( const IndexType &seed,
203  OutputPixelType requiredLabel,
204  OutputPixelType outputLabel,
205  std::vector<IndexType> & indexStack);
206 
208  {
209  size_t count;
210  vnl_vector<ClusterComponentType> cluster;
211  };
212 
213  using UpdateClusterMap = std::map<size_t, UpdateCluster>;
214 
216 
217  std::vector<UpdateClusterMap> m_UpdateClusterPerThread;
218 
221 
223 
225 
227  std::mutex m_Mutex;
228 };
229 } // end namespace itk
230 
231 #ifndef ITK_MANUAL_INSTANTIATION
232 #include "itkSLICImageFilter.hxx"
233 #endif
234 
235 #endif //itkSLICImageFilter_h
void AfterThreadedGenerateData() override
std::map< vcl_size_t, UpdateCluster > UpdateClusterMap
TDistancePixel DistanceType
vnl_vector< ClusterComponentType > cluster
void BeforeThreadedGenerateData() override
DistanceType Distance(const ClusterType &cluster1, const ClusterType &cluster2)
std::vector< ClusterComponentType > m_OldClusters
void EnlargeOutputRequestedRegion(DataObject *output) override
MarkerImageType::Pointer m_MarkerImage
static constexpr unsigned int ImageDimension
std::vector< ClusterComponentType > m_Clusters
typename OutputImageType::PixelType OutputPixelType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
vnl_vector_ref< ClusterComponentType > ClusterType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
typename InputImageType::PointType PointType
typename InputImageType::IndexType IndexType
virtual void SetSuperGridSize(SuperGridSizeType _arg)
The expected superpixel size and shape.
Base class for all process objects that output image data.
void PrintSelf(std::ostream &os, Indent indent) const override
void SingleThreadedConnectivity()
unsigned int m_MaximumNumberOfIterations
SuperGridSizeType m_SuperGridSize
void ThreadedPerturbClusters(SizeValueType idx)
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
DistanceImageType::Pointer m_DistanceImage
FixedArray< double, ImageDimension > m_DistanceScales
void RelabelConnectedRegion(const IndexType &seed, OutputPixelType requiredLabel, OutputPixelType outputLabel, std::vector< IndexType > &indexStack)
std::vector< UpdateClusterMap > m_UpdateClusterPerThread
A templated class holding a point in n-Dimensional image space.
void ThreadedUpdateClusters(const OutputImageRegionType &outputRegionForThread)
Simple Linear Iterative Clustering (SLIC) super-pixel segmentation.
~SLICImageFilter() override=default
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void GenerateData() override
typename InputImageType::PixelType InputPixelType
void ThreadedUpdateDistanceAndLabel(const OutputImageRegionType &outputRegionForThread)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void ThreadedConnectivity(SizeValueType idx)