go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itk::StandardGradientDescentOptimizer Class Reference

#include <itkStandardGradientDescentOptimizer.h>

Detailed Description

This class implements a gradient descent optimizer with a decaying gain.

If $C(x)$ is a costfunction that has to be minimised, the following iterative algorithm is used to find the optimal parameters $x$:

\[ x(k+1) = x(k) - a(k) dC/dx \]

The gain $a(k)$ at each iteration $k$ is defined by:

\[ a(k) =  a / (A + k + 1)^\alpha \]

.

It is very suitable to be used in combination with a stochastic estimate of the gradient $dC/dx$. For example, in image registration problems it is often advantageous to compute the metric derivative ( $dC/dx$) on a new set of randomly selected image samples in each iteration. You may set the parameter NewSamplesEveryIteration to "true" to achieve this effect. For more information on this strategy, you may have a look at:

S. Klein, M. Staring, J.P.W. Pluim, "Comparison of gradient approximation techniques for optimisation of mutual information in nonrigid registration", in: SPIE Medical Imaging: Image Processing, Editor(s): J.M. Fitzpatrick, J.M. Reinhardt, SPIE press, 2005, vol. 5747, Proceedings of SPIE, pp. 192-203.

Or:

S. Klein, M. Staring, J.P.W. Pluim, "Evaluation of Optimization Methods for Nonrigid Medical Image Registration using Mutual Information and B-Splines" IEEE Transactions on Image Processing, 2007, nr. 16(12), December.

This class also serves as a base class for other GradientDescent type algorithms, like the AcceleratedGradientDescentOptimizer.

See also
StandardGradientDescent, AcceleratedGradientDescentOptimizer

Definition at line 65 of file itkStandardGradientDescentOptimizer.h.

Inheritance diagram for itk::StandardGradientDescentOptimizer:

Public Types

using ConstPointer = SmartPointer<const Self>
 
using Pointer = SmartPointer<Self>
 
using ScaledCostFunctionPointer
 
using ScaledCostFunctionType
 
using ScalesType
 
using Self = StandardGradientDescentOptimizer
 
enum  StopConditionType
 
using Superclass = GradientDescentOptimizer2
 
- Public Types inherited from itk::GradientDescentOptimizer2
using ConstPointer = SmartPointer<const Self>
 
using Pointer = SmartPointer<Self>
 
using ScaledCostFunctionPointer
 
using ScaledCostFunctionType
 
using ScalesType
 
using Self = GradientDescentOptimizer2
 
enum  StopConditionType { MaximumNumberOfIterations , MetricError , MinimumStepSize }
 
using Superclass = ScaledSingleValuedNonLinearOptimizer
 
- Public Types inherited from itk::ScaledSingleValuedNonLinearOptimizer
using ConstPointer = SmartPointer<const Self>
 
using Pointer = SmartPointer<Self>
 
using ScaledCostFunctionPointer = ScaledCostFunctionType::Pointer
 
using ScaledCostFunctionType = ScaledSingleValuedCostFunction
 
using ScalesType = NonLinearOptimizer::ScalesType
 
using Self = ScaledSingleValuedNonLinearOptimizer
 
using Superclass = SingleValuedNonLinearOptimizer
 

Public Member Functions

void AdvanceOneStep () override
 
virtual const char * GetClassName () const
 
virtual double GetCurrentTime () const
 
virtual double GetInitialTime () const
 
virtual double GetParam_A () const
 
virtual double GetParam_a () const
 
virtual double GetParam_alpha () const
 
 ITK_DISALLOW_COPY_AND_MOVE (StandardGradientDescentOptimizer)
 
virtual void ResetCurrentTimeToInitialTime ()
 
virtual void SetInitialTime (double _arg)
 
virtual void SetParam_A (double _arg)
 
virtual void SetParam_a (double _arg)
 
virtual void SetParam_alpha (double _arg)
 
void StartOptimization () override
 
- Public Member Functions inherited from itk::GradientDescentOptimizer2
virtual unsigned int GetCurrentIteration () const
 
virtual const DerivativeType & GetGradient ()
 
virtual const doubleGetLearningRate ()
 
virtual const unsigned long & GetNumberOfIterations ()
 
virtual const DerivativeType & GetSearchDirection ()
 
