ITK  5.0.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:
67  ITK_DISALLOW_COPY_AND_ASSIGN(LBFGSBOptimizer);
68 
74 
76  itkNewMacro(Self);
77 
80 
85 
90 
92  using InternalBoundValueType = vnl_vector< double >;
93 
95  using InternalBoundSelectionType = vnl_vector< long >;
96 
99 
101  void StartOptimization() override;
102 
104  void SetCostFunction(SingleValuedCostFunction *costFunction) override;
105 
109  virtual void SetTrace(bool flag);
110 
111  itkGetMacro(Trace, bool);
112  itkBooleanMacro(Trace);
113 
115  virtual void SetLowerBound(const BoundValueType & value);
116  itkGetConstReferenceMacro(LowerBound,BoundValueType);
118 
120  virtual void SetUpperBound(const BoundValueType & value);
121  itkGetConstReferenceMacro(UpperBound,BoundValueType);
123 
130  virtual void SetBoundSelection(const BoundSelectionType & select);
131  itkGetConstReferenceMacro(BoundSelection,BoundSelectionType);
133 
140  virtual void SetCostFunctionConvergenceFactor(double);
141 
142  itkGetMacro(CostFunctionConvergenceFactor, double);
143 
148  virtual void SetProjectedGradientTolerance(double);
149 
150  itkGetMacro(ProjectedGradientTolerance, double);
151 
153  virtual void SetMaximumNumberOfIterations(unsigned int);
154 
155  itkGetMacro(MaximumNumberOfIterations, unsigned int);
156 
158  virtual void SetMaximumNumberOfEvaluations(unsigned int);
159 
160  itkGetMacro(MaximumNumberOfEvaluations, unsigned int);
161 
163  virtual void SetMaximumNumberOfCorrections(unsigned int);
164 
165  itkGetMacro(MaximumNumberOfCorrections, unsigned int);
166 
168  void SetScales(const ScalesType &)
169  {
170  itkExceptionMacro(<< "This optimizer does not support scales.");
171  }
172 
174  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
175 
177  MeasureType GetValue() const;
178 
181  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
182 
184  const std::string GetStopConditionDescription() const override;
185 
186 protected:
187  LBFGSBOptimizer();
188  ~LBFGSBOptimizer() override;
189  void PrintSelf(std::ostream & os, Indent indent) const override;
190 
191  using CostFunctionAdaptorType = Superclass::CostFunctionAdaptorType;
192 
193 private:
194  // give the helper access to member variables, to update iteration
195  // counts, etc.
196  friend class LBFGSBOptimizerHelper;
197 
198  bool m_Trace{false};
199  bool m_OptimizerInitialized{false};
200  double m_CostFunctionConvergenceFactor{1e+7};
201  double m_ProjectedGradientTolerance{1e-5};
202  unsigned int m_MaximumNumberOfIterations{500};
203  unsigned int m_MaximumNumberOfEvaluations{500};
204  unsigned int m_MaximumNumberOfCorrections{5};
205  unsigned int m_CurrentIteration{0};
206  double m_InfinityNormOfProjectedGradient{0.0};
207 
208  InternalOptimizerType * m_VnlOptimizer{nullptr};
212 };
213 } // end namespace itk
214 
215 #endif
This class is a base for the CostFunctions returning a single value.
Light weight base class for most itk classes.
BoundSelectionType m_BoundSelection
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
BoundValueType m_UpperBound
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds. ...
vnl_vector< double > InternalBoundValueType
This class is a base for the Optimization methods that optimize a single valued function.
BoundValueType m_LowerBound
Superclass::ScalesType ScalesType
static constexpr double e
The base of the natural logarithm or Euler&#39;s number
Definition: itkMath.h:53
class ITK_FORWARD_EXPORT LBFGSBOptimizerHelper
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void SetScales(const ScalesType &)
Wrapper helper around vnl_lbfgsb.
vnl_vector< long > InternalBoundSelectionType