ITK  5.2.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  * 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 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);
136 
140 
143 
145  itkNewMacro(Self);
146 
149 
151  using InputImageType = TInputImage;
152  using InputImagePointer = typename TInputImage::Pointer;
153  using InputImageConstPointer = typename TInputImage::ConstPointer;
154 
158 
160  using InputImagePixelType = typename TInputImage::PixelType;
161 
164 
166 
168  using MembershipFunctionPointer = typename TMembershipFunction::Pointer;
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 
194  itkGetConstMacro(Threshold, double);
195 
197  itkSetMacro(OffsetAdd, double);
198 
200  itkGetConstMacro(OffsetAdd, double);
201 
203  itkSetMacro(OffsetMultiply, double);
204 
206  itkGetConstMacro(OffsetMultiply, double);
207 
209  itkSetMacro(MaxSplitAttempts, int);
210 
212  itkGetConstMacro(MaxSplitAttempts, int);
213 
215  CodebookMatrixOfDoubleType
217  {
218  return m_Centroid;
219  }
220 
221 protected:
223  ~ImageKmeansModelEstimator() override = default;
224  void
225  PrintSelf(std::ostream & os, Indent indent) const override;
226 
228  void
229  GenerateData() override;
230 
232  void
233  Allocate();
234 
236  void
237  PrintKmeansAlgorithmResults();
238 
239 private:
246  void
247  EstimateModels() override;
248 
249  void
250  EstimateKmeansModelParameters();
251 
253 
256 
257  void
258  Reallocate(int oldSize, int newSize);
259 
260  // Local functions
261  int
262  WithCodebookUseGLA(); // GLA stands for the Generalized Lloyd Algorithm
263 
264  int
265  WithoutCodebookUseLBG(); // LBG stands for the Lindo Buzo Gray Algorithm
266 
267  void
268  NearestNeighborSearchBasic(double * distortion);
269 
270  void
271  SplitCodewords(int currentSize, int numDesired, int scale);
272 
273  void
274  Perturb(double * oldCodeword, int scale, double * newCodeword);
275 
277 
278  // Buffer for K-means calculations
280 
281  double m_Threshold;
282  double m_OffsetAdd;
285 
290 
294 
297 }; // class ImageKmeansModelEstimator
298 
299 } // namespace itk
300 
301 #ifndef ITK_MANUAL_INSTANTIATION
302 # include "itkImageKmeansModelEstimator.hxx"
303 #endif
304 
305 #endif
itk::ImageKmeansModelEstimator::GetOutCodebook
CodebookMatrixOfDoubleType GetOutCodebook()
Definition: itkImageKmeansModelEstimator.h:185
itk::ImageKmeansModelEstimator::CodebookMatrixOfIntegerType
vnl_matrix< int > CodebookMatrixOfIntegerType
Definition: itkImageKmeansModelEstimator.h:174
itk::ImageKmeansModelEstimator
Base class for ImageKmeansModelEstimator object.
Definition: itkImageKmeansModelEstimator.h:131
itk::ImageKmeansModelEstimator::m_OutputDistortion
double m_OutputDistortion
Definition: itkImageKmeansModelEstimator.h:288
itk::ImageKmeansModelEstimator::m_NumberOfCodewords
SizeValueType m_NumberOfCodewords
Definition: itkImageKmeansModelEstimator.h:292
itk::GTest::TypedefsAndConstructors::Dimension2::VectorType
ImageBaseType::SpacingType VectorType
Definition: itkGTestTypedefsAndConstructors.h:53
itk::ImageKmeansModelEstimator::m_CodewordDistortion
CodebookMatrixOfDoubleType m_CodewordDistortion
Definition: itkImageKmeansModelEstimator.h:296
itk::ImageKmeansModelEstimator::m_CodewordHistogram
CodebookMatrixOfIntegerType m_CodewordHistogram
Definition: itkImageKmeansModelEstimator.h:295
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:252
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::ImageKmeansModelEstimator::m_OffsetMultiply
double m_OffsetMultiply
Definition: itkImageKmeansModelEstimator.h:283
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::ImageKmeansModelEstimator::m_Threshold
double m_Threshold
Definition: itkImageKmeansModelEstimator.h:281
itkMacro.h
itk::ImageKmeansModelEstimator::m_CurrentNumberOfCodewords
SizeValueType m_CurrentNumberOfCodewords
Definition: itkImageKmeansModelEstimator.h:293
itk::ImageModelEstimatorBase
Base class for model estimation from images used for classification.
Definition: itkImageModelEstimatorBase.h:64
itk::ImageKmeansModelEstimator::m_OffsetAdd
double m_OffsetAdd
Definition: itkImageKmeansModelEstimator.h:282
itk::ImageKmeansModelEstimator::InputImageConstPointer
typename TInputImage::ConstPointer InputImageConstPointer
Definition: itkImageKmeansModelEstimator.h:153
itk::ImageKmeansModelEstimator::m_Centroid
CodebookMatrixOfDoubleType m_Centroid
Definition: itkImageKmeansModelEstimator.h:279
itk::ImageKmeansModelEstimator::InputImageVectorType
typename TInputImage::PixelType::VectorType InputImageVectorType
Definition: itkImageKmeansModelEstimator.h:157
itk::ImageKmeansModelEstimator::m_Codebook
CodebookMatrixOfDoubleType m_Codebook
Definition: itkImageKmeansModelEstimator.h:276
itk::ImageKmeansModelEstimator::m_DoubleMaximum
double m_DoubleMaximum
Definition: itkImageKmeansModelEstimator.h:287
itk::ImageKmeansModelEstimator::m_VectorDimension
SizeValueType m_VectorDimension
Definition: itkImageKmeansModelEstimator.h:291
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageKmeansModelEstimator::m_OutputNumberOfEmptyCells
int m_OutputNumberOfEmptyCells
Definition: itkImageKmeansModelEstimator.h:289
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
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:216
itk::ImageKmeansModelEstimator::MembershipFunctionPointer
typename TMembershipFunction::Pointer MembershipFunctionPointer
Definition: itkImageKmeansModelEstimator.h:168
itk::ImageKmeansModelEstimator::m_ValidInCodebook
bool m_ValidInCodebook
Definition: itkImageKmeansModelEstimator.h:286
itk::ImageKmeansModelEstimator::InputPixelVectorType
typename TInputImage::PixelType::VectorType InputPixelVectorType
Definition: itkImageKmeansModelEstimator.h:255
itkMath.h
itk::ImageKmeansModelEstimator::m_MaxSplitAttempts
int m_MaxSplitAttempts
Definition: itkImageKmeansModelEstimator.h:284
itk::ImageKmeansModelEstimator::InputImageType
TInputImage InputImageType
Definition: itkImageKmeansModelEstimator.h:151
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83