HFitImpl.cxx

Go to the documentation of this file.
00001 // new HFit function 
00002 //______________________________________________________________________________
00003 
00004 
00005 #include "TH1.h"
00006 #include "TH2.h"
00007 #include "TF1.h"
00008 #include "TF2.h"
00009 #include "TF3.h"
00010 #include "TError.h"
00011 #include "TGraph.h"
00012 #include "TMultiGraph.h"
00013 #include "TGraph2D.h"
00014 #include "THnSparse.h"
00015 
00016 #include "Fit/Fitter.h"
00017 #include "Fit/BinData.h"
00018 #include "Fit/UnBinData.h"
00019 #include "HFitInterface.h"
00020 #include "Math/MinimizerOptions.h"
00021 
00022 #include "Math/WrappedTF1.h"
00023 #include "Math/WrappedMultiTF1.h"
00024 
00025 #include "TList.h"
00026 #include "TMath.h"
00027 
00028 #include "TClass.h"
00029 #include "TVirtualPad.h" // for gPad
00030 
00031 #include "TBackCompFitter.h"
00032 #include "TFitResultPtr.h"
00033 #include "TFitResult.h"
00034 
00035 #include <stdlib.h>
00036 #include <cmath>
00037 #include <memory>
00038 #include <limits>
00039 
00040 //#define DEBUG
00041 
00042 // utility functions used in TH1::Fit
00043 
00044 namespace HFit { 
00045 
00046    int GetDimension(const TH1 * h1) { return h1->GetDimension(); }
00047    int GetDimension(const TGraph * ) { return 1; }
00048    int GetDimension(const TMultiGraph * ) { return 1; }
00049    int GetDimension(const TGraph2D * ) { return 2; }
00050    int GetDimension(const THnSparse * s1) { return s1->GetNdimensions(); }
00051 
00052    int CheckFitFunction(const TF1 * f1, int hdim);
00053 
00054    void FitOptionsMake(const char *option, Foption_t &fitOption);
00055 
00056    void CheckGraphFitOptions(Foption_t &fitOption);
00057 
00058 
00059    void GetDrawingRange(TH1 * h1, ROOT::Fit::DataRange & range);
00060    void GetDrawingRange(TGraph * gr, ROOT::Fit::DataRange & range);
00061    void GetDrawingRange(TMultiGraph * mg, ROOT::Fit::DataRange & range);
00062    void GetDrawingRange(TGraph2D * gr, ROOT::Fit::DataRange & range);
00063    void GetDrawingRange(THnSparse * s, ROOT::Fit::DataRange & range);
00064 
00065 
00066    template <class FitObject>
00067    TFitResultPtr Fit(FitObject * h1, TF1 *f1 , Foption_t & option , const ROOT::Math::MinimizerOptions & moption, const char *goption,  ROOT::Fit::DataRange & range); 
00068 
00069    template <class FitObject>
00070    void StoreAndDrawFitFunction(FitObject * h1, const TF1 * f1, const ROOT::Fit::DataRange & range, bool, bool, const char *goption);
00071 
00072 } 
00073 
00074 int HFit::CheckFitFunction(const TF1 * f1, int dim) { 
00075    // Check validity of fitted function
00076    if (!f1) {
00077       Error("Fit", "function may not be null pointer");
00078       return -1;
00079    }
00080    if (f1->IsZombie()) {
00081       Error("Fit", "function is zombie");
00082       return -2;
00083    }
00084 
00085    int npar = f1->GetNpar();
00086    if (npar <= 0) {
00087       Error("Fit", "function %s has illegal number of parameters = %d", f1->GetName(), npar);
00088       return -3;
00089    }
00090 
00091    // Check that function has same dimension as histogram
00092    if (f1->GetNdim() > dim) {
00093       Error("Fit","function %s dimension, %d, is greater than fit object dimension, %d",
00094             f1->GetName(), f1->GetNdim(), dim);
00095       return -4;
00096    }
00097    if (f1->GetNdim() < dim-1) {
00098       Error("Fit","function %s dimension, %d, is smaller than fit object dimension -1, %d",
00099             f1->GetName(), f1->GetNdim(), dim);
00100       return -5;
00101    }
00102 
00103    return 0; 
00104 
00105 }
00106 
00107 template<class FitObject>
00108 TFitResultPtr HFit::Fit(FitObject * h1, TF1 *f1 , Foption_t & fitOption , const ROOT::Math::MinimizerOptions & minOption, const char *goption, ROOT::Fit::DataRange & range)
00109 {
00110    // perform fit of histograms, or graphs using new fitting classes 
00111    // use same routines for fitting both graphs and histograms
00112 
00113 #ifdef DEBUG
00114    printf("fit function %s\n",f1->GetName() ); 
00115 #endif
00116 
00117    // replacement function using  new fitter
00118    int hdim = HFit::GetDimension(h1); 
00119    int iret = HFit::CheckFitFunction(f1, hdim);
00120    if (iret != 0) return iret; 
00121 
00122 
00123 
00124    // integral option is not supported in this case
00125    if (f1->GetNdim() < hdim ) { 
00126       if (fitOption.Integral) Info("Fit","Ignore Integral option. Model function dimension is less than the data object dimension");
00127       if (fitOption.Like) Info("Fit","Ignore Likelihood option. Model function dimension is less than the data object dimension");
00128       fitOption.Integral = 0; 
00129       fitOption.Like = 0; 
00130    }
00131 
00132    Int_t special = f1->GetNumber();
00133    Bool_t linear = f1->IsLinear();
00134    Int_t npar = f1->GetNpar();
00135    if (special==299+npar)  linear = kTRUE;
00136    // do not use linear fitter in these case
00137    if (fitOption.Bound || fitOption.Like || fitOption.Errors || fitOption.Gradient || fitOption.More || fitOption.User|| fitOption.Integral || fitOption.Minuit)
00138       linear = kFALSE;
00139 
00140 
00141    // create the fitter
00142    std::auto_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter() );
00143    ROOT::Fit::FitConfig & fitConfig = fitter->Config();
00144 
00145    // create options 
00146    ROOT::Fit::DataOptions opt; 
00147    opt.fIntegral = fitOption.Integral; 
00148    opt.fUseRange = fitOption.Range; 
00149    if (fitOption.Like) opt.fUseEmpty = true;  // use empty bins in log-likelihood fits 
00150    if (linear) opt.fCoordErrors = false; // cannot use coordinate errors in a linear fit
00151    if (fitOption.NoErrX) opt.fCoordErrors = false;  // do not use coordinate errors when requested
00152    if (fitOption.W1) opt.fErrors1 = true;
00153    if (fitOption.W1 > 1) opt.fUseEmpty = true; // use empty bins with weight=1
00154 
00155    //opt.fBinVolume = 1; // for testing
00156 
00157    if (opt.fUseRange) { 
00158 #ifdef DEBUG
00159       printf("use range \n" ); 
00160 #endif
00161       // check if function has range 
00162       Double_t fxmin, fymin, fzmin, fxmax, fymax, fzmax;
00163       f1->GetRange(fxmin, fymin, fzmin, fxmax, fymax, fzmax);
00164       // support only one range - so add only if was not set before
00165       if (range.Size(0) == 0) range.AddRange(0,fxmin,fxmax);
00166       if (range.Size(1) == 0) range.AddRange(1,fymin,fymax);
00167       if (range.Size(2) == 0) range.AddRange(2,fzmin,fzmax);
00168    }
00169 #ifdef DEBUG
00170    printf("range  size %d\n", range.Size(0) ); 
00171    if (range.Size(0)) {
00172       double x1; double x2; range.GetRange(0,x1,x2); 
00173       printf(" range in x = [%f,%f] \n",x1,x2);
00174    }
00175 #endif
00176 
00177    // fill data  
00178    std::auto_ptr<ROOT::Fit::BinData> fitdata(new ROOT::Fit::BinData(opt,range) );
00179    ROOT::Fit::FillData(*fitdata, h1, f1); 
00180    if (fitdata->Size() == 0 ) { 
00181       Warning("Fit","Fit data is empty ");
00182       return -1;
00183    }
00184 
00185 #ifdef DEBUG
00186    printf("HFit:: data size is %d \n",fitdata->Size());
00187    for (unsigned int i = 0; i < fitdata->Size(); ++i) { 
00188       if (fitdata->NDim() == 1) printf(" x[%d] = %f - value = %f \n", i,*(fitdata->Coords(i)),fitdata->Value(i) ); 
00189    }
00190 #endif   
00191 
00192    // this functions use the TVirtualFitter
00193    if (special != 0 && !fitOption.Bound && !linear) { 
00194       if      (special == 100)      ROOT::Fit::InitGaus  (*fitdata,f1); // gaussian
00195       else if (special == 110)      ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D gaussian
00196       else if (special == 400)      ROOT::Fit::InitGaus  (*fitdata,f1); // landau (use the same)
00197       else if (special == 410)      ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D landau (use the same)
00198 
00199       else if (special == 200)      ROOT::Fit::InitExpo  (*fitdata, f1); // exponential
00200 
00201    }
00202 
00203 
00204    // set the fit function 
00205    // if option grad is specified use gradient 
00206    if ( (linear || fitOption.Gradient) )  
00207       fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*f1) );
00208    else 
00209       fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*f1) ) );
00210 
00211    // error normalization in case of zero error in the data
00212    if (fitdata->GetErrorType() == ROOT::Fit::BinData::kNoError) fitConfig.SetNormErrors(true);
00213    // normalize errors also in case you are fitting a Ndim histo with a N-1 function
00214    if (int(fitdata->NDim())  == hdim -1 ) fitConfig.SetNormErrors(true);
00215 
00216    
00217    // here need to get some static extra information (like max iterations, error def, etc...)
00218 
00219 
00220    // parameter settings and transfer the parameters values, names and limits from the functions
00221    // is done automatically in the Fitter.cxx 
00222    for (int i = 0; i < npar; ++i) { 
00223       ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i); 
00224 
00225       // check limits
00226       double plow,pup; 
00227       f1->GetParLimits(i,plow,pup);  
00228       if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
00229          parSettings.Fix();
00230       }
00231       else if (plow < pup ) { 
00232          parSettings.SetLimits(plow,pup);
00233       }
00234 
00235       // set the parameter step size (by default are set to 0.3 of value)
00236       // if function provides meaningful error values
00237       double err = f1->GetParError(i); 
00238       if ( err > 0) 
00239          parSettings.SetStepSize(err); 
00240       else if (plow < pup) { // in case of limits improve step sizes 
00241          double step = 0.1 * (pup - plow); 
00242          // check if value is not too close to limit otherwise trim value
00243          if (  parSettings.Value() < pup && pup - parSettings.Value() < 2 * step  ) 
00244             step = (pup - parSettings.Value() ) / 2; 
00245          else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step ) 
00246             step = (parSettings.Value() - plow ) / 2; 
00247          
00248          parSettings.SetStepSize(step); 
00249       }
00250       
00251       
00252    }
00253 
00254    // needed for setting precision ? 
00255    //   - Compute sum of squares of errors in the bin range
00256    // should maybe use stat[1] ??
00257  //   Double_t ey, sumw2=0;
00258 //    for (i=hxfirst;i<=hxlast;i++) {
00259 //       ey = GetBinError(i);
00260 //       sumw2 += ey*ey;
00261 //    }
00262 
00263 
00264    // set all default minimizer options (tolerance, max iterations, etc..)
00265    fitConfig.SetMinimizerOptions(minOption); 
00266 
00267    // specific  print level options 
00268    if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3); 
00269    if (fitOption.Quiet)    fitConfig.MinimizerOptions().SetPrintLevel(0); 
00270 
00271    // specific minimizer options depending on minimizer 
00272    if (linear) { 
00273       if (fitOption.Robust  ) { 
00274          // robust fitting
00275          std::string type = "Robust";
00276          // if an h is specified print out the value adding to the type 
00277          if (fitOption.hRobust > 0 && fitOption.hRobust < 1.)
00278             type += " (h=" + ROOT::Math::Util::ToString(fitOption.hRobust) + ")";
00279          fitConfig.SetMinimizer("Linear",type.c_str());
00280          fitConfig.MinimizerOptions().SetTolerance(fitOption.hRobust); // use tolerance for passing robust parameter
00281       }
00282       else 
00283          fitConfig.SetMinimizer("Linear","");
00284    }
00285    else { 
00286       if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
00287    }
00288 
00289 
00290    // check if Error option (run Hesse and Minos) then 
00291    if (fitOption.Errors) { 
00292       // run Hesse and Minos
00293       fitConfig.SetParabErrors(true);
00294       fitConfig.SetMinosErrors(true);
00295    }
00296 
00297 
00298    // do fitting 
00299 
00300 #ifdef DEBUG
00301    if (fitOption.Like)   printf("do  likelihood fit...\n");
00302    if (linear)    printf("do linear fit...\n");
00303    printf("do now  fit...\n");
00304 #endif   
00305  
00306    bool fitok = false; 
00307 
00308 
00309    // check if can use option user 
00310    //typedef  void (* MinuitFCN_t )(int &npar, double *gin, double &f, double *u, int flag);
00311    TVirtualFitter::FCNFunc_t  userFcn = 0;  
00312    if (fitOption.User && TVirtualFitter::GetFitter() ) { 
00313       userFcn = (TVirtualFitter::GetFitter())->GetFCN(); 
00314       (TVirtualFitter::GetFitter())->SetUserFunc(f1); 
00315    }
00316    
00317 
00318    if (fitOption.User && userFcn) // user provided fit objective function
00319       fitok = fitter->FitFCN( userFcn );
00320    else if (fitOption.Like) // likelihood fit 
00321       fitok = fitter->LikelihoodFit(*fitdata);
00322    else // standard least square fit
00323       fitok = fitter->Fit(*fitdata); 
00324 
00325 
00326    if ( !fitok  && !fitOption.Quiet )
00327       Warning("Fit","Abnormal termination of minimization.");
00328    iret |= !fitok; 
00329         
00330 
00331    const ROOT::Fit::FitResult & fitResult = fitter->Result(); 
00332    // one could set directly the fit result in TF1
00333    iret = fitResult.Status(); 
00334    if (!fitResult.IsEmpty() ) { 
00335       // set in f1 the result of the fit      
00336       f1->SetChisquare(fitResult.Chi2() );
00337       f1->SetNDF(fitResult.Ndf() );
00338       f1->SetNumberFitPoints(fitdata->Size() );
00339 
00340       f1->SetParameters( &(fitResult.Parameters().front()) ); 
00341       if ( int( fitResult.Errors().size()) >= f1->GetNpar() ) 
00342          f1->SetParErrors( &(fitResult.Errors().front()) ); 
00343   
00344    }
00345 
00346 //   - Store fitted function in histogram functions list and draw
00347       if (!fitOption.Nostore) {
00348          HFit::GetDrawingRange(h1, range);
00349          HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption); 
00350       }
00351 
00352       // store result in the backward compatible VirtualFitter
00353       TVirtualFitter * lastFitter = TVirtualFitter::GetFitter(); 
00354       // pass ownership of Fitter and Fitdata to TBackCompFitter (fitter pointer cannot be used afterwards)
00355       // need to get the raw pointer due to the  missing template copy ctor of auto_ptr on solaris
00356       // reset fitdata(cannot use anymore , ownership is passed)
00357       TBackCompFitter * bcfitter = new TBackCompFitter(fitter, std::auto_ptr<ROOT::Fit::FitData>(fitdata.release()));
00358       bcfitter->SetFitOption(fitOption); 
00359       bcfitter->SetObjectFit(h1);
00360       bcfitter->SetUserFunc(f1);
00361       bcfitter->SetBit(TBackCompFitter::kCanDeleteLast);
00362       if (userFcn) { 
00363          bcfitter->SetFCN(userFcn); 
00364          // for interpreted FCN functions
00365          if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
00366       }
00367          
00368       // delete last fitter if it has been created here before
00369       if (lastFitter) {          
00370          TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter); 
00371          if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) ) 
00372             delete lastBCFitter; 
00373       }
00374       //N.B=  this might create a memory leak if user does not delete the fitter he creates
00375       TVirtualFitter::SetFitter( bcfitter ); 
00376 
00377       // print results
00378 //       if (!fitOption.Quiet) fitResult.Print(std::cout);
00379 //       if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout); 
00380 
00381       // use old-style for printing the results
00382       if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
00383       else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
00384 
00385       if (fitOption.StoreResult) 
00386       {
00387          TFitResult* fr = new TFitResult(fitResult);
00388          TString name = "TFitResult-";
00389          name = name + h1->GetName() + "-" + f1->GetName();
00390          TString title = "TFitResult-";
00391          title += h1->GetTitle();
00392          fr->SetName(name);
00393          fr->SetTitle(title);
00394          return TFitResultPtr(fr);
00395       }
00396       else 
00397          return TFitResultPtr(iret);
00398 }
00399 
00400 
00401 void HFit::GetDrawingRange(TH1 * h1, ROOT::Fit::DataRange & range) { 
00402    // get range from histogram and update the DataRange class  
00403    // if a ranges already exist in that dimension use that one
00404 
00405    Int_t ndim = GetDimension(h1);
00406 
00407    double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0; 
00408    if (range.Size(0) == 0) { 
00409       TAxis  & xaxis = *(h1->GetXaxis()); 
00410       Int_t hxfirst = xaxis.GetFirst();
00411       Int_t hxlast  = xaxis.GetLast();
00412       Double_t binwidx = xaxis.GetBinWidth(hxlast);
00413       xmin    = xaxis.GetBinLowEdge(hxfirst);
00414       xmax    = xaxis.GetBinLowEdge(hxlast) +binwidx;
00415       range.AddRange(xmin,xmax);
00416    } 
00417 
00418    if (ndim > 1) {
00419       if (range.Size(1) == 0) { 
00420          TAxis  & yaxis = *(h1->GetYaxis()); 
00421          Int_t hyfirst = yaxis.GetFirst();
00422          Int_t hylast  = yaxis.GetLast();
00423          Double_t binwidy = yaxis.GetBinWidth(hylast);
00424          ymin    = yaxis.GetBinLowEdge(hyfirst);
00425          ymax    = yaxis.GetBinLowEdge(hylast) +binwidy;
00426          range.AddRange(1,ymin,ymax);
00427       }
00428    }      
00429    if (ndim > 2) {
00430       if (range.Size(2) == 0) { 
00431          TAxis  & zaxis = *(h1->GetZaxis()); 
00432          Int_t hzfirst = zaxis.GetFirst();
00433          Int_t hzlast  = zaxis.GetLast();
00434          Double_t binwidz = zaxis.GetBinWidth(hzlast);
00435          zmin    = zaxis.GetBinLowEdge(hzfirst);
00436          zmax    = zaxis.GetBinLowEdge(hzlast) +binwidz;
00437          range.AddRange(2,zmin,zmax);
00438       }
00439    }      
00440 #ifdef DEBUG
00441    std::cout << "xmin,xmax" << xmin << "  " << xmax << std::endl;
00442 #endif
00443 
00444 }
00445 
00446 void HFit::GetDrawingRange(TGraph * gr,  ROOT::Fit::DataRange & range) { 
00447    // get range for graph (used sub-set histogram)
00448    // N.B. : this is different than in previous implementation of TGraph::Fit where range used was from xmin to xmax.
00449    TH1 * h1 = gr->GetHistogram();
00450    // an histogram is normally always returned for a TGraph
00451    if (h1) HFit::GetDrawingRange(h1, range);
00452 }
00453 void HFit::GetDrawingRange(TMultiGraph * mg,  ROOT::Fit::DataRange & range) { 
00454    // get range for multi-graph (used sub-set histogram)
00455    // N.B. : this is different than in previous implementation of TMultiGraph::Fit where range used was from data xmin to xmax.
00456    TH1 * h1 = mg->GetHistogram();
00457    if (h1) {
00458       HFit::GetDrawingRange(h1, range);
00459    } 
00460    else if (range.Size(0) == 0) { 
00461       // compute range from all the TGraph's belonging to the MultiGraph
00462       double xmin = std::numeric_limits<double>::infinity(); 
00463       double xmax = -std::numeric_limits<double>::infinity(); 
00464       TIter next(mg->GetListOfGraphs() );
00465       TGraph * g = 0; 
00466       while (  (g = (TGraph*) next() ) ) { 
00467          double x1 = 0, x2 = 0, y1 = 0, y2 = 0;
00468          g->ComputeRange(x1,y1,x2,y2); 
00469          if (x1 < xmin) xmin = x1; 
00470          if (x2 > xmax) xmax = x2; 
00471       }
00472       range.AddRange(xmin,xmax);
00473    }
00474 }
00475 void HFit::GetDrawingRange(TGraph2D * gr,  ROOT::Fit::DataRange & range) { 
00476    // get range for graph2D (used sub-set histogram)
00477    // N.B. : this is different than in previous implementation of TGraph2D::Fit. There range used was always(0,0)
00478    TH1 * h1 = gr->GetHistogram();
00479    if (h1) HFit::GetDrawingRange(h1, range);
00480 }
00481 
00482 void HFit::GetDrawingRange(THnSparse * s1, ROOT::Fit::DataRange & range) { 
00483    // get range from histogram and update the DataRange class  
00484    // if a ranges already exist in that dimension use that one
00485 
00486    Int_t ndim = GetDimension(s1);
00487 
00488    for ( int i = 0; i < ndim; ++i ) {
00489       if ( range.Size(i) == 0 ) {
00490          TAxis *axis = s1->GetAxis(i);
00491          range.AddRange(i, axis->GetXmin(), axis->GetXmax());
00492       }
00493    }
00494 }
00495 
00496 template<class FitObject>
00497 void HFit::StoreAndDrawFitFunction(FitObject * h1, const TF1 * f1, const ROOT::Fit::DataRange & range, bool delOldFunction, bool drawFunction, const char *goption) { 
00498 //   - Store fitted function in histogram functions list and draw
00499 // should have separate functions for 1,2,3d ? t.b.d in case
00500 
00501 #ifdef DEBUG
00502    std::cout <<"draw and store fit function " << f1->GetName() << std::endl;
00503 #endif
00504  
00505    TF1 *fnew1;
00506    TF2 *fnew2;
00507    TF3 *fnew3;
00508 
00509    Int_t ndim = GetDimension(h1);
00510    double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0; 
00511    if (range.Size(0) ) range.GetRange(0,xmin,xmax);
00512    if (range.Size(1) ) range.GetRange(1,ymin,ymax);
00513    if (range.Size(2) ) range.GetRange(2,zmin,zmax);
00514 
00515 
00516 #ifdef DEBUG
00517    std::cout <<"draw and store fit function " << f1->GetName() 
00518              << " Range in x = [ " << xmin << " , " << xmax << " ]" << std::endl;
00519 #endif
00520 
00521    TList * funcList = h1->GetListOfFunctions();
00522    if (funcList == 0){
00523       Error("StoreAndDrawFitFunction","Function list has not been created - cannot store the fitted function");
00524       return;
00525    } 
00526 
00527    if (delOldFunction) {
00528       TIter next(funcList, kIterBackward);
00529       TObject *obj;
00530       while ((obj = next())) {
00531          if (obj->InheritsFrom(TF1::Class())) {
00532             funcList->Remove(obj);
00533             delete obj;
00534          }
00535       }
00536    }
00537 
00538    // copy TF1 using TClass to avoid slicing in case of derived classes
00539    if (ndim < 2) {
00540       fnew1 = (TF1*)f1->IsA()->New();
00541       f1->Copy(*fnew1);
00542       funcList->Add(fnew1);
00543       fnew1->SetParent( h1 );
00544       fnew1->SetRange(xmin,xmax);
00545       fnew1->Save(xmin,xmax,0,0,0,0);
00546       if (!drawFunction) fnew1->SetBit(TF1::kNotDraw);
00547       fnew1->SetBit(TFormula::kNotGlobal);
00548    } else if (ndim < 3) {
00549       fnew2 = (TF2*)f1->IsA()->New();
00550       f1->Copy(*fnew2);
00551       funcList->Add(fnew2);
00552       fnew2->SetRange(xmin,ymin,xmax,ymax);
00553       fnew2->SetParent( h1 );
00554       fnew2->Save(xmin,xmax,ymin,ymax,0,0);
00555       if (!drawFunction) fnew2->SetBit(TF1::kNotDraw);
00556       fnew2->SetBit(TFormula::kNotGlobal);
00557    } else {
00558       // 3D- why f3d is not saved ???
00559       fnew3 = (TF3*)f1->IsA()->New();
00560       f1->Copy(*fnew3);
00561       funcList->Add(fnew3);
00562       fnew3->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
00563       fnew3->SetParent( h1 );
00564       fnew3->Save(xmin,xmax,ymin,ymax,zmin,zmax);
00565       if (!drawFunction) fnew3->SetBit(TF1::kNotDraw);
00566       fnew3->SetBit(TFormula::kNotGlobal);
00567    }
00568    if (h1->TestBit(kCanDelete)) return;
00569    // draw only in case of histograms
00570    if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) h1->Draw(goption);
00571    if (gPad) gPad->Modified(); // this is not in TH1 code (needed ??)
00572    
00573    return; 
00574 }
00575 
00576 
00577 void ROOT::Fit::FitOptionsMake(const char *option, Foption_t &fitOption) { 
00578    //   - Decode list of options into fitOption (used by the TGraph)
00579    Double_t h=0;
00580    TString opt = option;
00581    opt.ToUpper();
00582    opt.ReplaceAll("ROB", "H");
00583    opt.ReplaceAll("EX0", "T");
00584 
00585    //for robust fitting, see if # of good points is defined
00586    // decode parameters for robust fitting
00587    if (opt.Contains("H=0.")) {
00588       int start = opt.Index("H=0.");
00589       int numpos = start + strlen("H=0.");
00590       int numlen = 0;
00591       int len = opt.Length();
00592       while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
00593       TString num = opt(numpos,numlen);
00594       opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
00595       h = atof(num.Data());
00596       h*=TMath::Power(10, -numlen);
00597    }
00598 
00599    if (opt.Contains("U")) fitOption.User    = 1;
00600    if (opt.Contains("Q")) fitOption.Quiet   = 1;
00601    if (opt.Contains("V")){fitOption.Verbose = 1; fitOption.Quiet   = 0;}
00602    if (opt.Contains("L")) fitOption.Like    = 1;
00603    if (opt.Contains("X")) fitOption.Chi2    = 1;
00604    if (opt.Contains("I")) fitOption.Integral= 1;
00605    if (opt.Contains("LL")) fitOption.Like   = 2;
00606    if (opt.Contains("W")) fitOption.W1      = 1;
00607    if (opt.Contains("E")) fitOption.Errors  = 1;
00608    if (opt.Contains("R")) fitOption.Range   = 1;
00609    if (opt.Contains("G")) fitOption.Gradient= 1;
00610    if (opt.Contains("M")) fitOption.More    = 1;
00611    if (opt.Contains("N")) fitOption.Nostore = 1;
00612    if (opt.Contains("0")) fitOption.Nograph = 1;
00613    if (opt.Contains("+")) fitOption.Plus    = 1;
00614    if (opt.Contains("B")) fitOption.Bound   = 1;
00615    if (opt.Contains("C")) fitOption.Nochisq = 1;
00616    if (opt.Contains("F")) fitOption.Minuit  = 1;
00617    if (opt.Contains("T")) fitOption.NoErrX   = 1;
00618    if (opt.Contains("S")) fitOption.StoreResult   = 1;
00619    if (opt.Contains("H")) { fitOption.Robust  = 1;   fitOption.hRobust = h; } 
00620 
00621 }
00622 
00623 void HFit::CheckGraphFitOptions(Foption_t & foption) { 
00624    if (foption.Like) { 
00625       Info("CheckGraphFitOptions","L (Log Likelihood fit) is an invalid option when fitting a graph. It is ignored");
00626       foption.Like = 0; 
00627    }
00628    if (foption.Integral) { 
00629       Info("CheckGraphFitOptions","I (use function integral) is an invalid option when fitting a graph. It is ignored");
00630       foption.Integral = 0; 
00631    }
00632    return;
00633 } 
00634 
00635 // implementation of unbin fit function (defined in HFitInterface)
00636 
00637 TFitResultPtr ROOT::Fit::UnBinFit(ROOT::Fit::UnBinData * fitdata, TF1 * fitfunc, Foption_t & fitOption , const ROOT::Math::MinimizerOptions & minOption) { 
00638    // do unbin fit, ownership of fitdata is passed later to the TBackFitter class
00639 
00640 #ifdef DEBUG
00641    printf("tree data size is %d \n",fitdata->Size());
00642    for (unsigned int i = 0; i < fitdata->Size(); ++i) { 
00643       if (fitdata->NDim() == 1) printf(" x[%d] = %f \n", i,*(fitdata->Coords(i) ) ); 
00644    }
00645 #endif   
00646    if (fitdata->Size() == 0 ) { 
00647       Warning("Fit","Fit data is empty ");
00648       return -1;
00649    }
00650       
00651    // create the fitter
00652    std::auto_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter() );
00653    ROOT::Fit::FitConfig & fitConfig = fitter->Config();
00654 
00655    // dimension is given by data because TF1 pointer can have wrong one
00656    unsigned int dim = fitdata->NDim();    
00657 
00658    // set the fit function
00659    // if option grad is specified use gradient 
00660    // need to create a wrapper for an automatic  normalized TF1 ???
00661    if ( fitOption.Gradient ) {
00662       assert ( (int) dim == fitfunc->GetNdim() );
00663       fitter->SetFunction(ROOT::Math::WrappedTF1(*fitfunc) );
00664    }
00665    else 
00666       fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*fitfunc, dim) ) );
00667 
00668    // parameter setting is done automaticaly in the Fitter class 
00669    // need only to set limits
00670    int npar = fitfunc->GetNpar();
00671    for (int i = 0; i < npar; ++i) { 
00672       ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i); 
00673       double plow,pup; 
00674       fitfunc->GetParLimits(i,plow,pup);  
00675       // this is a limitation of TF1 interface - cannot fix a parameter to zero value
00676       if (plow*pup != 0 && plow >= pup) {
00677          parSettings.Fix();
00678       }
00679       else if (plow < pup ) 
00680          parSettings.SetLimits(plow,pup);
00681 
00682       // set the parameter step size (by default are set to 0.3 of value)
00683       // if function provides meaningful error values
00684       double err = fitfunc->GetParError(i); 
00685       if ( err > 0) 
00686          parSettings.SetStepSize(err); 
00687       else if (plow < pup) { // in case of limits improve step sizes 
00688          double step = 0.1 * (pup - plow); 
00689          // check if value is not too close to limit otherwise trim value
00690          if (  parSettings.Value() < pup && pup - parSettings.Value() < 2 * step  ) 
00691             step = (pup - parSettings.Value() ) / 2; 
00692          else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step ) 
00693             step = (parSettings.Value() - plow ) / 2; 
00694          
00695          parSettings.SetStepSize(step); 
00696       }
00697 
00698    }
00699 
00700    fitConfig.SetMinimizerOptions(minOption); 
00701 
00702    if (fitOption.Verbose)   fitConfig.MinimizerOptions().SetPrintLevel(3); 
00703    if (fitOption.Quiet)     fitConfig.MinimizerOptions().SetPrintLevel(0); 
00704   
00705    // more 
00706    if (fitOption.More)   fitConfig.SetMinimizer("Minuit","MigradImproved");
00707 
00708    // chech if Minos or more options
00709    if (fitOption.Errors) { 
00710       // run Hesse and Minos
00711       fitConfig.SetParabErrors(true);
00712       fitConfig.SetMinosErrors(true);
00713    }
00714    
00715    bool fitok = false; 
00716    fitok = fitter->Fit(*fitdata); 
00717  
00718    const ROOT::Fit::FitResult & fitResult = fitter->Result(); 
00719    // one could set directly the fit result in TF1
00720    int iret = fitResult.Status(); 
00721    if (!fitResult.IsEmpty() ) { 
00722       // set in fitfunc the result of the fit      
00723       fitfunc->SetNDF(fitResult.Ndf() );
00724       fitfunc->SetNumberFitPoints(fitdata->Size() );
00725 
00726       fitfunc->SetParameters( &(fitResult.Parameters().front()) ); 
00727       if ( int( fitResult.Errors().size()) >= fitfunc->GetNpar() ) 
00728          fitfunc->SetParErrors( &(fitResult.Errors().front()) ); 
00729   
00730    }
00731 
00732    // store result in the backward compatible VirtualFitter
00733    TVirtualFitter * lastFitter = TVirtualFitter::GetFitter(); 
00734    // pass ownership of Fitter and Fitdata to TBackCompFitter (fitter pointer cannot be used afterwards)
00735    TBackCompFitter * bcfitter = new TBackCompFitter(fitter, std::auto_ptr<ROOT::Fit::FitData>(fitdata));
00736  // cannot use anymore now fitdata (given away ownership)
00737    fitdata = 0;
00738    bcfitter->SetFitOption(fitOption); 
00739    //bcfitter->SetObjectFit(fTree);
00740    bcfitter->SetUserFunc(fitfunc);
00741    
00742    if (lastFitter) delete lastFitter; 
00743    TVirtualFitter::SetFitter( bcfitter ); 
00744    
00745    // print results
00746 //       if (!fitOption.Quiet) fitResult.Print(std::cout);
00747 //       if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout); 
00748    
00749    // use old-style for printing the results
00750    if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
00751    else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
00752 
00753    if (fitOption.StoreResult) 
00754    {
00755       TFitResult* fr = new TFitResult(fitResult);
00756       TString name = "TFitResult-";
00757       name = name + "UnBinData-" + fitfunc->GetName();
00758       TString title = "TFitResult-";
00759       title += name;
00760       fr->SetName(name);
00761       fr->SetTitle(title);
00762       return TFitResultPtr(fr);
00763    }
00764    else 
00765       return TFitResultPtr(iret);
00766 }
00767 
00768 
00769 // implementations of ROOT::Fit::FitObject functions (defined in HFitInterface) in terms of the template HFit::Fit
00770 
00771 TFitResultPtr ROOT::Fit::FitObject(TH1 * h1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) { 
00772    // histogram fitting
00773    return HFit::Fit(h1,f1,foption,moption,goption,range); 
00774 }
00775 
00776 TFitResultPtr ROOT::Fit::FitObject(TGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) { 
00777   // exclude options not valid for graphs
00778    HFit::CheckGraphFitOptions(foption);
00779     // TGraph fitting
00780    return HFit::Fit(gr,f1,foption,moption,goption,range); 
00781 }
00782 
00783 TFitResultPtr ROOT::Fit::FitObject(TMultiGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) { 
00784   // exclude options not valid for graphs
00785    HFit::CheckGraphFitOptions(foption);
00786     // TMultiGraph fitting
00787    return HFit::Fit(gr,f1,foption,moption,goption,range); 
00788 }
00789 
00790 TFitResultPtr ROOT::Fit::FitObject(TGraph2D * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) { 
00791   // exclude options not valid for graphs
00792    HFit::CheckGraphFitOptions(foption);
00793     // TGraph2D fitting
00794    return HFit::Fit(gr,f1,foption,moption,goption,range); 
00795 }
00796 
00797 TFitResultPtr ROOT::Fit::FitObject(THnSparse * s1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) { 
00798    // sparse histogram fitting
00799    return HFit::Fit(s1,f1,foption,moption,goption,range); 
00800 }
00801 
00802 
00803 
00804 // Int_t TGraph2D::DoFit(TF2 *f2 ,Option_t *option ,Option_t *goption) { 
00805 //    // internal graph2D fitting methods
00806 //    Foption_t fitOption;
00807 //    ROOT::Fit::FitOptionsMake(option,fitOption);
00808 
00809 //    // create range and minimizer options with default values 
00810 //    ROOT::Fit::DataRange range(2); 
00811 //    ROOT::Math::MinimizerOptions minOption; 
00812 //    return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range); 
00813 // }
00814 

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