TMinuitMinimizer.cxx

Go to the documentation of this file.
00001 // @(#)root/minuit:$Id: TMinuitMinimizer.cxx 37461 2010-12-09 22:19:18Z moneta $
00002 // Author: L. Moneta Wed Oct 25 16:28:55 2006
00003 
00004 /**********************************************************************
00005  *                                                                    *
00006  * Copyright (c) 2006  LCG ROOT Math Team, CERN/PH-SFT                *
00007  *                                                                    *
00008  *                                                                    *
00009  **********************************************************************/
00010 
00011 // Implementation file for class TMinuitMinimizer
00012 
00013 #include "TMinuitMinimizer.h"
00014 #include "Math/IFunction.h"
00015 
00016 #include "TMinuit.h"
00017 #include "TROOT.h"
00018 
00019 #include "TGraph.h" // needed for scan 
00020 #include "TError.h"
00021 
00022 #include <iostream>
00023 #include <cassert>
00024 #include <algorithm>
00025 #include <functional>
00026 #include <cmath>
00027 
00028 //______________________________________________________________________________
00029 //
00030 //  TMinuitMinimizer class implementing the ROOT::Math::Minimizer interface using 
00031 //  TMinuit. 
00032 //  This class is normally instantiates using the plug-in manager 
00033 //  (plug-in with name Minuit or TMinuit)
00034 //  In addition the user can choose the minimizer algorithm: Migrad (the default one), Simplex, or Minimize (combined Migrad + Simplex)
00035 //  
00036 //__________________________________________________________________________________________
00037 
00038 // initialize the static instances
00039 
00040 
00041 ROOT::Math::IMultiGenFunction * TMinuitMinimizer::fgFunc = 0; 
00042 TMinuit * TMinuitMinimizer::fgMinuit = 0;
00043 bool TMinuitMinimizer::fgUsed = false; 
00044 bool TMinuitMinimizer::fgUseStaticMinuit = true;   // default case use static Minuit instance
00045 
00046 ClassImp(TMinuitMinimizer)
00047 
00048 
00049 TMinuitMinimizer::TMinuitMinimizer(ROOT::Minuit::EMinimizerType type, unsigned int ndim ) : 
00050    fUsed(false),
00051    fMinosRun(false),
00052    fDim(ndim),
00053    fStrategy(1),
00054    fType(type), 
00055    fMinuit(0)
00056 {
00057    // Constructor for TMinuitMinimier class via an enumeration specifying the minimization 
00058    // algorithm type. Supported types are : kMigrad, kSimplex, kCombined (a combined 
00059    // Migrad + Simplex minimization) and kMigradImproved (a Migrad mininimization folloed by an 
00060    // improved search for global minima). The default type is Migrad (kMigrad). 
00061 
00062    // initialize if npar is given
00063    if (fDim > 0) InitTMinuit(fDim);   
00064 }
00065 
00066 TMinuitMinimizer::TMinuitMinimizer(const char *  type, unsigned int ndim ) : 
00067    fUsed(false),
00068    fMinosRun(false),
00069    fDim(ndim),
00070    fStrategy(1),
00071    fMinuit(0)
00072 {
00073    // constructor from a char * for the algorithm type, used by the plug-in manager
00074    // The names supported (case unsensitive) are: 
00075    //  Migrad (default), Simplex, Minimize (for the combined Migrad+ Simplex) and Migrad_imp 
00076 
00077    // select type from the string
00078    std::string algoname(type);
00079    std::transform(algoname.begin(), algoname.end(), algoname.begin(), (int(*)(int)) tolower ); 
00080 
00081    ROOT::Minuit::EMinimizerType algoType = ROOT::Minuit::kMigrad; 
00082    if (algoname == "simplex")   algoType = ROOT::Minuit::kSimplex; 
00083    if (algoname == "minimize" ) algoType = ROOT::Minuit::kCombined; 
00084    if (algoname == "migradimproved" ) algoType = ROOT::Minuit::kMigradImproved; 
00085    if (algoname == "scan" )           algoType = ROOT::Minuit::kScan; 
00086    if (algoname == "seek" )           algoType = ROOT::Minuit::kSeek; 
00087 
00088    fType = algoType; 
00089 
00090    // initialize if npar is given
00091    if (fDim > 0) InitTMinuit(fDim);
00092 
00093 }
00094 
00095 TMinuitMinimizer::~TMinuitMinimizer() 
00096 {
00097    // Destructor implementation.
00098    if (fMinuit && !fgUseStaticMinuit) { 
00099       delete fMinuit; 
00100       fgMinuit = 0;
00101    }
00102 }
00103 
00104 TMinuitMinimizer::TMinuitMinimizer(const TMinuitMinimizer &) : 
00105    Minimizer()
00106 {
00107    // Implementation of copy constructor (it is private).
00108 }
00109 
00110 TMinuitMinimizer & TMinuitMinimizer::operator = (const TMinuitMinimizer &rhs) 
00111 {
00112    // Implementation of assignment operator (private)
00113    if (this == &rhs) return *this;  // time saving self-test
00114    return *this;
00115 }
00116 
00117 bool TMinuitMinimizer::UseStaticMinuit(bool on ) { 
00118    // static method to control usage of global TMinuit instance
00119    bool prev = fgUseStaticMinuit; 
00120    fgUseStaticMinuit = on; 
00121    return prev; 
00122 }
00123 
00124 void TMinuitMinimizer::InitTMinuit(int dim) {
00125 
00126    // when called a second time check dimension - create only if needed
00127    // initialize the minuit instance - recreating a new one if needed 
00128    if (fMinuit ==0 ||  dim > fMinuit->fMaxpar) { 
00129 
00130       // case not using the global instance - recreate it all the time 
00131       if (fgUseStaticMinuit) { 
00132 
00133          // re-use gMinuit as static instance of TMinuit
00134          // which can be accessed by the user after minimization
00135          // check if fgMinuit is different than gMinuit
00136          // case 1: fgMinuit not zero but fgMinuit has been deleted (not in gROOT): set to zero
00137          // case 2: fgMinuit not zero and exists in global list  : set fgMinuit to gMinuit 
00138          // case 3: fgMinuit zero - and gMinuit not zero: create a new instance locally to avoid conflict
00139          if (fgMinuit != gMinuit) { 
00140             // if object exists in gROOT remove it to avoid a memory leak 
00141             if (fgMinuit ) { 
00142                if (gROOT->GetListOfSpecials()->FindObject(fgMinuit) == 0) { 
00143                   // case 1: object does not exists in gROOT - means it has been deleted
00144                   fgMinuit = 0; 
00145                }
00146                else {
00147                   // case 2: object exists - but gMinuit points to something else
00148                   // restore gMinuit to the one used before by TMinuitMinimizer
00149                   gMinuit = fgMinuit; 
00150                }
00151             }
00152             else {
00153                // case 3: avoid reusing existing one - mantain fgMinuit to zero
00154                // otherwise we will get a double delete if user deletes externally gMinuit
00155                // in this case we will loose gMinuit instance
00156 //                fgMinuit = gMinuit;
00157 //                fgUsed = true;  // need to reset in case  other gMinuit instance is later used
00158             }
00159          }
00160          
00161          // check if need to create a new TMinuit instance
00162          if (fgMinuit == 0) {
00163             fgUsed = false;
00164             fgMinuit =  new TMinuit(dim);
00165          }   
00166          else if (fgMinuit->GetNumPars() != int(dim) ) { 
00167             delete fgMinuit; 
00168             fgUsed = false;
00169             fgMinuit =  new TMinuit(dim);
00170          }
00171       
00172          fMinuit = fgMinuit; 
00173       }
00174       
00175       else { 
00176          // re- create all the time a new instance of TMinuit (fgUseStaticMinuit is false)
00177          if (fMinuit) delete fMinuit; 
00178          fMinuit =  new TMinuit(dim);
00179          fgMinuit = fMinuit; 
00180          fgUsed = false; 
00181       }
00182 
00183    }  // endif fMinuit ==0 || dim > fMaxpar
00184 
00185    fDim = dim; 
00186 
00187    R__ASSERT(fMinuit);
00188 
00189    // set print level in TMinuit
00190    double arglist[1];
00191    int ierr= 0; 
00192    // TMinuit level is shift by 1 -1 means 0;
00193    arglist[0] = PrintLevel() - 1;
00194    fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00195 }
00196 
00197 
00198 void TMinuitMinimizer::SetFunction(const  ROOT::Math::IMultiGenFunction & func) { 
00199    // Set the objective function to be minimized, by passing a function object implement the 
00200    // basic multi-dim Function interface. In this case the derivatives will be 
00201    // calculated by Minuit 
00202    // Here a TMinuit instance is created since only at this point we know the number of parameters 
00203 
00204 
00205    fDim = func.NDim(); 
00206 
00207    // create TMinuit if needed
00208    InitTMinuit(fDim); 
00209    
00210    // assign to the static pointer (NO Thread safety here)
00211    fgFunc = const_cast<ROOT::Math::IMultiGenFunction *>(&func); 
00212    fMinuit->SetFCN(&TMinuitMinimizer::Fcn);
00213 
00214    // switch off gradient calculations
00215    double arglist[1]; 
00216    int ierr = 0;
00217    fMinuit->mnexcm("SET NOGrad",arglist,0,ierr);
00218 }
00219 
00220 void TMinuitMinimizer::SetFunction(const  ROOT::Math::IMultiGradFunction & func) { 
00221    // Set the objective function to be minimized, by passing a function object implement the 
00222    // multi-dim gradient Function interface. In this case the function derivatives are provided  
00223    // by the user via this interface and there not calculated by Minuit. 
00224 
00225    fDim = func.NDim(); 
00226 
00227    // create TMinuit if needed
00228    InitTMinuit(fDim); 
00229    
00230    // assign to the static pointer (NO Thread safety here)
00231    fgFunc = const_cast<ROOT::Math::IMultiGradFunction *>(&func); 
00232    fMinuit->SetFCN(&TMinuitMinimizer::FcnGrad);
00233 
00234    // set gradient 
00235    // by default do not check gradient calculation 
00236    // it cannot be done here, check can be done only after having defined the parameters
00237    double arglist[1]; 
00238    int ierr = 0;
00239    arglist[0] = 1; 
00240    fMinuit->mnexcm("SET GRAD",arglist,1,ierr);
00241 }
00242 
00243 void TMinuitMinimizer::Fcn( int &, double * , double & f, double * x , int /* iflag */) { 
00244    // implementation of FCN static function used internally by TMinuit.
00245    // Adapt IMultiGenFunction interface to TMinuit FCN static function
00246    f = fgFunc->operator()(x);
00247 }
00248 
00249 void TMinuitMinimizer::FcnGrad( int &, double * g, double & f, double * x , int iflag ) { 
00250    // implementation of FCN static function used internally by TMinuit.
00251    // Adapt IMultiGradFunction interface to TMinuit FCN static function in the case of user 
00252    // provided gradient.
00253    ROOT::Math::IMultiGradFunction * gFunc = dynamic_cast<ROOT::Math::IMultiGradFunction *> ( fgFunc); 
00254 
00255    assert(gFunc != 0);
00256    f = gFunc->operator()(x);
00257 
00258    // calculates also derivatives 
00259    if (iflag == 2) gFunc->Gradient(x,g);
00260 }
00261 
00262 bool TMinuitMinimizer::SetVariable(unsigned int ivar, const std::string & name, double val, double step) { 
00263    // set a free variable.
00264    if (fMinuit == 0) { 
00265       Error("SetVariable","invalid TMinuit pointer. Need to call first SetFunction"); 
00266       return false; 
00267    }
00268 
00269    fUsed = fgUsed; 
00270 
00271    // clear after minimization when setting params
00272    if (fUsed) DoClear(); 
00273 
00274    // check if parameter was defined and in case it was fixed, release it 
00275    DoReleaseFixParameter(ivar);
00276 
00277    int iret = fMinuit->DefineParameter(ivar , name.c_str(), val, step, 0., 0. ); 
00278    return (iret == 0); 
00279 }
00280 
00281 bool TMinuitMinimizer::SetLimitedVariable(unsigned int ivar, const std::string & name, double val, double step, double lower, double upper) { 
00282    // set a limited variable.
00283    if (fMinuit == 0) { 
00284       Error("SetVariable","invalid TMinuit pointer. Need to call first SetFunction"); 
00285       return false; 
00286    }
00287 
00288    fUsed = fgUsed; 
00289 
00290    // clear after minimization when setting params
00291    if (fUsed) DoClear(); 
00292 
00293    // check if parameter was defined and in case it was fixed, release it 
00294    DoReleaseFixParameter(ivar);
00295 
00296    int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, lower, upper ); 
00297    return (iret == 0); 
00298 }
00299 #ifdef LATER
00300 bool Minuit2Minimizer::SetLowerLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double lower ) {
00301     // add a lower bounded variable as a double bound one, using a very large number for the upper limit
00302 
00303    if (fMinuit == 0) { 
00304       Error("SetVariable","invalid TMinuit pointer. Need to call first SetFunction"); 
00305       return false; 
00306    }
00307 
00308    fUsed = fgUsed; 
00309 
00310    // clear after minimization when setting params
00311    if (fUsed) DoClear(); 
00312 
00313    // check if parameter was defined and in case it was fixed, release it 
00314    DoReleaseFixParameter(ivar);
00315 
00316    double s = val-lower; 
00317    double upper = s*1.0E15; 
00318    if (s != 0)  upper = 1.0E15;
00319    return SetLimitedVariable(ivar, name, val, step, lower,upper);
00320 }
00321 #endif
00322 
00323 
00324 bool TMinuitMinimizer::SetFixedVariable(unsigned int ivar, const std::string & name, double val) { 
00325    // set a fixed variable.
00326    if (fMinuit == 0) { 
00327       Error("SetVariable","invalid TMinuit pointer. Need to call first SetFunction"); 
00328       return false; 
00329    }
00330 
00331    // clear after minimization when setting params
00332    fUsed = fgUsed; 
00333 
00334    // clear after minimization when setting params
00335    if (fUsed) DoClear(); 
00336 
00337    // put an arbitrary step (0.1*abs(value) otherwise TMinuit consider the parameter as constant
00338    // constant parameters are treated differently (they are ignored inside TMinuit and not considered in the
00339    // total list of parameters) 
00340    double step = ( val != 0) ? 0.1 * std::abs(val) : 0.1;
00341    int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, 0., 0. ); 
00342    if (iret == 0) iret = fMinuit->FixParameter(ivar);
00343    return (iret == 0); 
00344 }
00345 
00346 bool TMinuitMinimizer::SetVariableValue(unsigned int ivar, double val) { 
00347    // set the value of an existing variable
00348    // parameter must exist or return false
00349 
00350    if (fMinuit == 0) { 
00351       Error("SetVariable","invalid TMinuit pointer. Need to call first SetFunction"); 
00352       return false; 
00353    }
00354 
00355    double arglist[2]; 
00356    int ierr = 0; 
00357    arglist[0] = ivar+1;  // TMinuit starts from 1 
00358    arglist[1] = val;
00359    fMinuit->mnexcm("SET PAR",arglist,2,ierr);
00360    return (ierr==0);
00361 }
00362 
00363 std::string TMinuitMinimizer::VariableName(unsigned int ivar) const { 
00364    // return the variable name
00365    if (!fMinuit || (int) ivar > fMinuit->fNu) return std::string();
00366    return std::string(fMinuit->fCpnam[ivar]);
00367 }
00368 
00369 int TMinuitMinimizer::VariableIndex(const std::string & ) const { 
00370    // return variable index
00371    Error("TMinuitMinimizer::VariableIndex"," find index of a variable from its name  is not implemented in TMinuit");
00372    return -1;
00373 }
00374 
00375 bool TMinuitMinimizer::Minimize() { 
00376    // perform the minimization using the algorithm chosen previously by the user  
00377    // By default Migrad is used. 
00378    // Return true if the found minimum is valid and update internal chached values of 
00379    // minimum values, errors and covariance matrix. 
00380    // Status of minimizer is set to: 
00381    // migradResult + 10*minosResult + 100*hesseResult + 1000*improveResult
00382 
00383 
00384    if (fMinuit == 0) { 
00385       Error("Minimize","invalid TMinuit pointer. Need to call first SetFunction and SetVariable"); 
00386       return false; 
00387    }
00388 
00389 
00390    // total number of parameter defined in Minuit is fNu
00391    if (fMinuit->fNu <  static_cast<int>(fDim) ) { 
00392       Error("Minimize","The total number of defined parameters is different than the function dimension, npar = %d, dim = %d",fMinuit->fNu, fDim);
00393       return false; 
00394    }
00395 
00396    double arglist[10]; 
00397    int ierr = 0; 
00398 
00399 
00400    // set error and print level 
00401    arglist[0] = ErrorDef(); 
00402    fMinuit->mnexcm("SET Err",arglist,1,ierr);
00403 
00404    int printlevel = PrintLevel(); 
00405    arglist[0] = printlevel - 1;
00406    fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00407 
00408    // suppress warning in case Printlevel() == 0 
00409    if (printlevel == 0)    fMinuit->mnexcm("SET NOW",arglist,0,ierr);
00410 
00411    // set precision if needed
00412    if (Precision() > 0)  { 
00413       arglist[0] = Precision();
00414       fMinuit->mnexcm("SET EPS",arglist,1,ierr);
00415    }
00416 
00417    arglist[0] = MaxFunctionCalls(); 
00418    arglist[1] = Tolerance(); 
00419    
00420    int nargs = 2; 
00421 
00422    switch (fType){  
00423    case ROOT::Minuit::kMigrad: 
00424       // case of Migrad 
00425       fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
00426       break; 
00427    case ROOT::Minuit::kCombined: 
00428       // case of combined (Migrad+ simplex)
00429       fMinuit->mnexcm("MINIMIZE",arglist,nargs,ierr);
00430       break; 
00431    case ROOT::Minuit::kSimplex: 
00432       // case of Simlex
00433       fMinuit->mnexcm("SIMPLEX",arglist,nargs,ierr);
00434       break; 
00435    case ROOT::Minuit::kScan: 
00436       // case of Scan (scan all parameters with default values)
00437       nargs = 0; 
00438       fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
00439       break; 
00440    case ROOT::Minuit::kSeek: 
00441       // case of Seek (random find minimum in a hypercube around current parameter values
00442       // use Tolerance as measures for standard deviation (if < 1) used default value in Minuit ( supposed to be  3)
00443       nargs = 1; 
00444       if (arglist[1] >= 1.) nargs = 2; 
00445       fMinuit->mnexcm("SEEK",arglist,nargs,ierr);
00446       break; 
00447    default: 
00448       // default: use Migrad 
00449       fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
00450 
00451    }
00452 
00453    fgUsed = true; 
00454    fUsed = true;
00455 
00456    fStatus = ierr; 
00457    int minErrStatus = ierr;
00458 
00459    if (printlevel>2) Info("Minimize","Finished to run MIGRAD - status %d",ierr);
00460 
00461    // run improved if needed
00462    if (ierr == 0 && fType == ROOT::Minuit::kMigradImproved) {
00463       fMinuit->mnexcm("IMPROVE",arglist,1,ierr);
00464       fStatus += 1000*ierr; 
00465       if (printlevel>2) Info("Minimize","Finished to run IMPROVE - status %d",ierr);
00466    }
00467 
00468 
00469    // check if Hesse needs to be run 
00470    if (ierr == 0 && IsValidError() ) { 
00471       fMinuit->mnexcm("HESSE",arglist,1,ierr);
00472       fStatus += 100*ierr; 
00473       if (printlevel>2) Info("Minimize","Finished to run HESSE - status %d",ierr);
00474    }
00475 
00476    // retrieve parameters and errors  from TMinuit
00477    RetrieveParams(); 
00478 
00479    if (minErrStatus == 0) { 
00480 
00481       // store global min results (only if minimization is OK)
00482       // ignore cases when Hesse or IMprove return error different than zero
00483       RetrieveErrorMatrix(); 
00484 
00485       // need to re-run Minos again if requested
00486       fMinosRun = false; 
00487 
00488       return true;
00489 
00490    }
00491    return false; 
00492 
00493 }
00494 
00495 void TMinuitMinimizer::RetrieveParams() {
00496    // retrieve from TMinuit minimum parameter values 
00497    // and errors
00498 
00499    assert(fMinuit != 0); 
00500 
00501    // get parameter values 
00502    if (fParams.size() != fDim) fParams.resize( fDim); 
00503    if (fErrors.size() != fDim) fErrors.resize( fDim); 
00504    for (unsigned int i = 0; i < fDim; ++i) { 
00505       fMinuit->GetParameter( i, fParams[i], fErrors[i]);
00506    }
00507 }
00508 
00509 void TMinuitMinimizer::RetrieveErrorMatrix() {
00510    // get covariance error matrix from TMinuit
00511    // when some parameters are fixed filled the corresponding rows and column with zero's
00512 
00513    assert(fMinuit != 0); 
00514 
00515    unsigned int nfree = NFree();
00516 
00517    unsigned int ndim2 = fDim*fDim; 
00518    if (fCovar.size() != ndim2 )  fCovar.resize(fDim*fDim); 
00519    if (nfree >= fDim) { // no fixed parameters 
00520       fMinuit->mnemat(&fCovar.front(), fDim); 
00521    } 
00522    else { 
00523       // case of fixed params need to take care 
00524       std::vector<double> tmpMat(nfree*nfree); 
00525       fMinuit->mnemat(&tmpMat.front(), nfree); 
00526       
00527       unsigned int l = 0; 
00528       for (unsigned int i = 0; i < fDim; ++i) { 
00529          
00530          if ( fMinuit->fNiofex[i] > 0 ) {  // not fixed ?
00531             unsigned int m = 0; 
00532             for (unsigned int j = 0; j <= i; ++j) { 
00533                if ( fMinuit->fNiofex[j] > 0 ) {  //not fixed
00534                   fCovar[i*fDim + j] = tmpMat[l*nfree + m];
00535                   fCovar[j*fDim + i] = fCovar[i*fDim + j]; 
00536                   m++;
00537                }
00538             }
00539             l++;
00540          }
00541       }
00542       
00543    }
00544 }
00545 
00546 unsigned int TMinuitMinimizer::NCalls() const { 
00547    // return total number of function calls
00548    if (fMinuit == 0) return 0; 
00549    return fMinuit->fNfcn;
00550 }
00551 
00552 double TMinuitMinimizer::MinValue() const { 
00553    // return minimum function value
00554 
00555    // use part of code from mnstat
00556    if (!fMinuit) return 0; 
00557    double minval = fMinuit->fAmin; 
00558    if (minval == fMinuit->fUndefi) return 0; 
00559    return minval; 
00560 }
00561 
00562 double TMinuitMinimizer::Edm() const { 
00563    // return expected distance from the minimum
00564 
00565    // use part of code from mnstat
00566    if (!fMinuit) return -1; 
00567    if (fMinuit->fAmin == fMinuit->fUndefi || fMinuit->fEDM == fMinuit->fBigedm) return fMinuit->fUp; 
00568    return fMinuit->fEDM; 
00569 }
00570 
00571 unsigned int TMinuitMinimizer::NFree() const { 
00572     // return number of free parameters 
00573    if (!fMinuit) return 0; 
00574    if (fMinuit->fNpar < 0) return 0; 
00575    return fMinuit->fNpar; 
00576 }
00577 
00578 int TMinuitMinimizer::CovMatrixStatus() const { 
00579    // return status of covariance matrix 
00580    //           status:  0= not calculated at all
00581    //                    1= approximation only, not accurate
00582    //                    2= full matrix, but forced positive-definite
00583    //                    3= full accurate covariance matrix
00584 
00585    // use part of code from mnstat
00586    if (!fMinuit) return 0; 
00587    if (fMinuit->fAmin == fMinuit->fUndefi) return 0; 
00588    return fMinuit->fISW[1];
00589 }
00590 
00591 double TMinuitMinimizer::GlobalCC(unsigned int i) const { 
00592    // global correlation coefficient for parameter i 
00593    if (!fMinuit) return 0; 
00594    if (!fMinuit->fGlobcc) return 0;
00595    if (int(i) >= fMinuit->fNu) return 0; 
00596    // get internal number in Minuit
00597    int iin = fMinuit->fNiofex[i];  
00598    // index in TMinuit starts from 1 
00599    if (iin < 1) return 0; 
00600    return fMinuit->fGlobcc[iin-1];   
00601 }
00602 
00603 bool TMinuitMinimizer::GetMinosError(unsigned int i, double & errLow, double & errUp, int ) { 
00604    // Perform Minos analysis for the given parameter  i 
00605 
00606    if (fMinuit == 0) { 
00607       Error("GetMinosError","invalid TMinuit pointer. Need to call first SetFunction and SetVariable"); 
00608       return false; 
00609    }
00610 
00611    double arglist[2];
00612    int ierr = 0; 
00613 
00614    // if Minos is not run run it 
00615    if (!fMinosRun) { 
00616 
00617       // set error and print level 
00618       arglist[0] = ErrorDef(); 
00619       fMinuit->mnexcm("SET Err",arglist,1,ierr);
00620       
00621       arglist[0] = PrintLevel()-1; 
00622       fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00623 
00624       // suppress warning in case Printlevel() == 0 
00625       if (PrintLevel() == 0)    fMinuit->mnexcm("SET NOW",arglist,0,ierr);
00626 
00627       // set precision if needed
00628       if (Precision() > 0)  { 
00629          arglist[0] = Precision();
00630          fMinuit->mnexcm("SET EPS",arglist,1,ierr);
00631       }
00632 
00633    }
00634 
00635    // syntax of MINOS is MINOS [maxcalls] [parno]
00636    // if parno = 0 all parameters are done 
00637    arglist[0] = MaxFunctionCalls(); 
00638    arglist[1] = i+1;  // par number starts from 1 in TMInuit
00639    
00640    int nargs = 2; 
00641    fMinuit->mnexcm("MINOS",arglist,nargs,ierr);
00642    fStatus += 10*ierr;
00643 
00644    fMinosRun = true; 
00645 
00646    double errParab = 0; 
00647    double gcor = 0; 
00648    // what returns if parameter fixed or constant or at limit ? 
00649    fMinuit->mnerrs(i,errUp,errLow, errParab, gcor); 
00650 
00651    if (fStatus%100 != 0 ) return false; 
00652    return true; 
00653 
00654 }
00655 
00656 void TMinuitMinimizer::DoClear() { 
00657    // reset TMinuit
00658 
00659    fMinuit->mncler();
00660    
00661    //reset the internal Minuit random generator to its initial state
00662    double val = 3;
00663    int inseed = 12345;
00664    fMinuit->mnrn15(val,inseed);
00665 
00666    fUsed = false; 
00667    fgUsed = false; 
00668 
00669 }
00670 
00671 void TMinuitMinimizer::DoReleaseFixParameter(int ivar) { 
00672    // check if a parameter is defined and in case it was fixed released
00673    // TMinuit is not able to release free parameters by redefining them
00674    // so we need to force the release
00675    if (fMinuit == 0) return; 
00676    if (fMinuit->GetNumFixedPars() == 0) return; 
00677    // check if parameter has already been defined
00678    if (int(ivar) >= fMinuit->GetNumPars() ) return;  
00679     
00680    // check if parameter is fixed
00681    for (int i = 0; i < fMinuit->fNpfix; ++i) { 
00682       if (fMinuit->fIpfix[i] == ivar+1 ) { 
00683          // parameter is fixed
00684          fMinuit->Release(ivar);
00685          return; 
00686       }
00687    }
00688     
00689 }
00690 
00691 
00692 void TMinuitMinimizer::PrintResults() { 
00693    // print-out results using classic Minuit format (mnprin)
00694    if (fMinuit == 0) return; 
00695 
00696    // print minimizer result
00697    if (PrintLevel() > 2) 
00698       fMinuit->mnprin(4,fMinuit->fAmin);
00699    else
00700       fMinuit->mnprin(3,fMinuit->fAmin);
00701 }
00702 
00703 
00704 bool TMinuitMinimizer::Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double * x, double * y) {
00705    // contour plot for parameter i and j
00706    // need a valid FunctionMinimum otherwise exits
00707    if (fMinuit == 0) { 
00708       Error("TMinuitMinimizer::Contour"," invalid TMinuit instance");
00709       return false;
00710    }
00711 
00712    // set error and print level 
00713    double arglist[1];
00714    int ierr = 0; 
00715    arglist[0] = ErrorDef(); 
00716    fMinuit->mnexcm("SET Err",arglist,1,ierr);
00717       
00718    arglist[0] = PrintLevel()-1; 
00719    fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00720 
00721    // suppress warning in case Printlevel() == 0 
00722    if (PrintLevel() == 0)    fMinuit->mnexcm("SET NOW",arglist,0,ierr);
00723 
00724    // set precision if needed
00725    if (Precision() > 0)  { 
00726       arglist[0] = Precision();
00727       fMinuit->mnexcm("SET EPS",arglist,1,ierr);
00728    }
00729 
00730 
00731    if (npoints < 4) { 
00732       Error("Contour","Cannot make contour with so few points");
00733       return false; 
00734    }
00735    int npfound = 0; 
00736    npoints -= 1;   // remove always one point in TMinuit
00737    // parameter numbers in mncont start from zero
00738    fMinuit->mncont( ipar,jpar,npoints, x, y,npfound); 
00739    if (npfound<4) {
00740       // mncont did go wrong
00741       Error("Contour","Cannot find more than 4 points");
00742       return false;
00743    }
00744    if (npfound!=(int)npoints) {
00745       // mncont did go wrong
00746       Warning("Contour","Returning only %d points ",npfound);
00747       npoints = npfound;
00748    }
00749    return true; 
00750    
00751 }
00752 
00753 bool TMinuitMinimizer::Scan(unsigned int ipar, unsigned int & nstep, double * x, double * y, double xmin, double xmax) { 
00754    // scan a parameter (variable) around the minimum value
00755    // the parameters must have been set before 
00756    // if xmin=0 && xmax == 0  by default scan around 2 sigma of the error
00757    // if the errors  are also zero then scan from min and max of parameter range
00758    // (if parameters are limited Minuit scan from min and max instead of 2 sigma by default)
00759    // (force in that case to use errors)
00760 
00761    // scan is not implemented for TMinuit, the way to return the array is only via the graph 
00762    if (fMinuit == 0) { 
00763       Error("TMinuitMinimizer::Scan"," invalid TMinuit instance");
00764       return false;
00765    }
00766 
00767    // case of default xmin and xmax
00768    if (xmin >= xmax && (int) ipar < fMinuit->GetNumPars() ) { 
00769       double val = 0; double err = 0;
00770       TString name; 
00771       double xlow = 0; double xup = 0 ;
00772       int iuint = 0; 
00773       // in mnpout index starts from ze
00774       fMinuit->mnpout( ipar, name, val, err, xlow, xup, iuint);
00775       // redefine 2 sigma for all parameters by default (TMinuit does 1 sigma and range if limited)
00776       if (iuint > 0 && err > 0) { 
00777          xmin = val - 2.*err; 
00778          xmax = val + 2 * err; 
00779       }
00780    }
00781    
00782    double arglist[4]; 
00783    int ierr = 0; 
00784 
00785    arglist[0] = PrintLevel()-1; 
00786    fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00787    // suppress warning in case Printlevel() == 0 
00788    if (PrintLevel() == 0)    fMinuit->mnexcm("SET NOW",arglist,0,ierr);
00789 
00790    // set precision if needed
00791    if (Precision() > 0)  { 
00792       arglist[0] = Precision();
00793       fMinuit->mnexcm("SET EPS",arglist,1,ierr);
00794    }
00795 
00796    if (nstep == 0) return false; 
00797    arglist[0] = ipar+1;  // TMinuit starts from 1 
00798    arglist[1] = nstep+2; // TMinuit deletes two points
00799    int nargs = 2; 
00800    if (xmax > xmin ) { 
00801       arglist[2] = xmin; 
00802       arglist[3] = xmax; 
00803       nargs = 4; 
00804    }
00805    fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
00806    if (ierr) { 
00807       Error("TMinuitMinimizer::Scan"," Error executing command SCAN");
00808       return false;      
00809    }
00810    // get TGraph object
00811    TGraph * gr = dynamic_cast<TGraph *>(fMinuit->GetPlot() );
00812    if (!gr) {
00813       Error("TMinuitMinimizer::Scan"," Error in returned graph object");
00814       return false;      
00815    }
00816    nstep = std::min(gr->GetN(), (int) nstep); 
00817 
00818 
00819    std::copy(gr->GetX(), gr->GetX()+nstep, x);
00820    std::copy(gr->GetY(), gr->GetY()+nstep, y);
00821    nstep = gr->GetN(); 
00822    return true; 
00823 }
00824 
00825 bool TMinuitMinimizer::Hesse() { 
00826    // perform calculation of Hessian
00827 
00828    if (fMinuit == 0) { 
00829       Error("Hesse","invalid TMinuit pointer. Need to call first SetFunction and SetVariable"); 
00830       return false; 
00831    }
00832 
00833 
00834    double arglist[10]; 
00835    int ierr = 0; 
00836 
00837    // set error and print level 
00838    arglist[0] = ErrorDef(); 
00839    fMinuit->mnexcm("SET ERR",arglist,1,ierr);
00840 
00841    int printlevel = PrintLevel(); 
00842    arglist[0] = printlevel - 1;
00843    fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
00844 
00845    // suppress warning in case Printlevel() == 0 
00846    if (printlevel == 0)    fMinuit->mnexcm("SET NOW",arglist,0,ierr);
00847 
00848    // set precision if needed
00849    if (Precision() > 0)  { 
00850       arglist[0] = Precision();
00851       fMinuit->mnexcm("SET EPS",arglist,1,ierr);
00852    }
00853 
00854    arglist[0] = MaxFunctionCalls(); 
00855 
00856    fMinuit->mnexcm("HESSE",arglist,1,ierr);
00857    fStatus += 100*ierr; 
00858 
00859    if (ierr != 0) return false;    
00860 
00861    // retrieve results (parameter and error matrix)
00862    // only if result is OK
00863 
00864    RetrieveParams(); 
00865    RetrieveErrorMatrix(); 
00866 
00867    return true;
00868 }
00869 
00870 
00871 //    } // end namespace Fit
00872 
00873 // } // end namespace ROOT
00874 

Generated on Tue Jul 5 14:36:48 2011 for ROOT_528-00b_version by  doxygen 1.5.1