ITK  5.4.0
Insight Toolkit
itkImageGaussianModelEstimator.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 itkImageGaussianModelEstimator_h
19 #define itkImageGaussianModelEstimator_h
20 
21 #include <cmath>
22 #include <cfloat>
23 #include <memory> // For unique_ptr.
24 
25 #include "vnl/vnl_vector.h"
26 #include "vnl/vnl_matrix.h"
27 #include "vnl/vnl_matrix_fixed.h"
28 #include "itkMath.h"
29 #include "vnl/algo/vnl_matrix_inverse.h"
30 
31 #include "itkImageRegionIterator.h"
32 #include "itkMacro.h"
33 
35 
36 namespace itk
37 {
76 template <typename TInputImage, typename TMembershipFunction, typename TTrainingImage>
77 class ITK_TEMPLATE_EXPORT ImageGaussianModelEstimator : public ImageModelEstimatorBase<TInputImage, TMembershipFunction>
78 {
79 public:
80  ITK_DISALLOW_COPY_AND_MOVE(ImageGaussianModelEstimator);
81 
87 
89  itkNewMacro(Self);
90 
92  itkOverrideGetNameOfClassMacro(ImageGaussianModelEstimator);
93 
95  using InputImageType = TInputImage;
98 
100  using TrainingImageType = TTrainingImage;
103 
106  using InputImagePixelType = typename TInputImage::PixelType;
107 
110  using TrainingImagePixelType = typename TTrainingImage::PixelType;
111 
117 
119  using MembershipFunctionType = TMembershipFunction;
121 
123  itkSetObjectMacro(TrainingImage, TrainingImageType);
124  itkGetModifiableObjectMacro(TrainingImage, TrainingImageType);
127 protected:
128  ImageGaussianModelEstimator() = default;
129  ~ImageGaussianModelEstimator() override = default;
130  void
131  PrintSelf(std::ostream & os, Indent indent) const override;
132 
134  void
135  GenerateData() override;
136 
137 private:
138  using MatrixType = vnl_matrix<double>;
139 
141 
143  static constexpr unsigned int VectorDimension = InputImagePixelType::Dimension;
144 
150  void
151  EstimateModels() override;
152 
153  void
154  EstimateGaussianModelParameters();
155 
156  MatrixType m_NumberOfSamples{};
157  MatrixType m_Means{};
158  std::unique_ptr<MatrixType[]> m_Covariance{ nullptr };
159 
160  TrainingImagePointer m_TrainingImage{};
161 };
162 } // end namespace itk
163 
164 #ifndef ITK_MANUAL_INSTANTIATION
165 # include "itkImageGaussianModelEstimator.hxx"
166 #endif
167 
168 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ImageGaussianModelEstimator::TrainingImageConstPointer
typename TTrainingImage::ConstPointer TrainingImageConstPointer
Definition: itkImageGaussianModelEstimator.h:102
itk::ImageGaussianModelEstimator::InputImagePointer
typename TInputImage::Pointer InputImagePointer
Definition: itkImageGaussianModelEstimator.h:96
itkImageModelEstimatorBase.h
itk::ImageGaussianModelEstimator::InputImageSizeType
typename TInputImage::SizeType InputImageSizeType
Definition: itkImageGaussianModelEstimator.h:140
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageGaussianModelEstimator::MembershipFunctionType
TMembershipFunction MembershipFunctionType
Definition: itkImageGaussianModelEstimator.h:119
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::ImageGaussianModelEstimator::InputImageConstPointer
typename TInputImage::ConstPointer InputImageConstPointer
Definition: itkImageGaussianModelEstimator.h:97
itk::ImageModelEstimatorBase
Base class for model estimation from images used for classification.
Definition: itkImageModelEstimatorBase.h:64
itk::ImageGaussianModelEstimator::MatrixType
vnl_matrix< double > MatrixType
Definition: itkImageGaussianModelEstimator.h:138
itk::ImageGaussianModelEstimator::MembershipFunctionPointer
typename TMembershipFunction::Pointer MembershipFunctionPointer
Definition: itkImageGaussianModelEstimator.h:120
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageGaussianModelEstimator::TrainingImagePointer
typename TTrainingImage::Pointer TrainingImagePointer
Definition: itkImageGaussianModelEstimator.h:101
itk::ImageGaussianModelEstimator::TrainingImagePixelType
typename TTrainingImage::PixelType TrainingImagePixelType
Definition: itkImageGaussianModelEstimator.h:110
itk::ImageGaussianModelEstimator::TrainingImageType
TTrainingImage TrainingImageType
Definition: itkImageGaussianModelEstimator.h:100
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::ImageGaussianModelEstimator::InputImageType
TInputImage InputImageType
Definition: itkImageGaussianModelEstimator.h:95
itk::ImageGaussianModelEstimator
Base class for ImageGaussianModelEstimator object.
Definition: itkImageGaussianModelEstimator.h:77
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::ImageGaussianModelEstimator::InputImagePixelType
typename TInputImage::PixelType InputImagePixelType
Definition: itkImageGaussianModelEstimator.h:106
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itkMath.h