ITK  6.0.0
Insight Toolkit
itkMahalanobisDistanceThresholdImageFunction.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 itkMahalanobisDistanceThresholdImageFunction_h
19#define itkMahalanobisDistanceThresholdImageFunction_h
20
21#include "itkImageFunction.h"
23
24namespace itk
25{
49template <typename TInputImage, typename TCoordRep = float>
50class ITK_TEMPLATE_EXPORT MahalanobisDistanceThresholdImageFunction : public ImageFunction<TInputImage, bool, TCoordRep>
51{
52public:
53 ITK_DISALLOW_COPY_AND_MOVE(MahalanobisDistanceThresholdImageFunction);
54
60
62 itkOverrideGetNameOfClassMacro(MahalanobisDistanceThresholdImageFunction);
63
65 itkNewMacro(Self);
66
68 using typename Superclass::InputImageType;
69
71 using PixelType = typename TInputImage::PixelType;
72
74 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
75
77 using typename Superclass::PointType;
78
80 using typename Superclass::IndexType;
81
83 using typename Superclass::ContinuousIndexType;
84
86 using CovarianceMatrixType = vnl_matrix<double>;
87
89 using MeanVectorType = vnl_vector<double>;
90
99 bool
100 Evaluate(const PointType & point) const override;
101
110 bool
111 EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override;
112
121 bool
122 EvaluateAtIndex(const IndexType & index) const override;
123
128 virtual double
130
135 virtual double
136 EvaluateDistanceAtIndex(const IndexType & index) const;
137
139 itkGetConstReferenceMacro(Threshold, double);
140 itkSetMacro(Threshold, double);
145 void
146 SetMean(const MeanVectorType & mean);
147
150 itkGetConstReferenceMacro(Mean, MeanVectorType);
151
154 void
156
159 itkGetConstReferenceMacro(Covariance, CovarianceMatrixType);
160
161protected:
164 void
165 PrintSelf(std::ostream & os, Indent indent) const override;
166
167private:
168 double m_Threshold{};
169
170 // This is intended only for Image of Vector pixel type.
172
174 MahalanobisDistanceFunctionPointer m_MahalanobisDistanceMembershipFunction{};
175
176 // Cached versions of the mean and covariance to manage the
177 // difference in vector/matrix types between this class and the
178 // membership function used internally.
180 CovarianceMatrixType m_Covariance{};
181};
182} // end namespace itk
183
184#ifndef ITK_MANUAL_INSTANTIATION
185# include "itkMahalanobisDistanceThresholdImageFunction.hxx"
186#endif
187
188#endif
Evaluates a function of an image at specified position.
typename InputImageType::IndexType IndexType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Returns true if the pixel value of a vector image has a Mahalanobis distance below the value specifie...
bool EvaluateAtIndex(const IndexType &index) const override
void PrintSelf(std::ostream &os, Indent indent) const override
virtual double EvaluateDistanceAtIndex(const IndexType &index) const
~MahalanobisDistanceThresholdImageFunction() override=default
virtual double EvaluateDistance(const PointType &point) const
void SetMean(const MeanVectorType &mean)
bool Evaluate(const PointType &point) const override
typename MahalanobisDistanceFunctionType::Pointer MahalanobisDistanceFunctionPointer
bool EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
void SetCovariance(const CovarianceMatrixType &covariance)
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:54
MahalanobisDistanceMembershipFunction models class membership using Mahalanobis distance.
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents