ITK  4.13.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:
74 
75  typedef Superclass::MetricType MetricType;
76  typedef Superclass::ParametersType ParametersType;
77  typedef Superclass::ScalesType ScalesType;
78 
80  itkNewMacro(Self);
81 
83  itkTypeMacro(LBFGSBOptimizerv4, Superclass);
84 
86  UNBOUNDED = 0,
87  LOWERBOUNDED = 1,
88  BOTHBOUNDED = 2,
89  UPPERBOUNDED = 3
90  };
91 
96 
101 
103  void SetInitialPosition(const ParametersType & param);
104 
107  {
108  return m_InitialPosition;
109  }
110 
112  virtual void StartOptimization(bool doOnlyInitialization = false) ITK_OVERRIDE;
113 
115  virtual void SetMetric(MetricType *metric) ITK_OVERRIDE;
116 
118  void SetLowerBound(const BoundValueType & value);
119 
120  itkGetConstReferenceMacro(LowerBound,BoundValueType);
121 
123  void SetUpperBound(const BoundValueType & value);
124 
125  itkGetConstReferenceMacro(UpperBound,BoundValueType);
126 
133  void SetBoundSelection(const BoundSelectionType & select);
134 
135  itkGetConstReferenceMacro(BoundSelection,BoundSelectionType);
136 
143  virtual void SetCostFunctionConvergenceFactor(double);
144 
145  itkGetConstMacro(CostFunctionConvergenceFactor, double);
146 
148  virtual void SetMaximumNumberOfCorrections(unsigned int);
149 
150  itkGetConstMacro(MaximumNumberOfCorrections, unsigned int);
151 
153  virtual void SetScales(const ScalesType &) ITK_OVERRIDE;
154 
157  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
158 
159 protected:
161  virtual ~LBFGSBOptimizerv4() ITK_OVERRIDE;
162  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
163 
165 
168 
170 
171 private:
172  ITK_DISALLOW_COPY_AND_ASSIGN(LBFGSBOptimizerv4);
173 
174  unsigned int m_MaximumNumberOfCorrections;
175 
176  ParametersType m_InitialPosition;
177  BoundValueType m_LowerBound;
178  BoundValueType m_UpperBound;
179  BoundSelectionType m_BoundSelection;
180 };
181 } // end namespace itk
182 #endif
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
Light weight base class for most itk classes.
Superclass::MetricType MetricType
SmartPointer< const Self > ConstPointer
Superclass::ScalesType ScalesType
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds. ...
Wrapper helper around vnl_lbfgsb.
ParametersType & GetInitialPosition(void)
Array< double > BoundValueType
LBFGSOptimizerBasev4< vnl_lbfgsb > Superclass
Superclass::ParametersType ParametersType
SmartPointer< Self > Pointer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Abstract base for vnl lbfgs algorithm optimizers in ITKv4 registration framework. ...