ROOT  6.06/09
Reference Guide
GSLMultiMinimizer.h
Go to the documentation of this file.
1 // @(#)root/mathmore:$Id$
2 // Author: L. Moneta Tue Dec 19 14:09:15 2006
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
7  * *
8  * This library is free software; you can redistribute it and/or *
9  * modify it under the terms of the GNU General Public License *
10  * as published by the Free Software Foundation; either version 2 *
11  * of the License, or (at your option) any later version. *
12  * *
13  * This library is distributed in the hope that it will be useful, *
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
16  * General Public License for more details. *
17  * *
18  * You should have received a copy of the GNU General Public License *
19  * along with this library (see file COPYING); if not, write *
20  * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
21  * 330, Boston, MA 02111-1307 USA, or contact the author. *
22  * *
23  **********************************************************************/
24 
25 // Header file for class GSLMultiMinimizer
26 
27 #ifndef ROOT_Math_GSLMultiMinimizer
28 #define ROOT_Math_GSLMultiMinimizer
29 
30 #include "gsl/gsl_vector.h"
31 #include "gsl/gsl_multimin.h"
32 #include "gsl/gsl_version.h"
34 
35 #include "Math/Error.h"
36 
37 #include "Math/IFunction.h"
38 
39 #include <cassert>
40 
41 namespace ROOT {
42 
43  namespace Math {
44 
45 
46 /**
47  GSLMultiMinimizer class , for minimizing multi-dimensional function
48  using derivatives
49 
50  @ingroup MultiMin
51 
52 */
54 
55 public:
56 
57  /**
58  Default constructor
59  */
61  fMinimizer(0),
62  fType(0),
63  fVec(0)
64  {
65  switch(type)
66  {
68  fType = gsl_multimin_fdfminimizer_conjugate_fr;
69  break;
71  fType = gsl_multimin_fdfminimizer_conjugate_pr;
72  break;
74  fType = gsl_multimin_fdfminimizer_vector_bfgs;
75  break;
77 #if (GSL_MAJOR_VERSION > 1) || ((GSL_MAJOR_VERSION == 1) && (GSL_MINOR_VERSION >= 9))
78  // bfgs2 is available only for v>= 1.9
79  fType = gsl_multimin_fdfminimizer_vector_bfgs2;
80 #else
81  MATH_INFO_MSG("GSLMultiMinimizer","minimizer BFSG2 does not exist with this GSL version , use BFGS");
82  fType = gsl_multimin_fdfminimizer_vector_bfgs;
83 #endif
84  break;
86  fType = gsl_multimin_fdfminimizer_steepest_descent;
87  break;
88  default:
89  fType = gsl_multimin_fdfminimizer_conjugate_fr;
90  break;
91  }
92 
93  }
94 
95  /**
96  Destructor
97  */
99  if (fMinimizer != 0 ) gsl_multimin_fdfminimizer_free(fMinimizer);
100  // can free vector (is copied inside)
101  if (fVec != 0) gsl_vector_free(fVec);
102  }
103 
104 private:
105  // usually copying is non trivial, so we make this unaccessible
106 
107  /**
108  Copy constructor
109  */
111 
112  /**
113  Assignment operator
114  */
116  if (this == &rhs) return *this; // time saving self-test
117  return *this;
118  }
119 
120 public:
121 
122  /**
123  set the function to be minimize the initial minimizer parameters,
124  step size and tolerance in the line search
125  */
126  int Set(const ROOT::Math::IMultiGradFunction & func, const double * x, double stepSize, double tol) {
127  // create function wrapper
128  fFunc.SetFunction(func);
129  // create minimizer object (free previous one if already existing)
130  unsigned int ndim = func.NDim();
131  CreateMinimizer( ndim );
132  // set initial values
133  if (fVec != 0) gsl_vector_free(fVec);
134  fVec = gsl_vector_alloc( ndim );
135  std::copy(x,x+ndim, fVec->data);
136  assert(fMinimizer != 0);
137  return gsl_multimin_fdfminimizer_set(fMinimizer, fFunc.GetFunc(), fVec, stepSize, tol);
138  }
139 
140  /// create the minimizer from the type and size
141  void CreateMinimizer(unsigned int n) {
142  if (fMinimizer) gsl_multimin_fdfminimizer_free(fMinimizer);
143  fMinimizer = gsl_multimin_fdfminimizer_alloc(fType, n);
144  }
145 
146  std::string Name() const {
147  if (fMinimizer == 0) return "undefined";
148  return std::string(gsl_multimin_fdfminimizer_name(fMinimizer) );
149  }
150 
151  int Iterate() {
152  if (fMinimizer == 0) return -1;
153  return gsl_multimin_fdfminimizer_iterate(fMinimizer);
154  }
155 
156  /// x values at the minimum
157  double * X() const {
158  if (fMinimizer == 0) return 0;
159  gsl_vector * x = gsl_multimin_fdfminimizer_x(fMinimizer);
160  return x->data;
161  }
162 
163  /// function value at the minimum
164  double Minimum() const {
165  if (fMinimizer == 0) return 0;
166  return gsl_multimin_fdfminimizer_minimum(fMinimizer);
167  }
168 
169  /// gradient value at the minimum
170  double * Gradient() const {
171  if (fMinimizer == 0) return 0;
172  gsl_vector * g = gsl_multimin_fdfminimizer_gradient(fMinimizer);
173  return g->data;
174  }
175 
176  /// restart minimization from current point
177  int Restart() {
178  if (fMinimizer == 0) return -1;
179  return gsl_multimin_fdfminimizer_restart(fMinimizer);
180  }
181 
182  /// test gradient (ask from minimizer gradient vector)
183  int TestGradient(double absTol) const {
184  if (fMinimizer == 0) return -1;
185  gsl_vector * g = gsl_multimin_fdfminimizer_gradient(fMinimizer);
186  return gsl_multimin_test_gradient( g, absTol);
187  }
188 
189  /// test gradient (require a vector gradient)
190  int TestGradient(const double * g, double absTol) const {
191  if (fVec == 0 ) return -1;
192  unsigned int n = fVec->size;
193  if (n == 0 ) return -1;
194  std::copy(g,g+n, fVec->data);
195  return gsl_multimin_test_gradient( fVec, absTol);
196  }
197 
198 
199 private:
200 
201  gsl_multimin_fdfminimizer * fMinimizer;
203  const gsl_multimin_fdfminimizer_type * fType;
204  // cached vector to avoid re-allocating every time a new one
205  mutable gsl_vector * fVec;
206 
207 };
208 
209  } // end namespace Math
210 
211 } // end namespace ROOT
212 
213 
214 #endif /* ROOT_Math_GSLMultiMinimizer */
GSLMultiMinimizer(ROOT::Math::EGSLMinimizerType type)
Default constructor.
Interface (abstract class) for multi-dimensional functions providing a gradient calculation.
Definition: IFunction.h:322
const double absTol
void SetFunction(const FuncType &f)
Fill gsl function structure from a C++ Function class.
Namespace for new ROOT classes and functions.
Definition: ROOT.py:1
const gsl_multimin_fdfminimizer_type * fType
void CreateMinimizer(unsigned int n)
create the minimizer from the type and size
#define assert(cond)
Definition: unittest.h:542
EGSLMinimizerType
enumeration specifying the types of GSL minimizers
Definition: GSLMinimizer.h:64
GSLMultiMinDerivFunctionWrapper fFunc
gsl_multimin_fdfminimizer * fMinimizer
Double_t x[n]
Definition: legend1.C:17
GSLMultiMinimizer(const GSLMultiMinimizer &)
Copy constructor.
int TestGradient(const double *g, double absTol) const
test gradient (require a vector gradient)
#define MATH_INFO_MSG(loc, str)
Definition: Error.h:44
const double tol
virtual unsigned int NDim() const =0
Retrieve the dimension of the function.
int Set(const ROOT::Math::IMultiGradFunction &func, const double *x, double stepSize, double tol)
set the function to be minimize the initial minimizer parameters, step size and tolerance in the line...
double * Gradient() const
gradient value at the minimum
int type
Definition: TGX11.cxx:120
int Restart()
restart minimization from current point
double * X() const
x values at the minimum
double func(double *x, double *p)
Definition: stressTF1.cxx:213
GSLMultiMinimizer class , for minimizing multi-dimensional function using derivatives.
Namespace for new Math classes and functions.
Wrapper for a multi-dimensional function with derivatives used in GSL multidim minimization algorithm...
double Minimum() const
function value at the minimum
GSLMultiMinimizer & operator=(const GSLMultiMinimizer &rhs)
Assignment operator.
const Int_t n
Definition: legend1.C:16
int TestGradient(double absTol) const
test gradient (ask from minimizer gradient vector)