00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #ifndef __itkPowellOptimizer_h
00018 #define __itkPowellOptimizer_h
00019
00020 #include <itkVector.h>
00021 #include <itkMatrix.h>
00022 #include <itkSingleValuedNonLinearOptimizer.h>
00023
00024 namespace itk
00025 {
00026
00054 class ITK_EXPORT PowellOptimizer:
00055 public SingleValuedNonLinearOptimizer
00056 {
00057 public:
00059 typedef PowellOptimizer Self ;
00060 typedef SingleValuedNonLinearOptimizer Superclass;
00061 typedef SmartPointer<Self> Pointer;
00062 typedef SmartPointer<const Self> ConstPointer;
00063
00064 typedef SingleValuedNonLinearOptimizer::ParametersType
00065 ParametersType;
00066
00068 itkNewMacro(Self);
00069
00071 itkTypeMacro(PowellOptimizer, SingleValuedNonLinearOptimizer );
00072
00074 typedef SingleValuedCostFunction CostFunctionType;
00075 typedef CostFunctionType::Pointer CostFunctionPointer;
00076
00078 itkSetMacro( Maximize, bool );
00079 itkGetConstReferenceMacro( Maximize, bool );
00080
00082 itkSetMacro( MaximumIteration, unsigned int );
00083 itkGetConstReferenceMacro( MaximumIteration, unsigned int );
00084
00086 itkSetMacro(MaximumLineIteration, unsigned int);
00087 itkGetConstMacro(MaximumLineIteration, unsigned int);
00088
00091 itkSetMacro( StepLength, double ) ;
00092 itkGetConstReferenceMacro( StepLength, double );
00093
00096 itkSetMacro( StepTolerance, double ) ;
00097 itkGetConstReferenceMacro( StepTolerance, double );
00098
00102 itkSetMacro( ValueTolerance, double ) ;
00103 itkGetConstReferenceMacro( ValueTolerance, double );
00104
00106 itkGetConstReferenceMacro( CurrentCost, MeasureType );
00107 MeasureType GetValue() const { return this->GetCurrentCost(); }
00108
00110 itkGetConstReferenceMacro( CurrentIteration, unsigned int);
00111
00113 itkGetConstReferenceMacro( CurrentLineIteration, unsigned int);
00114
00116 void StartOptimization() ;
00117
00121 void StopOptimization()
00122 { m_Stop = true ; }
00123
00124 protected:
00125 PowellOptimizer() ;
00126 PowellOptimizer(const PowellOptimizer&) ;
00127 virtual ~PowellOptimizer() ;
00128 void PrintSelf(std::ostream& os, Indent indent) const;
00129
00130 itkSetMacro(CurrentCost, double);
00131
00134 void SetLine(const ParametersType & origin,
00135 const vnl_vector<double> & direction);
00136
00140 double GetLineValue(double x) const;
00141
00144 void SetCurrentLinePoint(double x, double fx);
00145
00148 void Swap(double *a, double *b) const;
00149
00152 void Shift(double *a, double *b, double *c, double d) const;
00153
00160 virtual void LineBracket(double *ax, double *bx, double *cx,
00161 double *fa, double *fb, double *fc);
00162
00168 virtual void BracketedLineOptimize(double ax, double bx, double cx,
00169 double fa, double fb, double fc,
00170 double * extX, double * extVal);
00171
00172 itkGetMacro(SpaceDimension, unsigned int);
00173 itkSetMacro(SpaceDimension, unsigned int);
00174
00175 itkSetMacro(CurrentIteration, unsigned int);
00176
00177 itkGetMacro(Stop, bool);
00178 itkSetMacro(Stop, bool);
00179
00180 private:
00181
00182 unsigned int m_SpaceDimension;
00183
00185 unsigned int m_CurrentIteration ;
00186 unsigned int m_CurrentLineIteration ;
00187
00189 unsigned int m_MaximumIteration ;
00190 unsigned int m_MaximumLineIteration ;
00191
00193 bool m_Maximize;
00194
00196 double m_StepLength ;
00197 double m_StepTolerance ;
00198
00199 ParametersType m_LineOrigin;
00200 vnl_vector<double> m_LineDirection;
00201
00202 double m_ValueTolerance;
00203
00205 MeasureType m_CurrentCost;
00206
00211 bool m_Stop ;
00212
00213 } ;
00214
00215 }
00216
00217 #endif