ITK  5.2.0
Insight Toolkit
itkPowellOptimizer.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  * http://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 itkPowellOptimizer_h
19 #define itkPowellOptimizer_h
20 
21 #include "itkVector.h"
22 #include "itkMatrix.h"
24 #include "ITKOptimizersExport.h"
25 
26 namespace itk
27 {
62 class ITKOptimizers_EXPORT PowellOptimizer : public SingleValuedNonLinearOptimizer
63 {
64 public:
70 
72 
74  itkNewMacro(Self);
75 
78 
82 
84  itkSetMacro(Maximize, bool);
85  itkBooleanMacro(Maximize);
86  itkGetConstReferenceMacro(Maximize, bool);
88 
90  itkSetMacro(MaximumIteration, unsigned int);
91  itkGetConstReferenceMacro(MaximumIteration, unsigned int);
93 
95  itkSetMacro(MaximumLineIteration, unsigned int);
96  itkGetConstMacro(MaximumLineIteration, unsigned int);
98 
101  itkSetMacro(StepLength, double);
102  itkGetConstReferenceMacro(StepLength, double);
104 
107  itkSetMacro(StepTolerance, double);
108  itkGetConstReferenceMacro(StepTolerance, double);
110 
114  itkSetMacro(ValueTolerance, double);
115  itkGetConstReferenceMacro(ValueTolerance, double);
117 
119  itkGetConstReferenceMacro(CurrentCost, MeasureType);
120  MeasureType
121  GetValue() const
122  {
123  return this->GetCurrentCost();
124  }
126 
128  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
129 
131  itkGetConstReferenceMacro(CurrentLineIteration, unsigned int);
132 
134  void
135  StartOptimization() override;
136 
140  void
142  {
143  m_Stop = true;
144  }
145 
146  itkGetConstReferenceMacro(CatchGetValueException, bool);
147  itkSetMacro(CatchGetValueException, bool);
148 
149  itkGetConstReferenceMacro(MetricWorstPossibleValue, double);
150  itkSetMacro(MetricWorstPossibleValue, double);
151 
152  const std::string
153  GetStopConditionDescription() const override;
154 
155 protected:
156  PowellOptimizer();
158  ~PowellOptimizer() override;
159  void
160  PrintSelf(std::ostream & os, Indent indent) const override;
161 
162  itkSetMacro(CurrentCost, double);
163 
166  void
167  SetLine(const ParametersType & origin, const vnl_vector<double> & direction);
168 
172  double
173  GetLineValue(double x) const;
174 
175  double
176  GetLineValue(double x, ParametersType & tempCoord) const;
177 
180  void
181  SetCurrentLinePoint(double x, double fx);
182 
185  void
186  Swap(double * a, double * b) const;
187 
190  void
191  Shift(double * a, double * b, double * c, double d) const;
192 
202  virtual void
203  LineBracket(double * x1, double * x2, double * x3, double * f1, double * f2, double * f3);
204 
205  virtual void
206  LineBracket(double * x1, double * x2, double * x3, double * f1, double * f2, double * f3, ParametersType & tempCoord);
207 
213  virtual void
214  BracketedLineOptimize(double ax,
215  double bx,
216  double cx,
217  double fa,
218  double functionValueOfb,
219  double fc,
220  double * extX,
221  double * extVal);
222 
223  virtual void
224  BracketedLineOptimize(double ax,
225  double bx,
226  double cx,
227  double fa,
228  double functionValueOfb,
229  double fc,
230  double * extX,
231  double * extVal,
232  ParametersType & tempCoord);
233 
234  itkGetMacro(SpaceDimension, unsigned int);
235  void
236  SetSpaceDimension(unsigned int dim)
237  {
238  this->m_SpaceDimension = dim;
239  this->m_LineDirection.set_size(dim);
240  this->m_LineOrigin.set_size(dim);
241  this->m_CurrentPosition.set_size(dim);
242  this->Modified();
243  }
244 
245  itkSetMacro(CurrentIteration, unsigned int);
246 
247  itkGetMacro(Stop, bool);
248  itkSetMacro(Stop, bool);
249 
250 private:
251  unsigned int m_SpaceDimension;
252 
254  unsigned int m_CurrentIteration;
256 
258  unsigned int m_MaximumIteration;
260 
263 
266 
268  double m_StepLength;
270 
272  vnl_vector<double> m_LineDirection;
273 
275 
277  MeasureType m_CurrentCost;
278 
283  bool m_Stop;
284 
285  std::ostringstream m_StopConditionDescription;
286 }; // end of class
287 } // end of namespace itk
288 
289 #endif
itk::PowellOptimizer::m_LineOrigin
ParametersType m_LineOrigin
Definition: itkPowellOptimizer.h:271
itk::OptimizerParameters< double >
itk::SingleValuedNonLinearOptimizer
This class is a base for the Optimization methods that optimize a single valued function.
Definition: itkSingleValuedNonLinearOptimizer.h:35
itkMatrix.h
itk::PowellOptimizer::m_MetricWorstPossibleValue
double m_MetricWorstPossibleValue
Definition: itkPowellOptimizer.h:262
itk::PowellOptimizer::m_CurrentIteration
unsigned int m_CurrentIteration
Definition: itkPowellOptimizer.h:254
itk::PowellOptimizer::SetSpaceDimension
void SetSpaceDimension(unsigned int dim)
Definition: itkPowellOptimizer.h:236
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::PowellOptimizer::m_CurrentCost
MeasureType m_CurrentCost
Definition: itkPowellOptimizer.h:277
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::SingleValuedCostFunction
This class is a base for the CostFunctions returning a single value.
Definition: itkSingleValuedCostFunction.h:34
itk::PowellOptimizer::m_MaximumLineIteration
unsigned int m_MaximumLineIteration
Definition: itkPowellOptimizer.h:259
itk::PowellOptimizer::m_MaximumIteration
unsigned int m_MaximumIteration
Definition: itkPowellOptimizer.h:258
itk::PowellOptimizer::m_Stop
bool m_Stop
Definition: itkPowellOptimizer.h:283
itk::PowellOptimizer
Implements Powell optimization using Brent line search.
Definition: itkPowellOptimizer.h:62
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::PowellOptimizer::m_CurrentLineIteration
unsigned int m_CurrentLineIteration
Definition: itkPowellOptimizer.h:255
itkVector.h
itk::PowellOptimizer::m_ValueTolerance
double m_ValueTolerance
Definition: itkPowellOptimizer.h:274
itk::PowellOptimizer::m_LineDirection
vnl_vector< double > m_LineDirection
Definition: itkPowellOptimizer.h:272
itk::PowellOptimizer::m_StepTolerance
double m_StepTolerance
Definition: itkPowellOptimizer.h:269
itk::PowellOptimizer::GetValue
MeasureType GetValue() const
Definition: itkPowellOptimizer.h:121
itk::PowellOptimizer::StopOptimization
void StopOptimization()
Definition: itkPowellOptimizer.h:141
itk::SingleValuedNonLinearOptimizer::ParametersType
Superclass::ParametersType ParametersType
Definition: itkSingleValuedNonLinearOptimizer.h:54
itk::PowellOptimizer::m_StepLength
double m_StepLength
Definition: itkPowellOptimizer.h:268
itkSingleValuedNonLinearOptimizer.h
itk::PowellOptimizer::m_CatchGetValueException
bool m_CatchGetValueException
Definition: itkPowellOptimizer.h:261
itk::PowellOptimizer::m_StopConditionDescription
std::ostringstream m_StopConditionDescription
Definition: itkPowellOptimizer.h:285
itk::PowellOptimizer::m_Maximize
bool m_Maximize
Definition: itkPowellOptimizer.h:265
itk::PowellOptimizer::m_SpaceDimension
unsigned int m_SpaceDimension
Definition: itkPowellOptimizer.h:248