virtual const StopConditionTypeGetStopCondition ()
 
virtual const doubleGetValue ()
 
 ITK_DISALLOW_COPY_AND_MOVE (GradientDescentOptimizer2)
 
virtual void MetricErrorResponse (ExceptionObject &err)
 
virtual void ResumeOptimization ()
 
virtual void SetLearningRate (double _arg)
 
virtual void SetNumberOfIterations (unsigned long _arg)
 
void StartOptimization () override
 
virtual void StopOptimization ()
 
- Public Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer
const ParametersType & GetCurrentPosition () const override
 
virtual bool GetMaximize () const
 
virtual const ScaledCostFunctionTypeGetScaledCostFunction ()
 
virtual const ParametersType & GetScaledCurrentPosition ()
 
bool GetUseScales () const
 
virtual void InitializeScales ()
 
 ITK_DISALLOW_COPY_AND_MOVE (ScaledSingleValuedNonLinearOptimizer)
 
virtual void MaximizeOff ()
 
virtual void MaximizeOn ()
 
void SetCostFunction (CostFunctionType *costFunction) override
 
virtual void SetMaximize (bool _arg)
 
virtual void SetUseScales (bool arg)
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::GradientDescentOptimizer2
static Pointer New ()
 
- Static Public Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer
static Pointer New ()
 

Protected Member Functions

virtual double Compute_a (double k) const
 
 StandardGradientDescentOptimizer ()
 
virtual void UpdateCurrentTime ()
 
 ~StandardGradientDescentOptimizer () override=default
 
- Protected Member Functions inherited from itk::GradientDescentOptimizer2
 GradientDescentOptimizer2 ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~GradientDescentOptimizer2 () override=default
 
- Protected Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer
virtual void GetScaledDerivative (const ParametersType &parameters, DerivativeType &derivative) const
 
virtual MeasureType GetScaledValue (const ParametersType &parameters) const
 
