ITK  5.4.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  * https://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 #include <memory> // For unique_ptr.
25 
26 namespace itk
27 {
28 /* Necessary forward declaration see below for definition */
36 class ITK_FORWARD_EXPORT LBFGSBOptimizerHelper;
37 
64 class ITKOptimizers_EXPORT LBFGSBOptimizer : public SingleValuedNonLinearVnlOptimizer
65 {
66 public:
67  ITK_DISALLOW_COPY_AND_MOVE(LBFGSBOptimizer);
68 
74 
76  itkNewMacro(Self);
77 
79  itkOverrideGetNameOfClassMacro(LBFGSBOptimizer);
80 
85 
90 
92  using InternalBoundValueType = vnl_vector<double>;
93 
95  using InternalBoundSelectionType = vnl_vector<long>;
96 
99 
101  void
102  StartOptimization() override;
103 
105  void
106  SetCostFunction(SingleValuedCostFunction * costFunction) override;
107 
111  virtual void
112  SetTrace(bool flag);
113 
114  itkGetMacro(Trace, bool);
115  itkBooleanMacro(Trace);
116 
118  virtual void
119  SetLowerBound(const BoundValueType & value);
120  itkGetConstReferenceMacro(LowerBound, BoundValueType);
124  virtual void
125  SetUpperBound(const BoundValueType & value);
126  itkGetConstReferenceMacro(UpperBound, BoundValueType);
135  virtual void
136  SetBoundSelection(const BoundSelectionType & value);
137  itkGetConstReferenceMacro(BoundSelection, BoundSelectionType);
146  virtual void
147  SetCostFunctionConvergenceFactor(double);
148 
149  itkGetMacro(CostFunctionConvergenceFactor, double);
150 
155  virtual void
156  SetProjectedGradientTolerance(double);
157 
158  itkGetMacro(ProjectedGradientTolerance, double);
159 
161  virtual void
162  SetMaximumNumberOfIterations(unsigned int);
163 
164  itkGetMacro(MaximumNumberOfIterations, unsigned int);
165 
167  virtual void
168  SetMaximumNumberOfEvaluations(unsigned int);
169 
170  itkGetMacro(MaximumNumberOfEvaluations, unsigned int);
171 
173  virtual void
174  SetMaximumNumberOfCorrections(unsigned int);
175 
176  itkGetMacro(MaximumNumberOfCorrections, unsigned int);
177 
179  void
181  {
182  itkExceptionMacro("This optimizer does not support scales.");
183  }
184 
186  itkGetConstReferenceMacro(CurrentIteration, unsigned int);
187 
189  MeasureType
190  GetValue() const;
191 
194  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
195 
197  const std::string
198  GetStopConditionDescription() const override;
199 
201  bool
202  CanUseScales() const override
203  {
204  return false;
205  }
206 
207 protected:
208  LBFGSBOptimizer();
209  ~LBFGSBOptimizer() override;
210  void
211  PrintSelf(std::ostream & os, Indent indent) const override;
212 
213  using CostFunctionAdaptorType = Superclass::CostFunctionAdaptorType;
214 
215 private:
216  // give the helper access to member variables, to update iteration
217  // counts, etc.
218  friend class LBFGSBOptimizerHelper;
219 
220  bool m_Trace{ false };
221  bool m_OptimizerInitialized{ false };
222  double m_CostFunctionConvergenceFactor{ 1e+7 };
223  double m_ProjectedGradientTolerance{ 1e-5 };
224  unsigned int m_MaximumNumberOfIterations{ 500 };
225  unsigned int m_MaximumNumberOfEvaluations{ 500 };
226  unsigned int m_MaximumNumberOfCorrections{ 5 };
227  unsigned int m_CurrentIteration{ 0 };
228  double m_InfinityNormOfProjectedGradient{ 0.0 };
229 
230  std::unique_ptr<InternalOptimizerType> m_VnlOptimizer;
231  BoundValueType m_LowerBound{};
232  BoundValueType m_UpperBound{};
233  BoundSelectionType m_BoundSelection{};
234 };
235 } // end namespace itk
236 
237 #endif
itk::LBFGSBOptimizer::CanUseScales
bool CanUseScales() const override
Definition: itkLBFGSBOptimizer.h:202
itk::LBFGSBOptimizer::InternalBoundValueType
vnl_vector< double > InternalBoundValueType
Definition: itkLBFGSBOptimizer.h:92
itk::LBFGSBOptimizer::m_VnlOptimizer
std::unique_ptr< InternalOptimizerType > m_VnlOptimizer
Definition: itkLBFGSBOptimizer.h:230
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:64
itk::LBFGSBOptimizer::SetScales
void SetScales(const ScalesType &)
Definition: itkLBFGSBOptimizer.h:180
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
LBFGSBOptimizerHelper
Wrapper helper around vnl_lbfgsb.
itk::LBFGSBOptimizer::InternalBoundSelectionType
vnl_vector< long > InternalBoundSelectionType
Definition: itkLBFGSBOptimizer.h:95
itk::LBFGSBOptimizerHelper
class ITK_FORWARD_EXPORT LBFGSBOptimizerHelper
Definition: itkLBFGSBOptimizer.h:36
itk::SingleValuedCostFunction
This class is a base for the CostFunctions returning a single value.
Definition: itkSingleValuedCostFunction.h:34
itkIntTypes.h
itkSingleValuedNonLinearVnlOptimizer.h
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:213
itk::NonLinearOptimizer::ScalesType
Superclass::ScalesType ScalesType
Definition: itkNonLinearOptimizer.h:55
itk::Array< double >
itk::Math::e
static constexpr double e
Definition: itkMath.h:56
itk::SingleValuedNonLinearVnlOptimizer
This class is a base for the Optimization methods that optimize a single valued function.
Definition: itkSingleValuedNonLinearVnlOptimizer.h:37