ITK
6.0.0
Insight Toolkit
|
#include <itkRegularStepGradientDescentOptimizerv4.h>
Regular Step Gradient descent optimizer.
This optimizer is a variant of gradient descent that attempts to prevent it from taking steps that are too large. At each iteration, this optimizer will take a step along the direction of the metric derivative. Each time the direction of the derivative abruptly changes, the optimizer assumes that a local extrema has been passed and reacts by reducing the step length by a relaxation factor that is set to 0.5 by default. The default value for the initial step length is 1, and this value can only be changed manually via SetLearningRate() since this optimizer does not use the ScaleEstimator to automatically estimate the learning rate. Also note that unlike the previous version of RegularStepGradientDescentOptimizer, ITKv4 does not have a "maximize/minimize" option to modify the effect of the metric derivative. The assigned metric is assumed to return a parameter derivative result that "improves" the optimization.
Definition at line 47 of file itkRegularStepGradientDescentOptimizerv4.h.
Public Member Functions | |
void | EstimateLearningRate () override |
double | GetCurrentStepLength () const |
const char * | GetNameOfClass () const override |
void | StartOptimization (bool doOnlyInitialization=false) override |
virtual void | SetMinimumStepLength (TInternalComputationValueType _arg) |
virtual const TInternalComputationValueType & | GetMinimumStepLength () const |
virtual void | SetRelaxationFactor (TInternalComputationValueType _arg) |
virtual const TInternalComputationValueType & | GetRelaxationFactor () const |
virtual void | SetGradientMagnitudeTolerance (TInternalComputationValueType _arg) |
virtual const TInternalComputationValueType & | GetGradientMagnitudeTolerance () const |
virtual void | SetCurrentLearningRateRelaxation (MeasureType _arg) |
virtual const MeasureType & | GetCurrentLearningRateRelaxation () const |
Public Member Functions inherited from itk::GradientDescentOptimizerv4Template< double > | |
virtual void | EstimateLearningRate () |
virtual const double & | GetConvergenceValue () const |
const char * | GetNameOfClass () const override |
void | ResumeOptimization () override |
virtual void | SetConvergenceWindowSize (SizeValueType _arg) |
virtual void | SetMinimumConvergenceValue (double _arg) |
void | StartOptimization (bool doOnlyInitialization=false) override |
void | StopOptimization () override |
virtual void | SetLearningRate (double _arg) |
virtual const double & | GetLearningRate () const |
virtual void | SetMaximumStepSizeInPhysicalUnits (double _arg) |
virtual const double & | GetMaximumStepSizeInPhysicalUnits () const |
virtual void | SetDoEstimateLearningRateAtEachIteration (bool _arg) |
virtual const bool & | GetDoEstimateLearningRateAtEachIteration () const |
virtual void | DoEstimateLearningRateAtEachIterationOn () |
virtual void | SetDoEstimateLearningRateOnce (bool _arg) |
virtual const bool & | GetDoEstimateLearningRateOnce () const |
virtual void | DoEstimateLearningRateOnceOn () |
virtual void | SetReturnBestParametersAndValue (bool _arg) |
virtual const bool & | GetReturnBestParametersAndValue () const |
virtual void | ReturnBestParametersAndValueOn () |
Public Member Functions inherited from itk::GradientDescentOptimizerBasev4Template< double > | |
virtual const DerivativeType & | GetGradient () const |
const char * | GetNameOfClass () const override |
virtual const StopConditionObjectToObjectOptimizerEnum & | GetStopCondition () const |
StopConditionReturnStringType | GetStopConditionDescription () const override |
virtual void | ModifyGradientByLearningRateOverSubRange (const IndexRangeType &subrange)=0 |
virtual void | ModifyGradientByScalesOverSubRange (const IndexRangeType &subrange)=0 |
virtual void | ResumeOptimization ()=0 |
void | StartOptimization (bool doOnlyInitialization=false) override |
virtual void | StopOptimization () |
virtual void | ModifyGradientByScales () |
virtual void | ModifyGradientByLearningRate () |
Public Member Functions inherited from itk::ObjectToObjectOptimizerBaseTemplate< double > | |
virtual bool | CanUseScales () const |
virtual SizeValueType | GetCurrentIteration () const |
virtual const MeasureType & | GetCurrentMetricValue () const |
virtual const ParametersType & | GetCurrentPosition () const |
const char * | GetNameOfClass () const override |
virtual SizeValueType | GetNumberOfIterations () const |
virtual const ThreadIdType & | GetNumberOfWorkUnits () const |
virtual const ScalesType & | GetScales () const |
virtual const bool & | GetScalesAreIdentity () const |
bool | GetScalesInitialized () const |
virtual StopConditionReturnStringType | GetStopConditionDescription () const=0 |
virtual const MeasureType & | GetValue () const |
virtual const ScalesType & | GetWeights () const |
virtual const bool & | GetWeightsAreIdentity () const |
virtual void | SetNumberOfIterations (SizeValueType _arg) |
virtual void | SetNumberOfWorkUnits (ThreadIdType number) |
virtual void | SetScalesEstimator (ScalesEstimatorType *_arg) |
virtual void | SetWeights (ScalesType _arg) |
virtual void | StartOptimization (bool doOnlyInitialization=false) |
virtual void | SetMetric (MetricType *_arg) |
virtual MetricType * | GetModifiableMetric () |
virtual void | SetScales (const ScalesType &scales) |
virtual void | SetDoEstimateScales (bool _arg) |
virtual const bool & | GetDoEstimateScales () const |
virtual void | DoEstimateScalesOn () |
Public Member Functions inherited from itk::Object | |
unsigned long | AddObserver (const EventObject &event, Command *cmd) const |
unsigned long | AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const |
LightObject::Pointer | CreateAnother () const override |
virtual void | DebugOff () const |
virtual void | DebugOn () const |
Command * | GetCommand (unsigned long tag) |
bool | GetDebug () const |
MetaDataDictionary & | GetMetaDataDictionary () |
const MetaDataDictionary & | GetMetaDataDictionary () const |
virtual ModifiedTimeType | GetMTime () const |
const char * | GetNameOfClass () const override |
virtual const TimeStamp & | GetTimeStamp () const |
bool | HasObserver (const EventObject &event) const |
void | InvokeEvent (const EventObject &) |
void | InvokeEvent (const EventObject &) const |
virtual void | Modified () const |
void | Register () const override |
void | RemoveAllObservers () |
void | RemoveObserver (unsigned long tag) const |
void | SetDebug (bool debugFlag) const |
void | SetReferenceCount (int) override |
void | UnRegister () const noexcept override |
void | SetMetaDataDictionary (const MetaDataDictionary &rhs) |
void | SetMetaDataDictionary (MetaDataDictionary &&rrhs) |
virtual void | SetObjectName (std::string _arg) |
virtual const std::string & | GetObjectName () const |
Public Member Functions inherited from itk::LightObject | |
Pointer | Clone () const |
virtual Pointer | CreateAnother () const |
virtual void | Delete () |
virtual const char * | GetNameOfClass () const |
virtual int | GetReferenceCount () const |
void | Print (std::ostream &os, Indent indent=0) const |
virtual void | Register () const |
virtual void | SetReferenceCount (int) |
virtual void | UnRegister () const noexcept |
Static Public Member Functions | |
static Pointer | New () |
Static Public Member Functions inherited from itk::GradientDescentOptimizerv4Template< double > | |
static Pointer | New () |
Static Public Member Functions inherited from itk::Object | |
static bool | GetGlobalWarningDisplay () |
static void | GlobalWarningDisplayOff () |
static void | GlobalWarningDisplayOn () |
static Pointer | New () |
static void | SetGlobalWarningDisplay (bool val) |
Static Public Member Functions inherited from itk::LightObject | |
static void | BreakOnError () |
static Pointer | New () |
Protected Member Functions | |
void | AdvanceOneStep () override |
void | PrintSelf (std::ostream &os, Indent indent) const override |
RegularStepGradientDescentOptimizerv4 () | |
~RegularStepGradientDescentOptimizerv4 () override=default | |
void | ModifyGradientByScalesOverSubRange (const IndexRangeType &subrange) override |
void | ModifyGradientByLearningRateOverSubRange (const IndexRangeType &subrange) override |
Protected Member Functions inherited from itk::GradientDescentOptimizerv4Template< double > | |
virtual void | AdvanceOneStep () |
GradientDescentOptimizerv4Template () | |
void | ModifyGradientByLearningRateOverSubRange (const IndexRangeType &subrange) override |
void | ModifyGradientByScalesOverSubRange (const IndexRangeType &subrange) override |
void | PrintSelf (std::ostream &os, Indent indent) const override |
~GradientDescentOptimizerv4Template () override=default | |
Protected Member Functions inherited from itk::GradientDescentOptimizerBasev4Template< double > | |
void | PrintSelf (std::ostream &os, Indent indent) const override |
GradientDescentOptimizerBasev4Template () | |
~GradientDescentOptimizerBasev4Template () override=default | |
Protected Member Functions inherited from itk::ObjectToObjectOptimizerBaseTemplate< double > | |
void | PrintSelf (std::ostream &os, Indent indent) const override |
ObjectToObjectOptimizerBaseTemplate () | |
~ObjectToObjectOptimizerBaseTemplate () override | |
Protected Member Functions inherited from itk::Object | |
Object () | |
bool | PrintObservers (std::ostream &os, Indent indent) const |
void | PrintSelf (std::ostream &os, Indent indent) const override |
virtual void | SetTimeStamp (const TimeStamp &timeStamp) |
~Object () override | |
Protected Member Functions inherited from itk::LightObject | |
virtual LightObject::Pointer | InternalClone () const |
LightObject () | |
virtual void | PrintHeader (std::ostream &os, Indent indent) const |
virtual void | PrintSelf (std::ostream &os, Indent indent) const |
virtual void | PrintTrailer (std::ostream &os, Indent indent) const |
virtual | ~LightObject () |
Private Attributes | |
MeasureType | m_CurrentLearningRateRelaxation {} |
TInternalComputationValueType | m_GradientMagnitudeTolerance {} |
TInternalComputationValueType | m_MinimumStepLength {} |
TInternalComputationValueType | m_RelaxationFactor {} |
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::CompensatedSummationType = CompensatedSummation<InternalComputationValueType> |
Compensated summation type.
Definition at line 79 of file itkRegularStepGradientDescentOptimizerv4.h.
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::ConstPointer = SmartPointer<const Self> |
Definition at line 57 of file itkRegularStepGradientDescentOptimizerv4.h.
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::InternalComputationValueType = TInternalComputationValueType |
It should be possible to derive the internal computation type from the class object.
Definition at line 67 of file itkRegularStepGradientDescentOptimizerv4.h.
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::Pointer = SmartPointer<Self> |
Definition at line 56 of file itkRegularStepGradientDescentOptimizerv4.h.
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::Self = RegularStepGradientDescentOptimizerv4 |
Standard class type aliases.
Definition at line 54 of file itkRegularStepGradientDescentOptimizerv4.h.
using itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::Superclass = GradientDescentOptimizerv4Template<TInternalComputationValueType> |
Definition at line 55 of file itkRegularStepGradientDescentOptimizerv4.h.
|
protected |
Default constructor.
|
overrideprotecteddefault |
Destructor.
|
overrideprotectedvirtual |
Advance one Step following the gradient direction. Includes transform update.
Reimplemented from itk::GradientDescentOptimizerv4Template< double >.
|
overridevirtual |
Estimate the learning rate based on the current gradient.
Reimplemented from itk::GradientDescentOptimizerv4Template< double >.
|
virtual |
Set/Get current scale for learning rate.
double itk::RegularStepGradientDescentOptimizerv4< TInternalComputationValueType >::GetCurrentStepLength | ( | ) | const |
Get current gradient step value.
|
virtual |
Set/Get gradient magnitude tolerance value for convergence checking.
|
virtual |
Minimum step length (learning rate) value for convergence checking. When the local minima is passed by taking a large step, the step length is adjusted (decreased) by the relaxation factor, so that smaller steps are taken towards the minimum point (convergence). When the step length value reaches a small value, it would be treated as converged.
The default value is set to 1e-4 to pass all tests.
|
overridevirtual |
Reimplemented from itk::Object.
|
virtual |
Set/Get relaxation factor value.
|
overrideprotectedvirtual |
Modify the input gradient over a given index range.
Implements itk::GradientDescentOptimizerBasev4Template< double >.
|
overrideprotectedvirtual |
Modify the input gradient over a given index range.
Implements itk::GradientDescentOptimizerBasev4Template< double >.
|
static |
New macro for creation of through a Smart Pointer.
|
overrideprotectedvirtual |
Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.
Reimplemented from itk::Object.
|
virtual |
Set/Get current scale for learning rate.
|
virtual |
Set/Get gradient magnitude tolerance value for convergence checking.
|
virtual |
Minimum step length (learning rate) value for convergence checking. When the local minima is passed by taking a large step, the step length is adjusted (decreased) by the relaxation factor, so that smaller steps are taken towards the minimum point (convergence). When the step length value reaches a small value, it would be treated as converged.
The default value is set to 1e-4 to pass all tests.
|
virtual |
Set/Get relaxation factor value.
|
overridevirtual |
Start and run the optimization.
Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.
|
private |
Definition at line 152 of file itkRegularStepGradientDescentOptimizerv4.h.
|
private |
Definition at line 150 of file itkRegularStepGradientDescentOptimizerv4.h.
|
private |
Definition at line 148 of file itkRegularStepGradientDescentOptimizerv4.h.
|
private |
Definition at line 146 of file itkRegularStepGradientDescentOptimizerv4.h.