ITK  5.1.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  * 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>
60 class SLICImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
61 {
62 public:
63  ITK_DISALLOW_COPY_AND_ASSIGN(SLICImageFilter);
64 
70 
72  itkNewMacro(Self);
73 
75  itkTypeMacro(SLICImageFilter, ImageToImageFilter);
76 
78  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
79 
80 
82  using InputImageType = TInputImage;
83  using InputPixelType = typename InputImageType::PixelType;
84  using OutputImageType = TOutputImage;
85  using OutputPixelType = typename OutputImageType::PixelType;
86  using DistanceType = TDistancePixel;
88 
92 
93  using ClusterComponentType = double;
94  using ClusterType = vnl_vector_ref<ClusterComponentType>;
95 
97 
99 
106  itkSetMacro(SpatialProximityWeight, double);
107  itkGetConstMacro(SpatialProximityWeight, double);
108 
113  itkSetMacro(MaximumNumberOfIterations, unsigned int);
114  itkGetConstMacro(MaximumNumberOfIterations, unsigned int);
115 
123  itkSetMacro(SuperGridSize, SuperGridSizeType);
124  itkGetConstMacro(SuperGridSize, SuperGridSizeType);
125  void
126  SetSuperGridSize(unsigned int factor);
127  void
128  SetSuperGridSize(unsigned int i, unsigned int factor);
129 
137  itkSetMacro(InitializationPerturbation, bool);
138  itkGetMacro(InitializationPerturbation, bool);
139  itkBooleanMacro(InitializationPerturbation);
140 
141 
149  itkSetMacro(EnforceConnectivity, bool);
150  itkGetMacro(EnforceConnectivity, bool);
151  itkBooleanMacro(EnforceConnectivity);
152 
153 
160  itkGetConstMacro(AverageResidual, double);
161 
162 protected:
163  SLICImageFilter();
164  ~SLICImageFilter() override = default;
165 
166  void
167  PrintSelf(std::ostream & os, Indent indent) const override;
168 
170  void
171  EnlargeOutputRequestedRegion(DataObject * output) override;
172 
173  void
174  BeforeThreadedGenerateData() override;
175 
176  void
177  ThreadedUpdateDistanceAndLabel(const OutputImageRegionType & outputRegionForThread);
178 
179  void
180  ThreadedUpdateClusters(const OutputImageRegionType & outputRegionForThread);
181 
182  void
184 
185  void
187 
188  void
190 
191  void
192  GenerateData() override;
193 
194 
195  void
196  AfterThreadedGenerateData() override;
197 
199  Distance(const ClusterType & cluster1, const ClusterType & cluster2);
200 
202  Distance(const ClusterType & cluster, const InputPixelType & v, const PointType & pt);
203 
204 private:
207  double m_SpatialProximityWeight{ 10.0 };
208 
210  std::vector<ClusterComponentType> m_Clusters;
211  std::vector<ClusterComponentType> m_OldClusters;
212 
213 
214  void
215  RelabelConnectedRegion(const IndexType & seed,
216  OutputPixelType requiredLabel,
217  OutputPixelType outputLabel,
218  std::vector<IndexType> & indexStack);
219 
221  {
222  size_t count;
223  vnl_vector<ClusterComponentType> cluster;
224  };
225 
226  using UpdateClusterMap = std::map<size_t, UpdateCluster>;
227 
229 
230  std::vector<UpdateClusterMap> m_UpdateClusterPerThread;
231 
234 
235  bool m_EnforceConnectivity{ true };
236 
238 
240  std::mutex m_Mutex;
241 };
242 } // end namespace itk
243 
244 #ifndef ITK_MANUAL_INSTANTIATION
245 # include "itkSLICImageFilter.hxx"
246 #endif
247 
248 #endif // itkSLICImageFilter_h
itk::SLICImageFilter::m_MaximumNumberOfIterations
unsigned int m_MaximumNumberOfIterations
Definition: itkSLICImageFilter.h:206
itk::SLICImageFilter::UpdateCluster
Definition: itkSLICImageFilter.h:220
itk::SLICImageFilter::m_Mutex
std::mutex m_Mutex
Definition: itkSLICImageFilter.h:240
itk::SLICImageFilter::m_UpdateClusterPerThread
std::vector< UpdateClusterMap > m_UpdateClusterPerThread
Definition: itkSLICImageFilter.h:230
itk::SLICImageFilter::UpdateCluster::count
vcl_size_t count
Definition: itkSLICImageFilter.h:222
itk::SLICImageFilter::SingleThreadedConnectivity
void SingleThreadedConnectivity()
itk::SLICImageFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
itk::SLICImageFilter::ThreadedConnectivity
void ThreadedConnectivity(SizeValueType idx)
itk::SLICImageFilter::ThreadedUpdateClusters
void ThreadedUpdateClusters(const OutputImageRegionType &outputRegionForThread)
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::SLICImageFilter::m_DistanceImage
DistanceImageType::Pointer m_DistanceImage
Definition: itkSLICImageFilter.h:232
itk::SLICImageFilter::m_AverageResidual
double m_AverageResidual
Definition: itkSLICImageFilter.h:239
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::SLICImageFilter::m_EnforceConnectivity
bool m_EnforceConnectivity
Definition: itkSLICImageFilter.h:235
itk::SLICImageFilter::AfterThreadedGenerateData
void AfterThreadedGenerateData() override
itk::SLICImageFilter::ThreadedUpdateDistanceAndLabel
void ThreadedUpdateDistanceAndLabel(const OutputImageRegionType &outputRegionForThread)
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::SLICImageFilter::ClusterType
vnl_vector_ref< ClusterComponentType > ClusterType
Definition: itkSLICImageFilter.h:94
itk::SLICImageFilter::m_InitializationPerturbation
bool m_InitializationPerturbation
Definition: itkSLICImageFilter.h:237
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::SLICImageFilter::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const override
itk::SLICImageFilter::DistanceType
TDistancePixel DistanceType
Definition: itkSLICImageFilter.h:86
itk::SLICImageFilter::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkSLICImageFilter.h:83
itk::SLICImageFilter::IndexType
typename InputImageType::IndexType IndexType
Definition: itkSLICImageFilter.h:89
itk::SLICImageFilter::Distance
DistanceType Distance(const ClusterType &cluster1, const ClusterType &cluster2)
itk::SLICImageFilter::m_MarkerImage
MarkerImageType::Pointer m_MarkerImage
Definition: itkSLICImageFilter.h:233
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::SLICImageFilter::~SLICImageFilter
~SLICImageFilter() override=default
itk::SLICImageFilter::OutputPixelType
typename OutputImageType::PixelType OutputPixelType
Definition: itkSLICImageFilter.h:85
itk::SLICImageFilter
Simple Linear Iterative Clustering (SLIC) super-pixel segmentation.
Definition: itkSLICImageFilter.h:60
itk::SLICImageFilter::m_Clusters
std::vector< ClusterComponentType > m_Clusters
Definition: itkSLICImageFilter.h:210
itk::SLICImageFilter::PointType
typename InputImageType::PointType PointType
Definition: itkSLICImageFilter.h:90
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::SLICImageFilter::SetSuperGridSize
virtual void SetSuperGridSize(SuperGridSizeType _arg)
The expected superpixel size and shape.
itk::SLICImageFilter::m_SpatialProximityWeight
double m_SpatialProximityWeight
Definition: itkSLICImageFilter.h:207
itk::FixedArray< unsigned int, ImageDimension >
itk::SLICImageFilter::EnlargeOutputRequestedRegion
void EnlargeOutputRequestedRegion(DataObject *output) override
itk::SLICImageFilter::UpdateClusterMap
std::map< vcl_size_t, UpdateCluster > UpdateClusterMap
Definition: itkSLICImageFilter.h:226
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkArray.h:26
itk::SLICImageFilter::GenerateData
void GenerateData() override
itk::SLICImageFilter::m_SuperGridSize
SuperGridSizeType m_SuperGridSize
Definition: itkSLICImageFilter.h:205
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::SLICImageFilter::UpdateCluster::cluster
vnl_vector< ClusterComponentType > cluster
Definition: itkSLICImageFilter.h:223
itk::SLICImageFilter::SLICImageFilter
SLICImageFilter()
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::SLICImageFilter::ImageDimension
static constexpr unsigned int ImageDimension
Definition: itkSLICImageFilter.h:78
itk::SLICImageFilter::RelabelConnectedRegion
void RelabelConnectedRegion(const IndexType &seed, OutputPixelType requiredLabel, OutputPixelType outputLabel, std::vector< IndexType > &indexStack)
itk::SLICImageFilter::m_DistanceScales
FixedArray< double, ImageDimension > m_DistanceScales
Definition: itkSLICImageFilter.h:209
itk::SLICImageFilter::ClusterComponentType
double ClusterComponentType
Definition: itkSLICImageFilter.h:93
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::SLICImageFilter::m_OldClusters
std::vector< ClusterComponentType > m_OldClusters
Definition: itkSLICImageFilter.h:211
itk::SLICImageFilter::ThreadedPerturbClusters
void ThreadedPerturbClusters(SizeValueType idx)
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293