ITK  5.2.0
Insight Toolkit
itkLBFGSBOptimizer.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 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 : public SingleValuedNonLinearVnlOptimizer
64 {
65 public:
66  ITK_DISALLOW_COPY_AND_MOVE(LBFGSBOptimizer);
67 
73 
75  itkNewMacro(Self);
76 
79 
84 
89 
91  using InternalBoundValueType = vnl_vector<double>;
92 
94  using InternalBoundSelectionType = vnl_vector<long>;
95 
98 
100  void
101  StartOptimization() override;
102 
104  void
105  SetCostFunction(SingleValuedCostFunction * costFunction) override;
106 
110  virtual void
111  SetTrace(bool flag);
112 
113  itkGetMacro(Trace, bool);
114  itkBooleanMacro(Trace);
115 
117  virtual void
118  SetLowerBound(const BoundValueType & value);
119  itkGetConstReferenceMacro(LowerBound, BoundValueType);
121 
123  virtual void
124  SetUpperBound(const BoundValueType & value);
125  itkGetConstReferenceMacro(UpperBound, BoundValueType);
127 
134  virtual void
135  SetBoundSelection(const BoundSelectionType & value);
136  itkGetConstReferenceMacro(BoundSelection, BoundSelectionType);
138 
145  virtual void
146  SetCostFunctionConvergenceFactor(double);
147 
148  itkGetMacro(CostFunctionConvergenceFactor, double);
149 
154  virtual void
155  SetProjectedGradientTolerance(double);
156 
157  itkGetMacro(ProjectedGradientTolerance, double);
158 
160  virtual void
161  SetMaximumNumberOfIterations(unsigned int);
162 
163  itkGetMacro(MaximumNumberOfIterations, unsigned int);
164 
166  virtual void
167  SetMaximumNumberOfEvaluations(unsigned int);
168 
169  itkGetMacro(MaximumNumberOfEvaluations, unsigned int);
170 
172  virtual void
173  SetMaximumNumberOfCorrections(unsigned int);
174 
175  itkGetMacro(MaximumNumberOfCorrections, unsigned int);
176 
178  void
180  {
181  itkExceptionMacro(<< "This optimizer does not support scales.");
182  }
183 
185  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
186 
188  MeasureType
189  GetValue() const;
190 
193  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
194 
196  const std::string
197  GetStopConditionDescription() const override;
198 
199 protected:
200  LBFGSBOptimizer();
201  ~LBFGSBOptimizer() override;
202  void
203  PrintSelf(std::ostream & os, Indent indent) const override;
204 
205  using CostFunctionAdaptorType = Superclass::CostFunctionAdaptorType;
206 
207 private:
208  // give the helper access to member variables, to update iteration
209  // counts, etc.
210  friend class LBFGSBOptimizerHelper;
211 
212  bool m_Trace{ false };
213  bool m_OptimizerInitialized{ false };
214  double m_CostFunctionConvergenceFactor{ 1e+7 };
215  double m_ProjectedGradientTolerance{ 1e-5 };
216  unsigned int m_MaximumNumberOfIterations{ 500 };
217  unsigned int m_MaximumNumberOfEvaluations{ 500 };
218  unsigned int m_MaximumNumberOfCorrections{ 5 };
219  unsigned int m_CurrentIteration{ 0 };
220  double m_InfinityNormOfProjectedGradient{ 0.0 };
221 
222  InternalOptimizerType * m_VnlOptimizer{ nullptr };
226 };
227 } // end namespace itk
228 
229 #endif
itk::LBFGSBOptimizer::InternalBoundValueType
vnl_vector< double > InternalBoundValueType
Definition: itkLBFGSBOptimizer.h:91
itk::LBFGSBOptimizer::m_LowerBound
BoundValueType m_LowerBound
Definition: itkLBFGSBOptimizer.h:223
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::LBFGSBOptimizer
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
Definition: itkLBFGSBOptimizer.h:63
itk::LBFGSBOptimizer::SetScales
void SetScales(const ScalesType &)
Definition: itkLBFGSBOptimizer.h:179
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
LBFGSBOptimizerHelper
Wrapper helper around vnl_lbfgsb.
itk::LBFGSBOptimizer::InternalBoundSelectionType
vnl_vector< long > InternalBoundSelectionType
Definition: itkLBFGSBOptimizer.h:94
itk::LBFGSBOptimizerHelper
class ITK_FORWARD_EXPORT LBFGSBOptimizerHelper
Definition: itkLBFGSBOptimizer.h:35
itk::SingleValuedCostFunction
This class is a base for the CostFunctions returning a single value.
Definition: itkSingleValuedCostFunction.h:34
itk::LBFGSBOptimizer::m_BoundSelection
BoundSelectionType m_BoundSelection
Definition: itkLBFGSBOptimizer.h:225
itkIntTypes.h
itkSingleValuedNonLinearVnlOptimizer.h
itk::LBFGSBOptimizer::m_UpperBound
BoundValueType m_UpperBound
Definition: itkLBFGSBOptimizer.h:224
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::LBFGSBOptimizer::CostFunctionAdaptorType
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
Definition: itkLBFGSBOptimizer.h:205
itk::NonLinearOptimizer::ScalesType
Superclass::ScalesType ScalesType
Definition: itkNonLinearOptimizer.h:55
itk::Array< double >
itk::Math::e
static constexpr double e
Definition: itkMath.h:54
itk::SingleValuedNonLinearVnlOptimizer
This class is a base for the Optimization methods that optimize a single valued function.
Definition: itkSingleValuedNonLinearVnlOptimizer.h:37