ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkObjectToObjectMetricBase.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 itkObjectToObjectMetricBase_h
19#define itkObjectToObjectMetricBase_h
20
21#include "itkTransformBase.h"
23#include "ITKOptimizersv4Export.h"
24
25namespace itk
26{
60// Define how to print enumeration
61extern ITKOptimizersv4_EXPORT std::ostream &
63extern ITKOptimizersv4_EXPORT std::ostream &
65
89template <typename TInternalComputationValueType = double>
90class ITK_TEMPLATE_EXPORT ObjectToObjectMetricBaseTemplate
91 : public SingleValuedCostFunctionv4Template<TInternalComputationValueType>
92{
93public:
94 ITK_DISALLOW_COPY_AND_MOVE(ObjectToObjectMetricBaseTemplate);
95
101
103 itkOverrideGetNameOfClassMacro(ObjectToObjectMetricBaseTemplate);
104
106 using CoordinateRepresentationType = TInternalComputationValueType;
107
109 using typename Superclass::MeasureType;
110
112 using typename Superclass::DerivativeType;
114
116 using typename Superclass::ParametersType;
117 using ParametersValueType = TInternalComputationValueType;
118
122
124 itkSetConstObjectMacro(FixedObject, ObjectType);
125 itkGetConstObjectMacro(FixedObject, ObjectType);
127
129 itkSetConstObjectMacro(MovingObject, ObjectType);
130 itkGetConstObjectMacro(MovingObject, ObjectType);
132
134#if !defined(ITK_LEGACY_REMOVE)
136 // We need to expose the enum values at the class level
137 // for backwards compatibility
138 static constexpr GradientSourceEnum GRADIENT_SOURCE_FIXED = GradientSourceEnum::GRADIENT_SOURCE_FIXED;
139 static constexpr GradientSourceEnum GRADIENT_SOURCE_MOVING = GradientSourceEnum::GRADIENT_SOURCE_MOVING;
140 static constexpr GradientSourceEnum GRADIENT_SOURCE_BOTH = GradientSourceEnum::GRADIENT_SOURCE_BOTH;
141#endif
142
150
156
159 bool
161
164 bool
166
171 virtual void
173
176 using NumberOfParametersType = unsigned int;
177
182 GetValue() const override = 0;
183
187 virtual void
189
192 void
193 GetValueAndDerivative(MeasureType & value, DerivativeType & derivative) const override = 0;
194
199 GetNumberOfParameters() const override = 0;
203
205 virtual void
207
209 virtual const ParametersType &
210 GetParameters() const = 0;
211
214 virtual bool
215 HasLocalSupport() const = 0;
216
223 virtual void
226
234
236#if !defined(ITK_LEGACY_REMOVE)
238 static constexpr MetricCategoryEnum UNKNOWN_METRIC = MetricCategoryEnum::UNKNOWN_METRIC;
239 static constexpr MetricCategoryEnum OBJECT_METRIC = MetricCategoryEnum::OBJECT_METRIC;
240 static constexpr MetricCategoryEnum IMAGE_METRIC = MetricCategoryEnum::IMAGE_METRIC;
241 static constexpr MetricCategoryEnum POINT_SET_METRIC = MetricCategoryEnum::POINT_SET_METRIC;
242 static constexpr MetricCategoryEnum MULTI_METRIC = MetricCategoryEnum::MULTI_METRIC;
243#endif
244
246 virtual MetricCategoryEnum
248 {
249 return MetricCategoryEnum::UNKNOWN_METRIC;
250 }
251
252protected:
255
256 void
257 PrintSelf(std::ostream & os, Indent indent) const override;
258
262
264
267};
268
271} // end namespace itk
272
273#ifndef ITK_MANUAL_INSTANTIATION
274# include "itkObjectToObjectMetricBase.hxx"
275#endif
276
277#endif
TParametersValueType ValueType
Definition itkArray.h:52
Control indentation during Print() invocation.
Definition itkIndent.h:50
This class contains all the enum classes used by the ObjectToObjectMetricBaseTemplate class.
Base class for all object-to-object similarity metrics added in ITKv4.
virtual void GetDerivative(DerivativeType &) const =0
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory MetricCategoryEnum
virtual bool HasLocalSupport() const =0
virtual MetricCategoryEnum GetMetricCategory() const
itk::ObjectToObjectMetricBaseTemplateEnums::GradientSource GradientSourceEnum
MeasureType GetValue() const override=0
NumberOfParametersType GetNumberOfParameters() const override=0
void PrintSelf(std::ostream &os, Indent indent) const override
SingleValuedCostFunctionv4Template< TParametersValueType > Superclass
void GetValueAndDerivative(MeasureType &value, DerivativeType &derivative) const override=0
virtual void SetParameters(ParametersType &params)=0
virtual const ParametersType & GetParameters() const =0
virtual NumberOfParametersType GetNumberOfLocalParameters() const =0
~ObjectToObjectMetricBaseTemplate() override=default
virtual void UpdateTransformParameters(const DerivativeType &derivative, ParametersValueType factor=NumericTraits< ParametersValueType >::OneValue())=0
Base class for most ITK classes.
Definition itkObject.h:62
SmartPointer< const Self > ConstPointer
Definition itkObject.h:70
Array< TInternalComputationValueType > DerivativeType
OptimizerParameters< TInternalComputationValueType > ParametersType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ObjectToObjectMetricBaseTemplate< double > ObjectToObjectMetricBase
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)