1 /*=========================================================================
2 Program: vv http://www.creatis.insa-lyon.fr/rio/vv
5 - University of LYON http://www.universite-lyon.fr/
6 - Léon Bérard cancer center http://www.centreleonberard.fr
7 - CREATIS CNRS laboratory http://www.creatis.insa-lyon.fr
9 This software is distributed WITHOUT ANY WARRANTY; without even
10 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11 PURPOSE. See the copyright notices for more information.
13 It is distributed under dual licence
15 - BSD See included LICENSE.txt file
16 - CeCILL-B http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ===========================================================================**/
18 #ifndef __clitkLBFGSBOptimizer_h
19 #define __clitkLBFGSBOptimizer_h
20 #include "itkSingleValuedNonLinearVnlOptimizer.h"
25 /** \class LBFGSBOptimizerHelper
26 * \brief Wrapper helper around vnl_lbfgsb.
28 * This class is used to translate iteration events, etc, from
29 * vnl_lbfgsb into iteration events in ITK.
31 class ITK_EXPORT LBFGSBOptimizerHelper;
34 /** \class LBFGSBOptimizer
35 * \brief Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
37 * This class is a wrapper for converted fortan code for performing limited
38 * memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
39 * The algorithm miminizes a nonlinear function f(x) of n variables subject to
40 * simple bound constraints of l <= x <= u.
42 * See also the documentation in Numerics/lbfgsb.c
46 * [1] R. H. Byrd, P. Lu and J. Nocedal.
47 * A Limited Memory Algorithm for Bound Constrained Optimization, (1995),
48 * SIAM Journal on Scientific and Statistical Computing ,
49 * 16, 5, pp. 1190-1208.
51 * [2] C. Zhu, R. H. Byrd and J. Nocedal.
52 * L-BFGS-B: Algorithm 778: L-BFGS-B, FORTRAN routines for large scale
53 * bound constrained optimization (1997),
54 * ACM Transactions on Mathematical Software,
55 * Vol 23, Num. 4, pp. 550 - 560.
57 * \ingroup Numerics Optimizers
59 class ITK_EXPORT LBFGSBOptimizer :
60 public itk::SingleValuedNonLinearVnlOptimizer
63 /** Standard "Self" typedef. */
64 typedef LBFGSBOptimizer Self;
65 typedef itk::SingleValuedNonLinearVnlOptimizer Superclass;
66 typedef itk::SmartPointer<Self> Pointer;
67 typedef itk::SmartPointer<const Self> ConstPointer;
69 /** Method for creation through the object factory. */
72 /** Run-time type information (and related methods). */
73 itkTypeMacro( LBFGSBOptimizer, SingleValuedNonLinearVnlOptimizer );
76 * Use for defining the lower and upper bounds on the variables.
78 typedef itk::Array<double> BoundValueType;
80 /** BoundSelection type
81 * Use for defining the boundary condition for each variables.
83 typedef itk::Array<long> BoundSelectionType;
85 /** Internal boundary value storage type */
86 typedef vnl_vector<double> InternalBoundValueType;
88 /** Internal boundary selection storage type */
89 typedef vnl_vector<long> InternalBoundSelectionType;
91 /** The vnl optimizer */
92 typedef LBFGSBOptimizerHelper InternalOptimizerType;
95 /** Start optimization with an initial value. */
96 void StartOptimization( void );
98 /** Plug in a Cost Function into the optimizer */
99 virtual void SetCostFunction( itk::SingleValuedCostFunction * costFunction );
101 /** Set the lower bound value for each variable. */
102 virtual void SetLowerBound( const BoundValueType & value );
103 virtual const BoundValueType & GetLowerBound();
105 /** Set the upper bound value for each variable. */
106 virtual void SetUpperBound( const BoundValueType & value );
107 virtual const BoundValueType & GetUpperBound();
109 /** Set the boundary condition for each variable, where
110 * select[i] = 0 if x[i] is unbounded,
111 * = 1 if x[i] has only a lower bound,
112 * = 2 if x[i] has both lower and upper bounds, and
113 * = 3 if x[1] has only an upper bound
115 virtual void SetBoundSelection( const BoundSelectionType & select );
116 virtual const BoundSelectionType & GetBoundSelection();
118 /** Set/Get the CostFunctionConvergenceFactor. Algorithm terminates
119 * when the reduction in cost function is less than factor * epsmcj
120 * where epsmch is the machine precision.
121 * Typical values for factor: 1e+12 for low accuracy;
122 * 1e+7 for moderate accuracy and 1e+1 for extremely high accuracy.
124 virtual void SetCostFunctionConvergenceFactor( double );
125 itkGetMacro( CostFunctionConvergenceFactor, double );
127 /** Set/Get the ProjectedGradientTolerance. Algorithm terminates
128 * when the project gradient is below the tolerance. Default value
131 virtual void SetProjectedGradientTolerance( double );
132 itkGetMacro( ProjectedGradientTolerance, double );
134 /** Set/Get the MaximumNumberOfIterations. Default is 500 */
135 virtual void SetMaximumNumberOfIterations( unsigned int );
136 itkGetMacro( MaximumNumberOfIterations, unsigned int );
138 /** Set/Get the MaximumNumberOfEvaluations. Default is 500 */
139 virtual void SetMaximumNumberOfEvaluations( unsigned int );
140 itkGetMacro( MaximumNumberOfEvaluations, unsigned int );
142 /** Set/Get the MaximumNumberOfCorrections. Default is 5 */
143 virtual void SetMaximumNumberOfCorrections( unsigned int );
144 itkGetMacro( MaximumNumberOfCorrections, unsigned int );
146 /** This optimizer does not support scaling of the derivatives. */
147 void SetScales( const ScalesType & )
149 itkExceptionMacro( << "This optimizer does not support scales." );
152 /** Get the current iteration number. */
153 itkGetConstReferenceMacro( CurrentIteration, unsigned int );
155 /** Get the current cost function value. */
156 itkGetConstReferenceMacro( Value, MeasureType );
158 /** Get the current infinity norm of the project gradient of the cost
160 itkGetConstReferenceMacro( InfinityNormOfProjectedGradient, double );
162 /** Get the reason for termination */
163 const std::string GetStopConditionDescription() const;
167 virtual ~LBFGSBOptimizer();
168 void PrintSelf(std::ostream& os, itk::Indent indent) const;
170 typedef Superclass::CostFunctionAdaptorType CostFunctionAdaptorType;
173 LBFGSBOptimizer(const Self&); //purposely not implemented
174 void operator=(const Self&); //purposely not implemented
176 // give the helper access to member variables, to update iteration
178 friend class LBFGSBOptimizerHelper;
180 bool m_OptimizerInitialized;
181 InternalOptimizerType * m_VnlOptimizer;
182 #if ITK_VERSION_MAJOR > 3
183 mutable std::ostringstream m_StopConditionDescription;
185 mutable itk::OStringStream m_StopConditionDescription;
187 BoundValueType m_LowerBound;
188 BoundValueType m_UpperBound;
189 BoundSelectionType m_BoundSelection;
191 double m_CostFunctionConvergenceFactor;
192 double m_ProjectedGradientTolerance;
193 unsigned int m_MaximumNumberOfIterations;
194 unsigned int m_MaximumNumberOfEvaluations;
195 unsigned int m_MaximumNumberOfCorrections;
197 unsigned int m_CurrentIteration;
199 double m_InfinityNormOfProjectedGradient;
203 } // end namespace clitk