ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
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
40namespace itk
41{
128template <typename TInputImage, typename TMembershipFunction>
129class ITK_TEMPLATE_EXPORT ImageKmeansModelEstimator : public ImageModelEstimatorBase<TInputImage, TMembershipFunction>
130{
131public:
132 ITK_DISALLOW_COPY_AND_MOVE(ImageKmeansModelEstimator);
134
138
141
143 itkNewMacro(Self);
144
146 itkOverrideGetNameOfClassMacro(ImageKmeansModelEstimator);
147
149 using InputImageType = TInputImage;
150 using InputImagePointer = typename TInputImage::Pointer;
151 using InputImageConstPointer = typename TInputImage::ConstPointer;
152
155 using InputImageVectorType = typename TInputImage::PixelType::VectorType;
156
158 using InputImagePixelType = typename TInputImage::PixelType;
159
162
164
166 using MembershipFunctionPointer = typename TMembershipFunction::Pointer;
167
169 using CodebookMatrixOfDoubleType = vnl_matrix<double>;
170
172 using CodebookMatrixOfIntegerType = vnl_matrix<int>;
173
175 void
177
179 itkGetConstMacro(Codebook, CodebookMatrixOfDoubleType);
180
184 {
185 return m_Codebook;
186 }
187
189 itkSetMacro(Threshold, double);
190 itkGetConstMacro(Threshold, double);
192
194 itkSetMacro(OffsetAdd, double);
195 itkGetConstMacro(OffsetAdd, double);
197
199 itkSetMacro(OffsetMultiply, double);
200 itkGetConstMacro(OffsetMultiply, double);
202
204 itkSetMacro(MaxSplitAttempts, int);
205 itkGetConstMacro(MaxSplitAttempts, int);
207
211 {
212 return m_Centroid;
213 }
214
215protected:
217 ~ImageKmeansModelEstimator() override = default;
218 void
219 PrintSelf(std::ostream & os, Indent indent) const override;
220
222 void
223 GenerateData() override;
224
226 void
228
230 void
232
233private:
242 void
243 EstimateModels() override;
244
246 void
248
249 using ImageSizeType = typename TInputImage::SizeType;
250
252 using InputPixelVectorType = typename TInputImage::PixelType::VectorType;
253
255 void
256 Reallocate(int oldSize, int newSize);
257
258 int
260
261 int
263
264 void
265 NearestNeighborSearchBasic(double * distortion);
266
267 void
268 SplitCodewords(int currentSize, int numDesired, int scale);
269
270 void
271 Perturb(double * oldCodeword, int scale, double * newCodeword);
272
274
275 // Buffer for K-means calculations
277
278 double m_Threshold{};
279 double m_OffsetAdd{};
282
287
291
294}; // class ImageKmeansModelEstimator
295
296} // end namespace itk
297
298#ifndef ITK_MANUAL_INSTANTIATION
299# include "itkImageKmeansModelEstimator.hxx"
300#endif
301
302#endif
void Reallocate(int oldSize, int newSize)
ImageRegionConstIterator< TInputImage > InputImageConstIterator
typename TInputImage::Pointer InputImagePointer
typename TInputImage::PixelType::VectorType InputPixelVectorType
typename TMembershipFunction::Pointer MembershipFunctionPointer
~ImageKmeansModelEstimator() override=default
typename TInputImage::ConstPointer InputImageConstPointer
void NearestNeighborSearchBasic(double *distortion)
void PrintSelf(std::ostream &os, Indent indent) const override
CodebookMatrixOfIntegerType m_CodewordHistogram
void SplitCodewords(int currentSize, int numDesired, int scale)
CodebookMatrixOfDoubleType m_CodewordDistortion
typename TInputImage::PixelType InputImagePixelType
void SetCodebook(CodebookMatrixOfDoubleType inCodebook)
ImageRegionIterator< TInputImage > InputImageIterator
CodebookMatrixOfDoubleType GetKmeansResults()
void Perturb(double *oldCodeword, int scale, double *newCodeword)
typename TInputImage::PixelType::VectorType InputImageVectorType
ImageModelEstimatorBase< TInputImage, TMembershipFunction > Superclass
typename TInputImage::SizeType ImageSizeType
CodebookMatrixOfDoubleType GetOutCodebook()
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86