BLMinimize.hh

Go to the documentation of this file.
00001 //      BLMinimize.hh
00002 #define USE_GSL_MULTIMIN
00003 //#define USE_TMINUIT
00004 /*
00005 This source file is part of G4beamline, http://g4beamline.muonsinc.com
00006 Copyright (C) 2003,2004,2005,2006 by Tom Roberts, all rights reserved.
00007 
00008 This program is free software; you can redistribute it and/or
00009 modify it under the terms of the GNU General Public License
00010 as published by the Free Software Foundation; either version 2
00011 of the License, or (at your option) any later version.
00012 
00013 This program is distributed in the hope that it will be useful,
00014 but WITHOUT ANY WARRANTY; without even the implied warranty of
00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016 GNU General Public License for more details.
00017 
00018 http://www.gnu.org/copyleft/gpl.html
00019 */
00020 
00021 
00022 #ifndef BLMINIMIZE_HH
00023 #define BLMINIMIZE_HH
00024 
00025 #include <stdio.h>
00026 #include <vector>
00027 #include <string>
00028 #include <assert.h>
00029 
00030 /**     class BLMinimizeFunction is the base class for a class defining a
00031  *      function to minimize.
00032  **/
00033 class BLMinimizeFunction {
00034 public:
00035         BLMinimizeFunction() { }
00036         virtual ~BLMinimizeFunction() { }
00037 
00038         /// operator() is the function to minimize. Note it is not const
00039         /// for flexibility in derived classes.
00040         /// Its argument is necessarily const.
00041         virtual double operator()(const std::vector<double> &x) = 0;
00042 };
00043 
00044 #ifdef USE_GSL_MULTIMIN
00045 
00046 #include <gsl/gsl_multimin.h>
00047 
00048 /**     class BLMinimize implements a general function minimizer using
00049  *      the GSL multimin functions.
00050  *
00051  *      Uses the the Simplex algorithm of Nelder and Mead.
00052  *
00053  *      Compilation for simple programs:
00054  *              g++ -I... test.cc -lgsl -lgslcblas -lm
00055  *
00056  *      NOTE: the parameters to the minimizer must all have scales that are
00057  *      comparable, or the minimizer will often converge to a wildly incorrect
00058  *      value. This applies to track fitting, where positions have scales of
00059  *      1 (millimeters), but angles have scales of 0.001 (milliradians).
00060  *      This class scales the parameters appropriately. it also permits the
00061  *      function to have arguments that are held constant during the 
00062  *      minimization by setting their scale to 0.0.
00063  *
00064  *      NOTE: only a single instance of BLMinimize can be active at a time
00065  *      (i.e. inside the minimize() function). Any number of instances can
00066  *      exist. This basically means your BLMinimizeFunction cannot itself
00067  *      use a minimizer internally.
00068  **/
00069 class BLMinimize {
00070         std::vector<double> scale;
00071         std::vector<double> arg;
00072         std::vector<std::string> name;
00073         double val;
00074         double size;
00075         int iter;
00076         int print;
00077         double chisqMin;
00078         BLMinimizeFunction *func;
00079         gsl_multimin_function ff;
00080         gsl_multimin_fminimizer *m;
00081 private:
00082         // f() is the function the gsl routines call (they are in C)
00083         static double f(const gsl_vector *x, void *param) {
00084                 BLMinimize *p = (BLMinimize*)param;
00085                 return p->call(x);
00086         }
00087         // call() will call the user's function, scaling the argument vector.
00088         double call(const gsl_vector *x) {
00089                 if(print >= 2) printf("BLMinimize call args:");
00090                 int j=0;
00091                 for(unsigned i=0; i<scale.size(); ++i) {
00092                         if(scale[i] != 0.0) // fixed args are already set
00093                                 arg[i] = gsl_vector_get(x,j++) * scale[i];
00094                         if(print >= 2)
00095                                 printf(" %.4g",arg[i]);
00096                 }
00097                 if(print >= 2) printf("\n");
00098                 double v = (*func)(arg);
00099                 if(print >= 2) printf("BLMinimize call returns value=%.4g\n",v);
00100                 return v;
00101         }
00102 public:
00103         /// Constructor.
00104         BLMinimize(int _print=2) : scale(), arg(), name()
00105                 { val=GSL_NAN; size=GSL_NAN; iter=0; print=_print;
00106                   chisqMin=GSL_NAN; func=0; ff.f=0; ff.n=0; ff.params=0;  m=0;
00107                 }
00108 
00109         virtual ~BLMinimize() { if(m) gsl_multimin_fminimizer_free(m); }
00110 
00111         /// setFunc() will set the function to minimize and the scale of its
00112         /// arguments. The number of non-zero entries in _scale[] determines
00113         /// the number of variables to minimize. names is a colon-separated
00114         /// list of parameter names.
00115         void setFunc(BLMinimizeFunction &f, const std::vector<double> &_scale,
00116                                                 const std::string &names="") {
00117                 func = &f;
00118                 scale.clear();
00119                 arg.clear();
00120                 name.clear();
00121                 std::vector<std::string> n = splitString(names);
00122                 int nx = 0;
00123                 for(unsigned i=0; i<_scale.size(); ++i) {
00124                         double s=_scale[i];
00125                         scale.push_back(s);
00126                         arg.push_back(s);
00127                         name.push_back(i<n.size() ? n[i] : "");
00128                         if(s != 0.0) ++nx;
00129                 }
00130                 ff.f=BLMinimize::f;
00131                 ff.n=nx;
00132                 ff.params=this;
00133                 if(m) gsl_multimin_fminimizer_free(m);
00134                 m = 0;
00135                 if(nx == 0) return;
00136                 m = gsl_multimin_fminimizer_alloc(
00137                         gsl_multimin_fminimizer_nmsimplex,nx);
00138         }
00139 
00140         /// setPrintDetail() sets the detail of printing per iteration:
00141         /// 0=no print, 1=print value and size only, 2=print value, size, and x.
00142         /// the default is 2.
00143         void setPrintDetail(int v=0) { print = v; }
00144 
00145         /// setChisqMin() sets an additional convergence criterium of chisq
00146         /// less than the given value. Can only be used when the function is
00147         /// positive definite, usually a chi-squared. Required for track-
00148         /// fitting of test tracks that have truly 0 chisq -- the normal 
00149         /// convergence often fails for such tracks.
00150         void setChisqMin(double v) { chisqMin = v; }
00151 
00152         /// minimize() finds the minimum, returning status (0=success).
00153         ///   x         the initial and then the final argument values
00154         ///   tolerance the convergence limit on the simplex size
00155         ///             NOTE: this is after scaling!
00156         ///   maxIter   limits the number of iterations
00157         /// The initial step for x[i] is scale[i] (x[i] is scaled, so the
00158         /// step to the gsl routine is 1.0).
00159         int minimize(std::vector<double> &x, double tolerance=1e-3,
00160                                                         int maxIter=100) {
00161                 if(ff.n <= 0) return GSL_EINVAL;
00162                 assert(x.size()==scale.size());
00163                 gsl_vector *init = gsl_vector_alloc(ff.n);
00164                 gsl_vector *istep = gsl_vector_alloc(ff.n);
00165                 unsigned j=0;
00166                 for(unsigned i=0; i<scale.size(); ++i) {
00167                         arg[i] = x[i]; // fixed args, variable ones overwritten
00168                         if(scale[i] != 0.0) {
00169                                 gsl_vector_set(init,j,x[i]/scale[i]); 
00170                                 gsl_vector_set(istep,j,1.0);
00171                                 ++j;
00172                         }
00173                 }
00174                 assert(j==ff.n); assert(m!=0);
00175                 gsl_multimin_fminimizer_set(m,&ff,init,istep);
00176                 for(iter=1; iter<=maxIter; ++iter) {
00177                         int status = gsl_multimin_fminimizer_iterate(m);
00178                         if(status) {
00179                                 if(print >= 1)
00180                                         printf("BLMinimize FAILED status=%d\n",
00181                                                                         status);
00182                                 gsl_vector_free(init);
00183                                 gsl_vector_free(istep);
00184                                 return status;
00185                         }
00186                         val = m->fval;
00187                         size = gsl_multimin_fminimizer_size(m);
00188                         if(print >= 1)
00189                             printf("BLMinimize Iter %d: val=%.4g size=%.4g\n",
00190                                                         iter,val,size);
00191                         status = gsl_multimin_test_size(size,tolerance);
00192                         if(status == GSL_SUCCESS || (!gsl_isnan(chisqMin) &&
00193                            val <= chisqMin)) {
00194                                 unsigned j=0;
00195                                 for(unsigned i=0; i<scale.size(); ++i) {
00196                                     if(scale[i] != 0.0)
00197                                         x[i] = gsl_vector_get(m->x,j++) *
00198                                                                 scale[i];
00199                                 }
00200                                 gsl_vector_free(init);
00201                                 gsl_vector_free(istep);
00202                                 if(print >= 1) printf("BLMinimize SUCCESS\n");
00203                                 return 0;
00204                         }
00205                 }
00206                 gsl_vector_free(init);
00207                 gsl_vector_free(istep);
00208                 if(print >= 1) printf("BLMinimize FAILED iter > %d\n",maxIter);
00209                 return GSL_EMAXITER;
00210         }
00211 
00212         /// value() returns the function value from the previous minimize().
00213         double value() const { return val; }
00214 
00215         /// getSize() returns the Simplex size from the previous minimize().
00216         double getSize() { return size; }
00217 
00218         /// getIterations() returns the # iterations of the previous minimize().
00219         int getIterations() { return iter; }
00220 
00221 private:
00222         std::vector<std::string> splitString(const std::string &s) {
00223                 std::string tmp(s);
00224                 std::vector<std::string> v;
00225                 while(tmp.size() > 0) {
00226                         unsigned i=tmp.find_first_of(":");
00227                         if(i > tmp.size()) i = tmp.size();
00228                         v.push_back(tmp.substr(0,i));
00229                         tmp.erase(0,i+1);
00230                 }
00231                 return v;
00232         }
00233 };
00234 
00235 #endif // USE_GSL_MULTIMIN
00236 
00237 #ifdef USE_TMINUIT
00238 
00239 #include "TMinuit.h"
00240 
00241 /**     class BLMinimize implements a general function minimizer using
00242  *      the TMinuit class from Root.
00243  *
00244  *      NOTE: only a single instance of BLMinimize can be active at a time
00245  *      (i.e. inside the minimize() function). Any number of instances can
00246  *      exist. This basically means your BLMinimizeFunction cannot itself
00247  *      use a minimizer internally.
00248  **/
00249 class BLMinimize {
00250         std::vector<double> scale;
00251         std::vector<std::string> name;
00252         std::vector<double> arg;
00253         double val;
00254         int iter;
00255         int print;
00256         double chisqMin;
00257         BLMinimizeFunction *func;
00258         TMinuit *tm;
00259         static BLMinimize *current;
00260 private:
00261         static void fcn(Int_t &npar, Double_t *grad, Double_t &val, Double_t *x,
00262                                                                 Int_t iflag) {
00263                 assert(current != 0);
00264                 val = current->call(npar,x);
00265         }
00266         double call(int npar, double *x) {
00267                 // npar seems to be 1 too small assert(npar == (int)arg.size());
00268                 for(unsigned i=0; i<arg.size(); ++i)
00269                         arg[i] = x[i];
00270                 return (*func)(arg);
00271         }
00272 public:
00273         /// Constructor.
00274         BLMinimize(int _print=2) : scale(), name(), arg()
00275                 { val=0.0; iter=0; print=_print; chisqMin=0.0; func=0; tm=0; }
00276 
00277         virtual ~BLMinimize() {  }
00278 
00279         /// setFunc() will set the function to minimize and the scale of its
00280         /// arguments. The number of non-zero entries in _scale[] determines
00281         /// the number of variables to minimize. names is a colon-separated
00282         /// list of parameter names.
00283         void setFunc(BLMinimizeFunction &f, const std::vector<double> &_scale,
00284                                                 const std::string &names="") {
00285                 func = &f;
00286                 scale.clear();
00287                 name.clear();
00288                 arg.clear();
00289                 std::vector<std::string> n = splitString(names);
00290                 for(unsigned i=0; i<_scale.size(); ++i) {
00291                         double s=_scale[i];
00292                         scale.push_back(s);
00293                         name.push_back(i<n.size() ? n[i] : "");
00294                         arg.push_back(0.0);
00295                 }
00296                 
00297                 if(tm) delete tm;
00298                 tm = new TMinuit();
00299                 tm->SetFCN(fcn);
00300                 tm->SetPrintLevel(-1);
00301                 tm->Command("SET PRINTOUT -1");
00302         }
00303 
00304         /// setPrintDetail() sets the detail of printing per iteration:
00305         /// 0=no print, 1=print value and size only, 2=print value, size, and x.
00306         /// the default is 2.
00307         void setPrintDetail(int v=0) { print = v; }
00308 
00309         /// setChisqMin() sets an additional convergence criterium of chisq
00310         /// less than the given value. Can only be used when the function is
00311         /// positive definite, usually a chi-squared. Required for track-
00312         /// fitting of test tracks that have truly 0 chisq -- the normal 
00313         /// convergence often fails for such tracks.
00314         void setChisqMin(double v) { chisqMin = v; }
00315 
00316         /// minimize() finds the minimum, returning status (0=success).
00317         ///   x         the initial and then the final argument values
00318         ///   tolerance the convergence limit on the simplex size
00319         ///             NOTE: this is after scaling!
00320         ///   maxIter   limits the number of iterations
00321         /// The initial step for x[i] is scale[i].
00322         int minimize(std::vector<double> &x, double tolerance=1e-3,
00323                                                         int maxIter=100) {
00324                 assert(tm != 0);
00325                 assert(x.size() == scale.size());
00326                 for(unsigned i=0; i<scale.size(); ++i) {
00327                         tm->DefineParameter(i,name[i].c_str(),x[i],scale[i],0.0,0.0);
00328                 }
00329                 //@ What do I do with tolerance????
00330                 tm->SetMaxIterations(maxIter);
00331                 current = this;
00332                 int retval = tm->Migrad();
00333                 current = 0;
00334                 for(unsigned i=0; i<scale.size(); ++i) {
00335                         double err;
00336                         tm->GetParameter(i,x[i],err);
00337                 }
00338                 val = (*func)(x);
00339                 return retval;
00340         }
00341 
00342         /// value() returns the function value from the previous minimize().
00343         double value() const { return val; }
00344 
00345         /// getSize() returns the Simplex size from the previous minimize().
00346         double getSize() { return 0.0; }
00347 
00348         /// getIterations() returns the # iterations of the previous minimize().
00349         int getIterations() { return iter; }
00350 
00351 private:
00352         std::vector<std::string> splitString(const std::string &s) {
00353                 std::string tmp(s);
00354                 std::vector<std::string> v;
00355                 while(tmp.size() > 0) {
00356                         unsigned i=tmp.find_first_of(":");
00357                         if(i > tmp.size()) i = tmp.size();
00358                         v.push_back(tmp.substr(0,i));
00359                         tmp.erase(0,i+1);
00360                 }
00361                 return v;
00362         }
00363 };
00364 
00365 
00366 #endif // USE_TMINUIT
00367 
00368 #endif // BLMINIMIZE_HH
g4beamline