00001
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"
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
00041
00042
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
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
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
00111
00112
00113 #ifdef DEBUG
00114 printf("fit function %s\n",f1->GetName() );
00115 #endif
00116
00117
00118 int hdim = HFit::GetDimension(h1);
00119 int iret = HFit::CheckFitFunction(f1, hdim);
00120 if (iret != 0) return iret;
00121
00122
00123
00124
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
00137 if (fitOption.Bound || fitOption.Like || fitOption.Errors || fitOption.Gradient || fitOption.More || fitOption.User|| fitOption.Integral || fitOption.Minuit)
00138 linear = kFALSE;
00139
00140
00141
00142 std::auto_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter() );
00143 ROOT::Fit::FitConfig & fitConfig = fitter->Config();
00144
00145
00146 ROOT::Fit::DataOptions opt;
00147 opt.fIntegral = fitOption.Integral;
00148 opt.fUseRange = fitOption.Range;
00149 if (fitOption.Like) opt.fUseEmpty = true;
00150 if (linear) opt.fCoordErrors = false;
00151 if (fitOption.NoErrX) opt.fCoordErrors = false;
00152 if (fitOption.W1) opt.fErrors1 = true;
00153 if (fitOption.W1 > 1) opt.fUseEmpty = true;
00154
00155
00156
00157 if (opt.fUseRange) {
00158 #ifdef DEBUG
00159 printf("use range \n" );
00160 #endif
00161
00162 Double_t fxmin, fymin, fzmin, fxmax, fymax, fzmax;
00163 f1->GetRange(fxmin, fymin, fzmin, fxmax, fymax, fzmax);
00164
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
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
00193 if (special != 0 && !fitOption.Bound && !linear) {
00194 if (special == 100) ROOT::Fit::InitGaus (*fitdata,f1);
00195 else if (special == 110) ROOT::Fit::Init2DGaus(*fitdata,f1);
00196 else if (special == 400) ROOT::Fit::InitGaus (*fitdata,f1);
00197 else if (special == 410) ROOT::Fit::Init2DGaus(*fitdata,f1);
00198
00199 else if (special == 200) ROOT::Fit::InitExpo (*fitdata, f1);
00200
00201 }
00202
00203
00204
00205
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
00212 if (fitdata->GetErrorType() == ROOT::Fit::BinData::kNoError) fitConfig.SetNormErrors(true);
00213
00214 if (int(fitdata->NDim()) == hdim -1 ) fitConfig.SetNormErrors(true);
00215
00216
00217
00218
00219
00220
00221
00222 for (int i = 0; i < npar; ++i) {
00223 ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
00224
00225
00226 double plow,pup;
00227 f1->GetParLimits(i,plow,pup);
00228 if (plow*pup != 0 && plow >= pup) {
00229 parSettings.Fix();
00230 }
00231 else if (plow < pup ) {
00232 parSettings.SetLimits(plow,pup);
00233 }
00234
00235
00236
00237 double err = f1->GetParError(i);
00238 if ( err > 0)
00239 parSettings.SetStepSize(err);
00240 else if (plow < pup) {
00241 double step = 0.1 * (pup - plow);
00242
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
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265 fitConfig.SetMinimizerOptions(minOption);
00266
00267
00268 if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
00269 if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
00270
00271
00272 if (linear) {
00273 if (fitOption.Robust ) {
00274
00275 std::string type = "Robust";
00276
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);
00281 }
00282 else
00283 fitConfig.SetMinimizer("Linear","");
00284 }
00285 else {
00286 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
00287 }
00288
00289
00290
00291 if (fitOption.Errors) {
00292
00293 fitConfig.SetParabErrors(true);
00294 fitConfig.SetMinosErrors(true);
00295 }
00296
00297
00298
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
00310
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)
00319 fitok = fitter->FitFCN( userFcn );
00320 else if (fitOption.Like)
00321 fitok = fitter->LikelihoodFit(*fitdata);
00322 else
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
00333 iret = fitResult.Status();
00334 if (!fitResult.IsEmpty() ) {
00335
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
00347 if (!fitOption.Nostore) {
00348 HFit::GetDrawingRange(h1, range);
00349 HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption);
00350 }
00351
00352
00353 TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
00354
00355
00356
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
00365 if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
00366 }
00367
00368
00369 if (lastFitter) {
00370 TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter);
00371 if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) )
00372 delete lastBCFitter;
00373 }
00374
00375 TVirtualFitter::SetFitter( bcfitter );
00376
00377
00378
00379
00380
00381
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
00403
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
00448
00449 TH1 * h1 = gr->GetHistogram();
00450
00451 if (h1) HFit::GetDrawingRange(h1, range);
00452 }
00453 void HFit::GetDrawingRange(TMultiGraph * mg, ROOT::Fit::DataRange & range) {
00454
00455
00456 TH1 * h1 = mg->GetHistogram();
00457 if (h1) {
00458 HFit::GetDrawingRange(h1, range);
00459 }
00460 else if (range.Size(0) == 0) {
00461
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
00477
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
00484
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
00499
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
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
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
00570 if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) h1->Draw(goption);
00571 if (gPad) gPad->Modified();
00572
00573 return;
00574 }
00575
00576
00577 void ROOT::Fit::FitOptionsMake(const char *option, Foption_t &fitOption) {
00578
00579 Double_t h=0;
00580 TString opt = option;
00581 opt.ToUpper();
00582 opt.ReplaceAll("ROB", "H");
00583 opt.ReplaceAll("EX0", "T");
00584
00585
00586
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
00636
00637 TFitResultPtr ROOT::Fit::UnBinFit(ROOT::Fit::UnBinData * fitdata, TF1 * fitfunc, Foption_t & fitOption , const ROOT::Math::MinimizerOptions & minOption) {
00638
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
00652 std::auto_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter() );
00653 ROOT::Fit::FitConfig & fitConfig = fitter->Config();
00654
00655
00656 unsigned int dim = fitdata->NDim();
00657
00658
00659
00660
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
00669
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
00676 if (plow*pup != 0 && plow >= pup) {
00677 parSettings.Fix();
00678 }
00679 else if (plow < pup )
00680 parSettings.SetLimits(plow,pup);
00681
00682
00683
00684 double err = fitfunc->GetParError(i);
00685 if ( err > 0)
00686 parSettings.SetStepSize(err);
00687 else if (plow < pup) {
00688 double step = 0.1 * (pup - plow);
00689
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
00706 if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
00707
00708
00709 if (fitOption.Errors) {
00710
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
00720 int iret = fitResult.Status();
00721 if (!fitResult.IsEmpty() ) {
00722
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
00733 TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
00734
00735 TBackCompFitter * bcfitter = new TBackCompFitter(fitter, std::auto_ptr<ROOT::Fit::FitData>(fitdata));
00736
00737 fitdata = 0;
00738 bcfitter->SetFitOption(fitOption);
00739
00740 bcfitter->SetUserFunc(fitfunc);
00741
00742 if (lastFitter) delete lastFitter;
00743 TVirtualFitter::SetFitter( bcfitter );
00744
00745
00746
00747
00748
00749
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
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
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
00778 HFit::CheckGraphFitOptions(foption);
00779
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
00785 HFit::CheckGraphFitOptions(foption);
00786
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
00792 HFit::CheckGraphFitOptions(foption);
00793
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
00799 return HFit::Fit(s1,f1,foption,moption,goption,range);
00800 }
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814