ITK  4.2.0
Insight Segmentation and Registration Toolkit
itkPowellOptimizer.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 
25 namespace itk
26 {
61 class ITK_EXPORT PowellOptimizer:
63 {
64 public:
70 
73 
75  itkNewMacro(Self);
76 
79 
83 
85  itkSetMacro(Maximize, bool);
86  itkBooleanMacro(Maximize);
87  itkGetConstReferenceMacro(Maximize, bool);
89 
91  itkSetMacro(MaximumIteration, unsigned int);
92  itkGetConstReferenceMacro(MaximumIteration, unsigned int);
94 
96  itkSetMacro(MaximumLineIteration, unsigned int);
97  itkGetConstMacro(MaximumLineIteration, unsigned int);
99 
102  itkSetMacro(StepLength, double);
103  itkGetConstReferenceMacro(StepLength, double);
105 
108  itkSetMacro(StepTolerance, double);
109  itkGetConstReferenceMacro(StepTolerance, double);
111 
115  itkSetMacro(ValueTolerance, double);
116  itkGetConstReferenceMacro(ValueTolerance, double);
118 
120  itkGetConstReferenceMacro(CurrentCost, MeasureType);
121  MeasureType GetValue() const { return this->GetCurrentCost(); }
123 
125  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
126 
128  itkGetConstReferenceMacro(CurrentLineIteration, unsigned int);
129 
131  void StartOptimization();
132 
136  void StopOptimization()
137  { m_Stop = true; }
138 
139  itkGetConstReferenceMacro(CatchGetValueException, bool);
140  itkSetMacro(CatchGetValueException, bool);
141 
142  itkGetConstReferenceMacro(MetricWorstPossibleValue, double);
143  itkSetMacro(MetricWorstPossibleValue, double);
144 
145  const std::string GetStopConditionDescription() const;
146 
147 protected:
148  PowellOptimizer();
150  virtual ~PowellOptimizer();
151  void PrintSelf(std::ostream & os, Indent indent) const;
152 
153  itkSetMacro(CurrentCost, double);
154 
157  void SetLine(const ParametersType & origin,
158  const vnl_vector< double > & direction);
159 
163  double GetLineValue(double x) const;
164 
165  double GetLineValue(double x, ParametersType & tempCoord) const;
166 
169  void SetCurrentLinePoint(double x, double fx);
170 
173  void Swap(double *a, double *b) const;
174 
177  void Shift(double *a, double *b, double *c, double d) const;
178 
188  virtual void LineBracket(double *ax, double *bx, double *cx,
189  double *fa, double *fb, double *fc);
190 
191  virtual void LineBracket(double *ax, double *bx, double *cx,
192  double *fa, double *fb, double *fc,
193  ParametersType & tempCoord);
194 
200  virtual void BracketedLineOptimize(double ax, double bx, double cx,
201  double fa, double fb, double fc,
202  double *extX, double *extVal);
203 
204  virtual void BracketedLineOptimize(double ax, double bx, double cx,
205  double fa, double fb, double fc,
206  double *extX, double *extVal,
207  ParametersType & tempCoord);
208 
209  itkGetMacro(SpaceDimension, unsigned int);
210  void SetSpaceDimension(unsigned int dim)
211  {
212  this->m_SpaceDimension = dim;
213  this->m_LineDirection.set_size(dim);
214  this->m_LineOrigin.set_size(dim);
215  this->m_CurrentPosition.set_size(dim);
216  this->Modified();
217  }
218 
219  itkSetMacro(CurrentIteration, unsigned int);
220 
221  itkGetMacro(Stop, bool);
222  itkSetMacro(Stop, bool);
223 private:
224  unsigned int m_SpaceDimension;
225 
227  unsigned int m_CurrentIteration;
229 
231  unsigned int m_MaximumIteration;
233 
236 
239 
241  double m_StepLength;
243 
246 
248 
251 
256  bool m_Stop;
257 
258  std::ostringstream m_StopConditionDescription;
259 }; // end of class
260 } // end of namespace itk
261 
262 #endif
263