ITK  6.0.0
Insight Toolkit
itkImageKmeansModelEstimator.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 itkImageKmeansModelEstimator_h
19 #define itkImageKmeansModelEstimator_h
20 
21 #include <ctime>
22 #include <cmath>
23 #include <cfloat>
24 
25 #include "vnl/vnl_vector.h"
26 #include "vnl/vnl_matrix.h"
27 #include "itkMath.h"
28 #include "vnl/algo/vnl_matrix_inverse.h"
29 
30 #include "itkImageRegionIterator.h"
31 #include "itkMacro.h"
32 
34 
35 #define ONEBAND 1
36 #define GLA_CONVERGED 1
37 #define GLA_NOT_CONVERGED 2
38 #define LBG_COMPLETED 3
39 
40 namespace itk
41 {
130 template <typename TInputImage, typename TMembershipFunction>
131 class ITK_TEMPLATE_EXPORT ImageKmeansModelEstimator : public ImageModelEstimatorBase<TInputImage, TMembershipFunction>
132 {
133 public:
134  ITK_DISALLOW_COPY_AND_MOVE(ImageKmeansModelEstimator);
140 
143 
145  itkNewMacro(Self);
146 
148  itkOverrideGetNameOfClassMacro(ImageKmeansModelEstimator);
149 
151  using InputImageType = TInputImage;
154 
158 
160  using InputImagePixelType = typename TInputImage::PixelType;
161 
164 
166 
169 
171  using CodebookMatrixOfDoubleType = vnl_matrix<double>;
172 
174  using CodebookMatrixOfIntegerType = vnl_matrix<int>;
175 
177  void
178  SetCodebook(CodebookMatrixOfDoubleType inCodebook);
179 
181  itkGetConstMacro(Codebook, CodebookMatrixOfDoubleType);
182 
186  {
187  return m_Codebook;
188  }
189 
191  itkSetMacro(Threshold, double);
192  itkGetConstMacro(Threshold, double);
196  itkSetMacro(OffsetAdd, double);
197  itkGetConstMacro(OffsetAdd, double);
201  itkSetMacro(OffsetMultiply, double);
202  itkGetConstMacro(OffsetMultiply, double);
206  itkSetMacro(MaxSplitAttempts, int);
207  itkGetConstMacro(MaxSplitAttempts, int);
211  CodebookMatrixOfDoubleType
213  {
214  return m_Centroid;
215  }
216 
217 protected:
219  ~ImageKmeansModelEstimator() override = default;
220  void
221  PrintSelf(std::ostream & os, Indent indent) const override;
222 
224  void
225  GenerateData() override;
226 
228  void
229  Allocate();
230 
232  void
233  PrintKmeansAlgorithmResults();
234 
235 private:
244  void
245  EstimateModels() override;
246 
248  void
249  EstimateKmeansModelParameters();
250 
252 
255 
257  void
258  Reallocate(int oldSize, int newSize);
259 
260  int
261  WithCodebookUseGLA();
262 
263  int
264  WithoutCodebookUseLBG();
265 
266  void
267  NearestNeighborSearchBasic(double * distortion);
268 
269  void
270  SplitCodewords(int currentSize, int numDesired, int scale);
271 
272  void
273  Perturb(double * oldCodeword, int scale, double * newCodeword);
274 
276 
277  // Buffer for K-means calculations
279 
280  double m_Threshold{};
281  double m_OffsetAdd{};
282  double m_OffsetMultiply{};
283  int m_MaxSplitAttempts{};
284 
285  bool m_ValidInCodebook{};
286  double m_DoubleMaximum{};
287  double m_OutputDistortion{};
288  int m_OutputNumberOfEmptyCells{};
289 
290  SizeValueType m_VectorDimension{};
291  SizeValueType m_NumberOfCodewords{};
292  SizeValueType m_CurrentNumberOfCodewords{};
293 
294  CodebookMatrixOfIntegerType m_CodewordHistogram{};
295  CodebookMatrixOfDoubleType m_CodewordDistortion{};
296 }; // class ImageKmeansModelEstimator
297 
298 } // end namespace itk
299 
300 #ifndef ITK_MANUAL_INSTANTIATION
301 # include "itkImageKmeansModelEstimator.hxx"
302 #endif
303 
304 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ImageKmeansModelEstimator::GetOutCodebook
CodebookMatrixOfDoubleType GetOutCodebook()
Definition: itkImageKmeansModelEstimator.h:185
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ImageKmeansModelEstimator::CodebookMatrixOfIntegerType
vnl_matrix< int > CodebookMatrixOfIntegerType
Definition: itkImageKmeansModelEstimator.h:174
itk::ImageKmeansModelEstimator
Base class for ImageKmeansModelEstimator object.
Definition: itkImageKmeansModelEstimator.h:131
itk::GTest::TypedefsAndConstructors::Dimension2::VectorType
ImageBaseType::SpacingType VectorType
Definition: itkGTestTypedefsAndConstructors.h:53
itkImageModelEstimatorBase.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ImageKmeansModelEstimator::InputImagePointer
typename TInputImage::Pointer InputImagePointer
Definition: itkImageKmeansModelEstimator.h:152
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageKmeansModelEstimator::ImageSizeType
typename TInputImage::SizeType ImageSizeType
Definition: itkImageKmeansModelEstimator.h:251
itk::ImageKmeansModelEstimator::InputImagePixelType
typename TInputImage::PixelType InputImagePixelType
Definition: itkImageKmeansModelEstimator.h:160
itk::ImageKmeansModelEstimator::CodebookMatrixOfDoubleType
vnl_matrix< double > CodebookMatrixOfDoubleType
Definition: itkImageKmeansModelEstimator.h:171
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itkMacro.h
itk::ImageModelEstimatorBase
Base class for model estimation from images used for classification.
Definition: itkImageModelEstimatorBase.h:64
itk::ImageKmeansModelEstimator::InputImageConstPointer
typename TInputImage::ConstPointer InputImageConstPointer
Definition: itkImageKmeansModelEstimator.h:153
itk::ImageKmeansModelEstimator::InputImageVectorType
typename TInputImage::PixelType::VectorType InputImageVectorType
Definition: itkImageKmeansModelEstimator.h:157
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::ImageKmeansModelEstimator::GetKmeansResults
CodebookMatrixOfDoubleType GetKmeansResults()
Definition: itkImageKmeansModelEstimator.h:212
itk::ImageKmeansModelEstimator::MembershipFunctionPointer
typename TMembershipFunction::Pointer MembershipFunctionPointer
Definition: itkImageKmeansModelEstimator.h:168
itk::ImageKmeansModelEstimator::InputPixelVectorType
typename TInputImage::PixelType::VectorType InputPixelVectorType
Definition: itkImageKmeansModelEstimator.h:254
itkMath.h
itk::ImageKmeansModelEstimator::InputImageType
TInputImage InputImageType
Definition: itkImageKmeansModelEstimator.h:151
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86