ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkLBFGSBOptimizerv4.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 itkLBFGSBOptimizerv4_h
19 #define itkLBFGSBOptimizerv4_h
20 
22 #include "vnl/algo/vnl_lbfgsb.h"
23 #include "ITKOptimizersv4Export.h"
24 
25 namespace itk
26 {
27 /* Necessary forward declaration */
36 // Forward reference because of private implementation
37 class ITK_FORWARD_EXPORT LBFGSBOptimizerHelperv4;
38 
65 class ITKOptimizersv4_EXPORT LBFGSBOptimizerv4:
66  public LBFGSOptimizerBasev4< vnl_lbfgsb >
67 {
68 public:
69  ITK_DISALLOW_COPY_AND_ASSIGN(LBFGSBOptimizerv4);
70 
76 
77  using MetricType = Superclass::MetricType;
78  using ParametersType = Superclass::ParametersType;
79  using ScalesType = Superclass::ScalesType;
80 
82  itkNewMacro(Self);
83 
85  itkTypeMacro(LBFGSBOptimizerv4, Superclass);
86 
88  UNBOUNDED = 0,
89  LOWERBOUNDED = 1,
90  BOTHBOUNDED = 2,
91  UPPERBOUNDED = 3
92  };
93 
98 
103 
105  void SetInitialPosition(const ParametersType & param);
106 
109  {
110  return m_InitialPosition;
111  }
112 
114  void StartOptimization(bool doOnlyInitialization = false) override;
115 
117  void SetMetric(MetricType *metric) override;
118 
120  void SetLowerBound(const BoundValueType & value);
121 
122  itkGetConstReferenceMacro(LowerBound,BoundValueType);
123 
125  void SetUpperBound(const BoundValueType & value);
126 
127  itkGetConstReferenceMacro(UpperBound,BoundValueType);
128 
135  void SetBoundSelection(const BoundSelectionType & select);
136 
137  itkGetConstReferenceMacro(BoundSelection,BoundSelectionType);
138 
145  virtual void SetCostFunctionConvergenceFactor(double);
146 
147  itkGetConstMacro(CostFunctionConvergenceFactor, double);
148 
150  virtual void SetMaximumNumberOfCorrections(unsigned int);
151 
152  itkGetConstMacro(MaximumNumberOfCorrections, unsigned int);
153 
155  void SetScales(const ScalesType &) override;
156 
159  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
160 
161 protected:
163  ~LBFGSBOptimizerv4() override;
164  void PrintSelf(std::ostream & os, Indent indent) const override;
165 
166  using CostFunctionAdaptorType = Superclass::CostFunctionAdaptorType;
167 
170 
172 
173 private:
174  unsigned int m_MaximumNumberOfCorrections{5};
175 
180 };
181 } // end namespace itk
182 #endif
Superclass::ParametersType ParametersType
Light weight base class for most itk classes.
class ITK_FORWARD_EXPORT LBFGSBOptimizerHelperv4
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds. ...
Superclass::MetricType MetricType
Wrapper helper around vnl_lbfgsb.
Superclass::ScalesType ScalesType
BoundSelectionType m_BoundSelection
Control indentation during Print() invocation.
Definition: itkIndent.h:49
ParametersType & GetInitialPosition()
Abstract base for vnl lbfgs algorithm optimizers in ITKv4 registration framework. ...