]> Creatis software - clitk.git/blob - registration/clitkLBFGSBOptimizer.h
Debug RTStruct conversion with empty struc
[clitk.git] / registration / clitkLBFGSBOptimizer.h
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to: 
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
8
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.
12
13   It is distributed under dual licence
14
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"
21
22 #include "clitkCommon.h"
23
24 namespace clitk
25 {
26
27 /** \class LBFGSBOptimizerHelper
28  * \brief Wrapper helper around vnl_lbfgsb.
29  *
30  * This class is used to translate iteration events, etc, from
31  * vnl_lbfgsb into iteration events in ITK.
32  */
33 class ITK_EXPORT LBFGSBOptimizerHelper;
34
35
36 /** \class LBFGSBOptimizer
37  * \brief Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
38  *
39  * This class is a wrapper for converted fortan code for performing limited
40  * memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
41  * The algorithm miminizes a nonlinear function f(x) of n variables subject to
42  * simple bound constraints of l <= x <= u.
43  *
44  * See also the documentation in Numerics/lbfgsb.c
45  *
46  * References:
47  *
48  * [1] R. H. Byrd, P. Lu and J. Nocedal. 
49  * A Limited Memory Algorithm for Bound Constrained Optimization, (1995), 
50  * SIAM Journal on Scientific and Statistical Computing , 
51  * 16, 5, pp. 1190-1208. 
52  *
53  * [2] C. Zhu, R. H. Byrd and J. Nocedal. 
54  * L-BFGS-B: Algorithm 778: L-BFGS-B, FORTRAN routines for large scale 
55  * bound constrained optimization (1997), 
56  * ACM Transactions on Mathematical Software, 
57  * Vol 23, Num. 4, pp. 550 - 560. 
58  *
59  * \ingroup Numerics Optimizers
60  */
61 class ITK_EXPORT LBFGSBOptimizer : 
62     public itk::SingleValuedNonLinearVnlOptimizer
63 {
64 public:
65   /** Standard "Self" typedef. */
66   typedef LBFGSBOptimizer                     Self;
67   typedef itk::SingleValuedNonLinearVnlOptimizer   Superclass;
68   typedef itk::SmartPointer<Self>                  Pointer;
69   typedef itk::SmartPointer<const Self>            ConstPointer;
70   
71   /** Method for creation through the object factory. */
72   itkNewMacro(Self);
73
74   /** Run-time type information (and related methods). */
75   itkTypeMacro( LBFGSBOptimizer, SingleValuedNonLinearVnlOptimizer );
76
77   /**  BoundValue type.
78    *  Use for defining the lower and upper bounds on the variables. 
79    */
80   typedef itk::Array<double>             BoundValueType;
81
82   /** BoundSelection type
83    * Use for defining the boundary condition for each variables. 
84    */
85   typedef itk::Array<long>                BoundSelectionType;
86
87   /** Internal boundary value storage type */
88   typedef vnl_vector<double>    InternalBoundValueType;
89
90   /** Internal boundary selection storage type */
91   typedef vnl_vector<long>      InternalBoundSelectionType;
92
93   /** The vnl optimizer */
94   typedef LBFGSBOptimizerHelper InternalOptimizerType;
95
96
97   /** Start optimization with an initial value. */
98   virtual void StartOptimization( void ) ITK_OVERRIDE;
99
100   /** Plug in a Cost Function into the optimizer  */
101   virtual void SetCostFunction( itk::SingleValuedCostFunction * costFunction ) ITK_OVERRIDE;
102
103   /** Set the lower bound value for each variable. */
104   virtual void SetLowerBound( const BoundValueType & value );
105   virtual const BoundValueType & GetLowerBound();
106
107   /** Set the upper bound value for each variable. */
108   virtual void SetUpperBound( const BoundValueType & value );
109   virtual const BoundValueType & GetUpperBound();
110
111   /** Set the boundary condition for each variable, where
112    * select[i] = 0 if x[i] is unbounded,
113    *           = 1 if x[i] has only a lower bound,
114    *           = 2 if x[i] has both lower and upper bounds, and
115    *           = 3 if x[1] has only an upper bound
116    */
117   virtual void SetBoundSelection( const BoundSelectionType & select );
118   virtual const BoundSelectionType & GetBoundSelection();
119
120   /** Set/Get the CostFunctionConvergenceFactor. Algorithm terminates
121    * when the reduction in cost function is less than factor * epsmcj
122    * where epsmch is the machine precision.
123    * Typical values for factor: 1e+12 for low accuracy; 
124    * 1e+7 for moderate accuracy and 1e+1 for extremely high accuracy.
125    */
126   virtual void SetCostFunctionConvergenceFactor( double );
127   itkGetMacro( CostFunctionConvergenceFactor, double );
128
129   /** Set/Get the ProjectedGradientTolerance. Algorithm terminates
130    * when the project gradient is below the tolerance. Default value
131    * is 1e-5.
132    */
133   virtual void SetProjectedGradientTolerance( double );
134   itkGetMacro( ProjectedGradientTolerance, double );
135
136   /** Set/Get the MaximumNumberOfIterations. Default is 500 */
137   virtual void SetMaximumNumberOfIterations( unsigned int );
138   itkGetMacro( MaximumNumberOfIterations, unsigned int );
139
140   /** Set/Get the MaximumNumberOfEvaluations. Default is 500 */
141   virtual void SetMaximumNumberOfEvaluations( unsigned int );
142   itkGetMacro( MaximumNumberOfEvaluations, unsigned int );
143
144   /** Set/Get the MaximumNumberOfCorrections. Default is 5 */
145   virtual void SetMaximumNumberOfCorrections( unsigned int );
146   itkGetMacro( MaximumNumberOfCorrections, unsigned int );
147
148   /** This optimizer does not support scaling of the derivatives. */
149   void SetScales( const ScalesType & )
150     {
151     itkExceptionMacro( << "This optimizer does not support scales." );
152     }
153
154   /** Get the current iteration number. */
155   itkGetConstReferenceMacro( CurrentIteration, unsigned int );
156
157   /** Get the current cost function value. */
158   itkGetConstReferenceMacro( Value, MeasureType );
159
160   /** Get the current infinity norm of the project gradient of the cost
161    * function. */
162   itkGetConstReferenceMacro( InfinityNormOfProjectedGradient, double );
163
164   /** Get the reason for termination */
165   virtual const std::string GetStopConditionDescription() const ITK_OVERRIDE;
166
167 protected:
168   LBFGSBOptimizer();
169   virtual ~LBFGSBOptimizer();
170   void PrintSelf(std::ostream& os, itk::Indent indent) const ITK_OVERRIDE;
171
172   typedef Superclass::CostFunctionAdaptorType   CostFunctionAdaptorType;
173
174 private:
175   LBFGSBOptimizer(const Self&); //purposely not implemented
176   void operator=(const Self&); //purposely not implemented
177
178   // give the helper access to member variables, to update iteration
179   // counts, etc.
180   friend class LBFGSBOptimizerHelper;
181
182   bool                     m_OptimizerInitialized;
183   InternalOptimizerType  * m_VnlOptimizer;
184   mutable std::ostringstream    m_StopConditionDescription;
185   BoundValueType           m_LowerBound;
186   BoundValueType           m_UpperBound;
187   BoundSelectionType       m_BoundSelection;
188
189   double                   m_CostFunctionConvergenceFactor;
190   double                   m_ProjectedGradientTolerance;
191   unsigned int             m_MaximumNumberOfIterations;
192   unsigned int             m_MaximumNumberOfEvaluations;
193   unsigned int             m_MaximumNumberOfCorrections;
194
195   unsigned int             m_CurrentIteration;
196   MeasureType              m_Value;
197   double                   m_InfinityNormOfProjectedGradient;
198
199 };
200
201 } // end namespace clitk
202
203 #endif