ITK  4.8.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 class ITKOptimizersv4_EXPORT LBFGSBOptimizerHelperv4;
37 
64 class ITKOptimizersv4_EXPORT LBFGSBOptimizerv4:
65  public LBFGSOptimizerBasev4< vnl_lbfgsb >
66 {
67 public:
73 
74  typedef Superclass::MetricType MetricType;
75  typedef Superclass::ParametersType ParametersType;
76  typedef Superclass::ScalesType ScalesType;
77 
79  itkNewMacro(Self);
80 
82  itkTypeMacro(LBFGSBOptimizerv4, Superclass);
83 
85  UNBOUNDED = 0,
86  LOWERBOUNDED = 1,
87  BOTHBOUNDED = 2,
88  UPPERBOUNDED = 3
89  };
90 
95 
100 
102  void SetInitialPosition(const ParametersType & param);
103 
106  {
107  return m_InitialPosition;
108  }
109 
111  virtual void StartOptimization(bool doOnlyInitialization = false) ITK_OVERRIDE;
112 
114  virtual void SetMetric(MetricType *metric) ITK_OVERRIDE;
115 
117  void SetLowerBound(const BoundValueType & value);
118 
119  itkGetConstReferenceMacro(LowerBound,BoundValueType);
120 
122  void SetUpperBound(const BoundValueType & value);
123 
124  itkGetConstReferenceMacro(UpperBound,BoundValueType);
125 
132  void SetBoundSelection(const BoundSelectionType & select);
133 
134  itkGetConstReferenceMacro(BoundSelection,BoundSelectionType);
135 
142  virtual void SetCostFunctionConvergenceFactor(double);
143 
144  itkGetConstMacro(CostFunctionConvergenceFactor, double);
145 
147  virtual void SetMaximumNumberOfCorrections(unsigned int);
148 
149  itkGetConstMacro(MaximumNumberOfCorrections, unsigned int);
150 
152  virtual void SetScales(const ScalesType &) ITK_OVERRIDE;
153 
156  itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
157 
158 protected:
160  virtual ~LBFGSBOptimizerv4();
161  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
162 
164 
167 
169 
170 private:
171  LBFGSBOptimizerv4(const Self &); //purposely not implemented
172  void operator=(const Self &); //purposely not implemented
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. ...