virtual void GetScaledValueAndDerivative (const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ScaledSingleValuedNonLinearOptimizer ()
 
void SetCurrentPosition (const ParametersType &param) override
 
virtual void SetScaledCurrentPosition (const ParametersType &parameters)
 
 ~ScaledSingleValuedNonLinearOptimizer () override=default
 

Protected Attributes

double m_CurrentTime { 0.0 }
 
bool m_UseConstantStep { false }
 
- Protected Attributes inherited from itk::GradientDescentOptimizer2
DerivativeType m_Gradient {}
 
DerivativeType m_SearchDirection {}
 
StopConditionType m_StopCondition { MaximumNumberOfIterations }
 
- Protected Attributes inherited from itk::ScaledSingleValuedNonLinearOptimizer
ScaledCostFunctionPointer m_ScaledCostFunction {}
 
ParametersType m_ScaledCurrentPosition {}
 

Private Attributes

double m_InitialTime { 0.0 }
 
double m_Param_A { 1.0 }
 
double m_Param_a { 1.0 }
 
double m_Param_alpha { 0.602 }
 

Member Typedef Documentation

◆ ConstPointer

◆ Pointer

◆ ScaledCostFunctionPointer

◆ ScaledCostFunctionType

◆ ScalesType

◆ Self

◆ Superclass

Member Enumeration Documentation

◆ StopConditionType

Codes of stopping conditions The MinimumStepSize stopcondition never occurs, but may be implemented in inheriting classes

Definition at line 83 of file itkGradientDescentOptimizer2.h.

Constructor & Destructor Documentation

◆ StandardGradientDescentOptimizer()

itk::StandardGradientDescentOptimizer::StandardGradientDescentOptimizer ( )
protected

◆ ~StandardGradientDescentOptimizer()

itk::StandardGradientDescentOptimizer::~StandardGradientDescentOptimizer ( )
overrideprotecteddefault

Member Function Documentation

◆ AdvanceOneStep()

void itk::StandardGradientDescentOptimizer::AdvanceOneStep ( )
overridevirtual

Sets a new LearningRate before calling the Superclass' implementation, and updates the current time.

Reimplemented from itk::GradientDescentOptimizer2.

◆ Compute_a()

virtual double itk::StandardGradientDescentOptimizer::Compute_a ( double k) const
protectedvirtual

Function to compute the parameter at time/iteration k.

◆ GetClassName()

◆ GetCurrentTime()

virtual double itk::StandardGradientDescentOptimizer::GetCurrentTime ( ) const
virtual

Get the current time. This equals the CurrentIteration in this base class but may be different in inheriting classes, such as the AccelerateGradientDescent

◆ GetInitialTime()

virtual double itk::StandardGradientDescentOptimizer::GetInitialTime ( ) const
virtual

◆ GetParam_A()

virtual double itk::StandardGradientDescentOptimizer::GetParam_A ( ) const
virtual

◆ GetParam_a()

virtual double itk::StandardGradientDescentOptimizer::GetParam_a ( ) const
virtual

◆ GetParam_alpha()

virtual double itk::StandardGradientDescentOptimizer::GetParam_alpha ( ) const
virtual

◆ ITK_DISALLOW_COPY_AND_MOVE()

itk::StandardGradientDescentOptimizer::ITK_DISALLOW_COPY_AND_MOVE ( StandardGradientDescentOptimizer )

◆ New()

static Pointer itk::StandardGradientDescentOptimizer::New ( )
static

Method for creation through the object factory.

◆ ResetCurrentTimeToInitialTime()

virtual void itk::StandardGradientDescentOptimizer::ResetCurrentTimeToInitialTime ( )
inlinevirtual

Set the current time to the initial time. This can be useful to 'reset' the optimisation, for example if you changed the cost function while optimisation. Be careful with this function.

Definition at line 130 of file itkStandardGradientDescentOptimizer.h.

◆ SetInitialTime()

virtual void itk::StandardGradientDescentOptimizer::SetInitialTime ( double _arg)
virtual

Set/Get the initial time. Should be >=0. This function is superfluous, since Param_A does effectively the same. However, in inheriting classes, like the AcceleratedGradientDescent the initial time may have a different function than Param_A. Default: 0.0

◆ SetParam_A()

virtual void itk::StandardGradientDescentOptimizer::SetParam_A ( double _arg)
virtual

Set/Get A.

◆ SetParam_a()

virtual void itk::StandardGradientDescentOptimizer::SetParam_a ( double _arg)
virtual

Set/Get a.

◆ SetParam_alpha()

virtual void itk::StandardGradientDescentOptimizer::SetParam_alpha ( double _arg)
virtual

Set/Get alpha.

◆ StartOptimization()

void itk::StandardGradientDescentOptimizer::StartOptimization ( )
override

Set current time to 0 and call superclass' implementation.

◆ UpdateCurrentTime()

virtual void itk::StandardGradientDescentOptimizer::UpdateCurrentTime ( )
protectedvirtual

Function to update the current time This function just increments the CurrentTime by 1. Inheriting functions may implement something smarter, for example, dependent on the progress

Reimplemented in itk::AdaptiveStepsizeOptimizer, itk::AdaptiveStochasticGradientDescentOptimizer, and itk::PreconditionedASGDOptimizer.

Field Documentation

◆ m_CurrentTime

double itk::StandardGradientDescentOptimizer::m_CurrentTime { 0.0 }
protected

The current time, which serves as input for Compute_a

Definition at line 152 of file itkStandardGradientDescentOptimizer.h.

◆ m_InitialTime

double itk::StandardGradientDescentOptimizer::m_InitialTime { 0.0 }
private

Settings

Definition at line 164 of file itkStandardGradientDescentOptimizer.h.

◆ m_Param_A

double itk::StandardGradientDescentOptimizer::m_Param_A { 1.0 }
private

Definition at line 160 of file itkStandardGradientDescentOptimizer.h.

◆ m_Param_a

double itk::StandardGradientDescentOptimizer::m_Param_a { 1.0 }
private

Parameters, as described by Spall.

Definition at line 159 of file itkStandardGradientDescentOptimizer.h.

◆ m_Param_alpha

double itk::StandardGradientDescentOptimizer::m_Param_alpha { 0.602 }
private

Definition at line 161 of file itkStandardGradientDescentOptimizer.h.

◆ m_UseConstantStep

bool itk::StandardGradientDescentOptimizer::m_UseConstantStep { false }
protected

Constant step size or others, different value of k.

Definition at line 155 of file itkStandardGradientDescentOptimizer.h.



Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo