ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkLBFGSBOptimizer.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 itkLBFGSBOptimizer_h
19 #define itkLBFGSBOptimizer_h
20 
21 #include "itkIntTypes.h"
23 #include "ITKOptimizersExport.h"
24 
25 namespace itk
26 {
27 /* Necessary forward declaration see below for definition */
35 class ITK_FORWARD_EXPORT LBFGSBOptimizerHelper;
36 
63 class ITKOptimizers_EXPORT LBFGSBOptimizer:
65 {
66 public:
72 
74  itkNewMacro(Self);
75 
78 
83 
88 
90  typedef vnl_vector< double > InternalBoundValueType;
91 
93  typedef vnl_vector< long > InternalBoundSelectionType;
94 
97 
99  virtual void StartOptimization(void) ITK_OVERRIDE;
100 
102  virtual void SetCostFunction(SingleValuedCostFunction *costFunction) ITK_OVERRIDE;
103 
107  virtual void SetTrace(bool flag);
108 
109  itkGetMacro(Trace, bool);
110  itkBooleanMacro(Trace);
111 
113  virtual void SetLowerBound(const BoundValueType & value);
114  itkGetConstReferenceMacro(LowerBound,BoundValueType);
116 
118  virtual void SetUpperBound(const BoundValueType & value);
119  itkGetConstReferenceMacro(UpperBound,BoundValueType);
121 
128  virtual void SetBoundSelection(const BoundSelectionType & select);
129  itkGetConstReferenceMacro(BoundSelection,BoundSelectionType);
131 
138  virtual void SetCostFunctionConvergenceFactor(double);
139 
140  itkGetMacro(CostFunctionConvergenceFactor, double);
141 
146  virtual void SetProjectedGradientTolerance(double);
147 
148  itkGetMacro(ProjectedGradientTolerance, double);
149 
151  virtual void SetMaximumNumberOfIterations(unsigned int);
152 
153  itkGetMacro(MaximumNumberOfIterations, unsigned int);
154 
156  virtual void SetMaximumNumberOfEvaluations(unsigned int);
157 
158  itkGetMacro(MaximumNumberOfEvaluations, unsigned int);
159 
161  virtual void SetMaximumNumberOfCorrections(unsigned int);
162 
163  itkGetMacro(MaximumNumberOfCorrections, unsigned int);
164 
166  void SetScales(const ScalesType &)
167  {
168  itkExceptionMacro(<< "This optimizer does not support scales.");
169  }
170 
172  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
173 
175  MeasureType GetValue() const;
176 
179  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
180 
182  virtual const std::string GetStopConditionDescription() const ITK_OVERRIDE;
183 
184 protected:
185  LBFGSBOptimizer();
186  virtual ~LBFGSBOptimizer() ITK_OVERRIDE;
187  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
188 
190 
191 private:
192  ITK_DISALLOW_COPY_AND_ASSIGN(LBFGSBOptimizer);
193 
194  // give the helper access to member variables, to update iteration
195  // counts, etc.
196  friend class LBFGSBOptimizerHelper;
197 
198  bool m_Trace;
199  bool m_OptimizerInitialized;
200  double m_CostFunctionConvergenceFactor;
201  double m_ProjectedGradientTolerance;
202  unsigned int m_MaximumNumberOfIterations;
203  unsigned int m_MaximumNumberOfEvaluations;
204  unsigned int m_MaximumNumberOfCorrections;
205  unsigned int m_CurrentIteration;
206  double m_InfinityNormOfProjectedGradient;
207 
208  InternalOptimizerType * m_VnlOptimizer;
209  BoundValueType m_LowerBound;
210  BoundValueType m_UpperBound;
211  BoundSelectionType m_BoundSelection;
212 };
213 } // end namespace itk
214 
215 #endif
vnl_vector< long > InternalBoundSelectionType
This class is a base for the CostFunctions returning a single value.
Light weight base class for most itk classes.
vnl_vector< double > InternalBoundValueType
Array< long > BoundSelectionType
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
Superclass::ScalesType ScalesType
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds. ...
Array< double > BoundValueType
This class is a base for the Optimization methods that optimize a single valued function.
SmartPointer< Self > Pointer
SingleValuedNonLinearVnlOptimizer Superclass
LBFGSBOptimizerHelper InternalOptimizerType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
SmartPointer< const Self > ConstPointer
Wrapper helper around vnl_lbfgsb.