ITK  6.0.0
Insight Toolkit
itkEuclideanDistancePointMetric.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 itkEuclideanDistancePointMetric_h
19#define itkEuclideanDistancePointMetric_h
20
22#include "itkCovariantVector.h"
23#include "itkPointSet.h"
24#include "itkImage.h"
25
26namespace itk
27{
44template <typename TFixedPointSet,
45 typename TMovingPointSet,
47class ITK_TEMPLATE_EXPORT EuclideanDistancePointMetric
48 : public PointSetToPointSetMetric<TFixedPointSet, TMovingPointSet>
49{
50public:
51 ITK_DISALLOW_COPY_AND_MOVE(EuclideanDistancePointMetric);
52
56
59
61 itkNewMacro(Self);
62
64 itkOverrideGetNameOfClassMacro(EuclideanDistancePointMetric);
65
67 using typename Superclass::TransformType;
68 using typename Superclass::TransformPointer;
69 using typename Superclass::TransformParametersType;
70 using typename Superclass::TransformJacobianType;
71
72 using typename Superclass::MeasureType;
73 using typename Superclass::DerivativeType;
74 using typename Superclass::FixedPointSetType;
75 using typename Superclass::MovingPointSetType;
76 using typename Superclass::FixedPointSetConstPointer;
77 using typename Superclass::MovingPointSetConstPointer;
78
79 using typename Superclass::FixedPointIterator;
80 using typename Superclass::FixedPointDataIterator;
81
82 using typename Superclass::MovingPointIterator;
83 using typename Superclass::MovingPointDataIterator;
84
85 using DistanceMapType = TDistanceMap;
87
89 unsigned int
90 GetNumberOfValues() const override;
91
93 void
94 GetDerivative(const TransformParametersType & parameters, DerivativeType & Derivative) const override;
95
98 GetValue(const TransformParametersType & parameters) const override;
99
101 void
103 MeasureType & value,
104 DerivativeType & derivative) const;
105
107 itkSetConstObjectMacro(DistanceMap, DistanceMapType);
108 itkGetConstObjectMacro(DistanceMap, DistanceMapType);
115 itkSetMacro(ComputeSquaredDistance, bool);
116 itkGetConstMacro(ComputeSquaredDistance, bool);
117 itkBooleanMacro(ComputeSquaredDistance);
120protected:
122 ~EuclideanDistancePointMetric() override = default;
123
124 void
125 PrintSelf(std::ostream & os, Indent indent) const override;
126
127private:
128 DistanceMapPointer m_DistanceMap{};
129 bool m_ComputeSquaredDistance{ false };
130};
131} // end namespace itk
132
133#ifndef ITK_MANUAL_INSTANTIATION
134# include "itkEuclideanDistancePointMetric.hxx"
135#endif
136
137#endif
Array2D class representing a 2D array.
Definition: itkArray2D.h:43
Computes the minimum distance between a moving point-set and a fixed point-set. A vector of minimum c...
void GetDerivative(const TransformParametersType &parameters, DerivativeType &Derivative) const override
void GetValueAndDerivative(const TransformParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
MeasureType GetValue(const TransformParametersType &parameters) const override
typename DistanceMapType::ConstPointer DistanceMapPointer
~EuclideanDistancePointMetric() override=default
unsigned int GetNumberOfValues() const override
void PrintSelf(std::ostream &os, Indent indent) const override
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Computes similarity between two point sets.
typename TransformType::ParametersType TransformParametersType
SmartPointer< const Self > ConstPointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....