ITK  6.0.0
Insight Toolkit
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 
25 namespace itk
26 {
32 {
33 public:
40  enum class GradientSource : uint8_t
41  {
45  };
46 
51  enum class MetricCategory : uint8_t
52  {
53  UNKNOWN_METRIC = 0,
54  OBJECT_METRIC = 1,
55  IMAGE_METRIC = 2,
56  POINT_SET_METRIC = 3,
57  MULTI_METRIC = 4
58  };
59 };
60 // Define how to print enumeration
61 extern ITKOptimizersv4_EXPORT std::ostream &
63 extern ITKOptimizersv4_EXPORT std::ostream &
65 
89 template <typename TInternalComputationValueType = double>
90 class ITK_TEMPLATE_EXPORT ObjectToObjectMetricBaseTemplate
91  : public SingleValuedCostFunctionv4Template<TInternalComputationValueType>
92 {
93 public:
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);
129  itkSetConstObjectMacro(MovingObject, ObjectType);
130  itkGetConstObjectMacro(MovingObject, ObjectType);
134 #if !defined(ITK_LEGACY_REMOVE)
135 
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 
149  itkSetMacro(GradientSource, GradientSourceEnum);
150 
155  itkGetConstMacro(GradientSource, GradientSourceEnum);
156 
159  bool
160  GetGradientSourceIncludesFixed() const;
161 
164  bool
165  GetGradientSourceIncludesMoving() const;
166 
171  virtual void
172  Initialize() = 0;
173 
176  using NumberOfParametersType = unsigned int;
177 
181  MeasureType
182  GetValue() const override = 0;
183 
187  virtual void
188  GetDerivative(DerivativeType &) const = 0;
189 
192  void
193  GetValueAndDerivative(MeasureType & value, DerivativeType & derivative) const override = 0;
194 
199  GetNumberOfParameters() const override = 0;
200  virtual NumberOfParametersType
201  GetNumberOfLocalParameters() const = 0;
205  virtual void
206  SetParameters(ParametersType & params) = 0;
207 
209  virtual const ParametersType &
210  GetParameters() const = 0;
211 
214  virtual bool
215  HasLocalSupport() const = 0;
216 
223  virtual void
224  UpdateTransformParameters(const DerivativeType & derivative,
226 
232  MeasureType
233  GetCurrentValue() const;
234 
236 #if !defined(ITK_LEGACY_REMOVE)
237 
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 
252 protected:
254  ~ObjectToObjectMetricBaseTemplate() override = default;
255 
256  void
257  PrintSelf(std::ostream & os, Indent indent) const override;
258 
260  ObjectConstPointer m_FixedObject{};
261  ObjectConstPointer m_MovingObject{};
262 
263  GradientSourceEnum m_GradientSource{};
264 
266  mutable MeasureType m_Value{};
267 };
268 
271 } // end namespace itk
272 
273 #ifndef ITK_MANUAL_INSTANTIATION
274 # include "itkObjectToObjectMetricBase.hxx"
275 #endif
276 
277 #endif
MetricCategory
itk::OptimizerParameters
Class to hold and manage different parameter types used during optimization.
Definition: itkOptimizerParameters.h:36
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ObjectToObjectMetricBaseTemplate< TParametersValueType >::NumberOfParametersType
unsigned int NumberOfParametersType
Definition: itkObjectToObjectMetricBase.h:176
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory::IMAGE_METRIC
itk::ObjectToObjectMetricBaseTemplate::GetMetricCategory
virtual MetricCategoryEnum GetMetricCategory() const
Definition: itkObjectToObjectMetricBase.h:247
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory::OBJECT_METRIC
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
GradientSource
itk::ObjectToObjectMetricBaseTemplateEnums::GradientSource::GRADIENT_SOURCE_FIXED
itk::ObjectToObjectMetricBaseTemplateEnums::GradientSource::GRADIENT_SOURCE_BOTH
itk::ObjectToObjectMetricBaseTemplate< TParametersValueType >::DerivativeValueType
typename DerivativeType::ValueType DerivativeValueType
Definition: itkObjectToObjectMetricBase.h:113
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory::UNKNOWN_METRIC
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory
MetricCategory
Definition: itkObjectToObjectMetricBase.h:51
itkSingleValuedCostFunctionv4.h
itk::ObjectToObjectMetricBaseTemplateEnums
This class contains all the enum classes used by the ObjectToObjectMetricBaseTemplate class.
Definition: itkObjectToObjectMetricBase.h:31
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::ObjectToObjectMetricBaseTemplateEnums::GradientSource::GRADIENT_SOURCE_MOVING
itk::ObjectToObjectMetricBaseTemplateEnums::GradientSource
GradientSource
Definition: itkObjectToObjectMetricBase.h:40
itk::ObjectToObjectMetricBaseTemplate< TParametersValueType >::CoordinateRepresentationType
TParametersValueType CoordinateRepresentationType
Definition: itkObjectToObjectMetricBase.h:106
itk::CostFunctionTemplate< TParametersValueType >::ParametersValueType
TParametersValueType ParametersValueType
Definition: itkCostFunction.h:52
itk::ObjectToObjectMetricBaseTemplate
Base class for all object-to-object similarity metrics added in ITKv4.
Definition: itkObjectToObjectMetricBase.h:90
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory::POINT_SET_METRIC
itk::SingleValuedCostFunctionv4Template
This class is a base for a CostFunction that returns a single value.
Definition: itkSingleValuedCostFunctionv4.h:50
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::ObjectToObjectMetricBaseTemplateEnums::MetricCategory::MULTI_METRIC
itk::Array< TParametersValueType >::ValueType
TParametersValueType ValueType
Definition: itkArray.h:52
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::Array< TParametersValueType >
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itkTransformBase.h
itk::ObjectToObjectMetricBaseTemplate< TParametersValueType >::ObjectConstPointer
typename ObjectType::ConstPointer ObjectConstPointer
Definition: itkObjectToObjectMetricBase.h:121