/*========================================================================= * * Copyright NumFOCUS * * Licensed under the Apache License, Version 2.0 (the "License"); * you may not use this file except in compliance with the License. * You may obtain a copy of the License at * * https://www.apache.org/licenses/LICENSE-2.0.txt * * Unless required by applicable law or agreed to in writing, software * distributed under the License is distributed on an "AS IS" BASIS, * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. * See the License for the specific language governing permissions and * limitations under the License. * *=========================================================================*/ #ifndef itkExhaustiveOptimizer_h #define itkExhaustiveOptimizer_h #include "itkIntTypes.h" #include "itkSingleValuedNonLinearOptimizer.h" #include "ITKOptimizersExport.h" namespace itk { /** \class ExhaustiveOptimizer * \brief Optimizer that fully samples a grid on the parametric space. * * This optimizer is equivalent to an exhaustive search in a discrete grid * defined over the parametric space. The grid is centered on the initial * position. The subdivisions of the grid along each one of the dimensions * of the parametric space is defined by an array of number of steps. * * A typical use is to plot the metric space to get an idea of how noisy it * is. An example is given below, where it is desired to plot the metric * space with respect to translations along x, y and z in a 3D registration * application: * Here it is assumed that the transform is TranslationTransform. * \code OptimizerType::StepsType steps( m_Transform->GetNumberOfParameters() ); steps[0] = 10; steps[1] = 10; steps[2] = 10; m_Optimizer->SetNumberOfSteps( steps ); m_Optimizer->SetStepLength( 2 ); \endcode * * The optimizer throws IterationEvents after every iteration. We use this to plot * the metric space in an image as follows: * \code if( itk::IterationEvent().CheckEvent(& event ) ) { IndexType index; index[0] = m_Optimizer->GetCurrentIndex()[0]; index[1] = m_Optimizer->GetCurrentIndex()[1]; index[2] = m_Optimizer->GetCurrentIndex()[2]; image->SetPixel( index, m_Optimizer->GetCurrentValue() ); } \endcode * * The image size is expected to be 11 x 11 x 11. * * If you wish to use different step lengths along each parametric axis, * you can use the SetScales() method. This accepts an array, each element * represents the number of subdivisions per step length. For instance scales * of [0.5 1 4] along with a step length of 2 will cause the optimizer * to search the metric space on a grid with x,y,z spacing of [1 2 8]. * * The number of samples for each dimension of the parameter grid are * influenced by both the scales and the number of steps along each * dimension: * * parameter_samples[d] = stepLength*(2*numberOfSteps[d]+1)*scaling[d] * * start_parameter[d] = - stepLength * scaling[d] * numberOfSteps[d] * end_parameter[d] = + stepLength * scaling[d] * numberOfSteps[d] * * \ingroup Numerics Optimizers * \ingroup ITKOptimizers * * \sphinx * \sphinxexample{Numerics/Optimizers/ExhaustiveOptimizer,Exhaustive Optimizer} * \endsphinx */ class ITKOptimizers_EXPORT ExhaustiveOptimizer : public SingleValuedNonLinearOptimizer { public: ITK_DISALLOW_COPY_AND_MOVE(ExhaustiveOptimizer); /** Standard "Self" type alias. */ using Self = ExhaustiveOptimizer; using Superclass = SingleValuedNonLinearOptimizer; using Pointer = SmartPointer; using ConstPointer = SmartPointer; using StepsType = Array; /** Method for creation through the object factory. */ itkNewMacro(Self); /** \see LightObject::GetNameOfClass() */ itkOverrideGetNameOfClassMacro(ExhaustiveOptimizer); /** Start optimization. */ void StartOptimization() override; void StartWalking(); /** Resume optimization. */ void ResumeWalking(); /** Stop optimization. */ void StopWalking(); itkSetMacro(StepLength, double); itkSetMacro(NumberOfSteps, StepsType); itkGetConstReferenceMacro(StepLength, double); itkGetConstReferenceMacro(NumberOfSteps, StepsType); itkGetConstReferenceMacro(CurrentValue, MeasureType); itkGetConstReferenceMacro(MaximumMetricValue, MeasureType); itkGetConstReferenceMacro(MinimumMetricValue, MeasureType); itkGetConstReferenceMacro(MinimumMetricValuePosition, ParametersType); itkGetConstReferenceMacro(MaximumMetricValuePosition, ParametersType); itkGetConstReferenceMacro(CurrentIndex, ParametersType); itkGetConstReferenceMacro(MaximumNumberOfIterations, SizeValueType); /** Get the reason for termination */ const std::string GetStopConditionDescription() const override; protected: ExhaustiveOptimizer(); ~ExhaustiveOptimizer() override = default; void PrintSelf(std::ostream & os, Indent indent) const override; /** Advance to the next grid position. */ void AdvanceOneStep(); void IncrementIndex(ParametersType & newPosition); protected: MeasureType m_CurrentValue{ 0 }; StepsType m_NumberOfSteps{}; SizeValueType m_CurrentIteration{ 0 }; bool m_Stop{ false }; unsigned int m_CurrentParameter{ 0 }; double m_StepLength{ 1.0 }; ParametersType m_CurrentIndex{}; SizeValueType m_MaximumNumberOfIterations{ 1 }; MeasureType m_MaximumMetricValue{ itk::NumericTraits::max() }; MeasureType m_MinimumMetricValue{ itk::NumericTraits::Zero }; ParametersType m_MinimumMetricValuePosition{}; ParametersType m_MaximumMetricValuePosition{}; private: std::ostringstream m_StopConditionDescription{}; }; } // end namespace itk #endif