00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "Riostream.h"
00013 #include "TROOT.h"
00014 #include "TMath.h"
00015 #include "TF1.h"
00016 #include "TH1.h"
00017 #include "TGraph.h"
00018 #include "TVirtualPad.h"
00019 #include "TStyle.h"
00020 #include "TRandom.h"
00021 #include "TInterpreter.h"
00022 #include "TPluginManager.h"
00023 #include "TBrowser.h"
00024 #include "TColor.h"
00025 #include "TClass.h"
00026 #include "TMethodCall.h"
00027 #include "TF1Helper.h"
00028 #include "Math/WrappedFunction.h"
00029 #include "Math/WrappedTF1.h"
00030 #include "Math/BrentRootFinder.h"
00031 #include "Math/BrentMinimizer1D.h"
00032 #include "Math/BrentMethods.h"
00033 #include "Math/Integrator.h"
00034 #include "Math/GaussIntegrator.h"
00035 #include "Math/GaussLegendreIntegrator.h"
00036 #include "Math/AdaptiveIntegratorMultiDim.h"
00037 #include "Math/RichardsonDerivator.h"
00038 #include "Math/Functor.h"
00039 #include "Fit/FitResult.h"
00040
00041
00042
00043 Bool_t TF1::fgAbsValue = kFALSE;
00044 Bool_t TF1::fgRejectPoint = kFALSE;
00045 static Double_t gErrorTF1 = 0;
00046
00047
00048 ClassImp(TF1)
00049
00050
00051 class GFunc {
00052 const TF1* fFunction;
00053 const double fY0;
00054 public:
00055 GFunc(const TF1* function , double y ):fFunction(function), fY0(y) {}
00056 double operator()(double x) const {
00057 return fFunction->Eval(x) - fY0;
00058 }
00059 };
00060
00061
00062 class GInverseFunc {
00063 const TF1* fFunction;
00064 public:
00065 GInverseFunc(const TF1* function):fFunction(function) {}
00066 double operator()(double x) const {
00067 return - fFunction->Eval(x);
00068 }
00069 };
00070
00071
00072
00073
00074 class TF1_EvalWrapper : public ROOT::Math::IGenFunction {
00075 public:
00076 TF1_EvalWrapper(TF1 * f, const Double_t * par, bool useAbsVal, Double_t n = 1, Double_t x0 = 0) :
00077 fFunc(f),
00078 fPar( ( (par) ? par : f->GetParameters() ) ),
00079 fAbsVal(useAbsVal),
00080 fN(n),
00081 fX0(x0)
00082 {
00083 fFunc->InitArgs(fX, fPar);
00084 }
00085
00086 ROOT::Math::IGenFunction * Clone() const {
00087
00088 TF1_EvalWrapper * f = new TF1_EvalWrapper( *this);
00089 f->fFunc->InitArgs(f->fX, f->fPar);
00090 return f;
00091 }
00092
00093 Double_t DoEval( Double_t x) const {
00094 fX[0] = x;
00095 Double_t fval = fFunc->EvalPar( fX, fPar);
00096 if (fAbsVal && fval < 0) return -fval;
00097 return fval;
00098 }
00099
00100 Double_t EvalFirstMom( Double_t x) {
00101 fX[0] = x;
00102 return fX[0] * TMath::Abs( fFunc->EvalPar( fX, fPar) );
00103 }
00104
00105 Double_t EvalNMom( Double_t x) const {
00106 fX[0] = x;
00107 return TMath::Power( fX[0] - fX0, fN) * TMath::Abs( fFunc->EvalPar( fX, fPar) );
00108 }
00109
00110 TF1 * fFunc;
00111 mutable Double_t fX[1];
00112 const double * fPar;
00113 Bool_t fAbsVal;
00114 Double_t fN;
00115 Double_t fX0;
00116 };
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346 TF1 *TF1::fgCurrent = 0;
00347
00348
00349
00350 TF1::TF1(): TFormula(), TAttLine(), TAttFill(), TAttMarker()
00351 {
00352
00353
00354 fXmin = 0;
00355 fXmax = 0;
00356 fNpx = 100;
00357 fType = 0;
00358 fNpfits = 0;
00359 fNDF = 0;
00360 fNsave = 0;
00361 fChisquare = 0;
00362 fIntegral = 0;
00363 fParErrors = 0;
00364 fParMin = 0;
00365 fParMax = 0;
00366 fAlpha = 0;
00367 fBeta = 0;
00368 fGamma = 0;
00369 fParent = 0;
00370 fSave = 0;
00371 fHistogram = 0;
00372 fMinimum = -1111;
00373 fMaximum = -1111;
00374 fMethodCall = 0;
00375 fCintFunc = 0;
00376 SetFillStyle(0);
00377 }
00378
00379
00380
00381 TF1::TF1(const char *name,const char *formula, Double_t xmin, Double_t xmax)
00382 :TFormula(name,formula), TAttLine(), TAttFill(), TAttMarker()
00383 {
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397 if (xmin < xmax ) {
00398 fXmin = xmin;
00399 fXmax = xmax;
00400 } else {
00401 fXmin = xmax;
00402 fXmax = xmin;
00403 }
00404 fNpx = 100;
00405 fType = 0;
00406 if (fNpar) {
00407 fParErrors = new Double_t[fNpar];
00408 fParMin = new Double_t[fNpar];
00409 fParMax = new Double_t[fNpar];
00410 for (int i = 0; i < fNpar; i++) {
00411 fParErrors[i] = 0;
00412 fParMin[i] = 0;
00413 fParMax[i] = 0;
00414 }
00415 } else {
00416 fParErrors = 0;
00417 fParMin = 0;
00418 fParMax = 0;
00419 }
00420 fChisquare = 0;
00421 fIntegral = 0;
00422 fAlpha = 0;
00423 fBeta = 0;
00424 fGamma = 0;
00425 fParent = 0;
00426 fNpfits = 0;
00427 fNDF = 0;
00428 fNsave = 0;
00429 fSave = 0;
00430 fHistogram = 0;
00431 fMinimum = -1111;
00432 fMaximum = -1111;
00433 fMethodCall = 0;
00434 fCintFunc = 0;
00435
00436 if (fNdim != 1 && xmin < xmax) {
00437 Error("TF1","function: %s/%s has %d parameters instead of 1",name,formula,fNdim);
00438 MakeZombie();
00439 }
00440
00441 if (!gStyle) return;
00442 SetLineColor(gStyle->GetFuncColor());
00443 SetLineWidth(gStyle->GetFuncWidth());
00444 SetLineStyle(gStyle->GetFuncStyle());
00445 SetFillStyle(0);
00446 }
00447
00448
00449
00450 TF1::TF1(const char *name, Double_t xmin, Double_t xmax, Int_t npar)
00451 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00452 {
00453
00454
00455
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465 fXmin = xmin;
00466 fXmax = xmax;
00467 fNpx = 100;
00468 fType = 2;
00469 if (npar > 0 ) fNpar = npar;
00470 if (fNpar) {
00471 fNames = new TString[fNpar];
00472 fParams = new Double_t[fNpar];
00473 fParErrors = new Double_t[fNpar];
00474 fParMin = new Double_t[fNpar];
00475 fParMax = new Double_t[fNpar];
00476 for (int i = 0; i < fNpar; i++) {
00477 fParams[i] = 0;
00478 fParErrors[i] = 0;
00479 fParMin[i] = 0;
00480 fParMax[i] = 0;
00481 }
00482 } else {
00483 fParErrors = 0;
00484 fParMin = 0;
00485 fParMax = 0;
00486 }
00487 fChisquare = 0;
00488 fIntegral = 0;
00489 fAlpha = 0;
00490 fBeta = 0;
00491 fGamma = 0;
00492 fParent = 0;
00493 fNpfits = 0;
00494 fNDF = 0;
00495 fNsave = 0;
00496 fSave = 0;
00497 fHistogram = 0;
00498 fMinimum = -1111;
00499 fMaximum = -1111;
00500 fMethodCall = 0;
00501 fCintFunc = 0;
00502 fNdim = 1;
00503
00504 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00505 gROOT->GetListOfFunctions()->Remove(f1old);
00506 SetName(name);
00507
00508 if (gStyle) {
00509 SetLineColor(gStyle->GetFuncColor());
00510 SetLineWidth(gStyle->GetFuncWidth());
00511 SetLineStyle(gStyle->GetFuncStyle());
00512 }
00513 SetFillStyle(0);
00514
00515 SetTitle(name);
00516 if (name) {
00517 if (*name == '*') return;
00518 fMethodCall = new TMethodCall();
00519 fMethodCall->InitWithPrototype(name,"Double_t*,Double_t*");
00520 fNumber = -1;
00521 gROOT->GetListOfFunctions()->Add(this);
00522 if (! fMethodCall->IsValid() ) {
00523 Error("TF1","No function found with the signature %s(Double_t*,Double_t*)",name);
00524 }
00525 } else {
00526 Error("TF1","requires a proper function name!");
00527 }
00528 }
00529
00530
00531
00532 TF1::TF1(const char *name,void *fcn, Double_t xmin, Double_t xmax, Int_t npar)
00533 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00534 {
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553 fXmin = xmin;
00554 fXmax = xmax;
00555 fNpx = 100;
00556 fType = 2;
00557
00558 if (npar > 0 ) fNpar = npar;
00559 if (fNpar) {
00560 fNames = new TString[fNpar];
00561 fParams = new Double_t[fNpar];
00562 fParErrors = new Double_t[fNpar];
00563 fParMin = new Double_t[fNpar];
00564 fParMax = new Double_t[fNpar];
00565 for (int i = 0; i < fNpar; i++) {
00566 fParams[i] = 0;
00567 fParErrors[i] = 0;
00568 fParMin[i] = 0;
00569 fParMax[i] = 0;
00570 }
00571 } else {
00572 fParErrors = 0;
00573 fParMin = 0;
00574 fParMax = 0;
00575 }
00576 fChisquare = 0;
00577 fIntegral = 0;
00578 fAlpha = 0;
00579 fBeta = 0;
00580 fGamma = 0;
00581 fParent = 0;
00582 fNpfits = 0;
00583 fNDF = 0;
00584 fNsave = 0;
00585 fSave = 0;
00586 fHistogram = 0;
00587 fMinimum = -1111;
00588 fMaximum = -1111;
00589 fMethodCall = 0;
00590 fCintFunc = 0;
00591 fNdim = 1;
00592
00593 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00594 gROOT->GetListOfFunctions()->Remove(f1old);
00595 SetName(name);
00596
00597 if (gStyle) {
00598 SetLineColor(gStyle->GetFuncColor());
00599 SetLineWidth(gStyle->GetFuncWidth());
00600 SetLineStyle(gStyle->GetFuncStyle());
00601 }
00602 SetFillStyle(0);
00603
00604 if (!fcn) return;
00605 const char *funcname = gCint->Getp2f2funcname(fcn);
00606 SetTitle(funcname);
00607 if (funcname) {
00608 fMethodCall = new TMethodCall();
00609 fMethodCall->InitWithPrototype(funcname,"Double_t*,Double_t*");
00610 fNumber = -1;
00611 gROOT->GetListOfFunctions()->Add(this);
00612 if (! fMethodCall->IsValid() ) {
00613 Error("TF1","No function found with the signature %s(Double_t*,Double_t*)",funcname);
00614 }
00615 } else {
00616 Error("TF1","can not find any function at the address 0x%lx. This function requested for %s",(Long_t)fcn,name);
00617 }
00618
00619
00620 }
00621
00622
00623
00624 TF1::TF1(const char *name,Double_t (*fcn)(Double_t *, Double_t *), Double_t xmin, Double_t xmax, Int_t npar)
00625 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00626 {
00627
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639 fXmin = xmin;
00640 fXmax = xmax;
00641 fNpx = 100;
00642
00643 fType = 1;
00644 fMethodCall = 0;
00645 fCintFunc = 0;
00646 fFunctor = ROOT::Math::ParamFunctor(fcn);
00647
00648 if (npar > 0 ) fNpar = npar;
00649 if (fNpar) {
00650 fNames = new TString[fNpar];
00651 fParams = new Double_t[fNpar];
00652 fParErrors = new Double_t[fNpar];
00653 fParMin = new Double_t[fNpar];
00654 fParMax = new Double_t[fNpar];
00655 for (int i = 0; i < fNpar; i++) {
00656 fParams[i] = 0;
00657 fParErrors[i] = 0;
00658 fParMin[i] = 0;
00659 fParMax[i] = 0;
00660 }
00661 } else {
00662 fParErrors = 0;
00663 fParMin = 0;
00664 fParMax = 0;
00665 }
00666 fChisquare = 0;
00667 fIntegral = 0;
00668 fAlpha = 0;
00669 fBeta = 0;
00670 fGamma = 0;
00671 fNsave = 0;
00672 fSave = 0;
00673 fParent = 0;
00674 fNpfits = 0;
00675 fNDF = 0;
00676 fHistogram = 0;
00677 fMinimum = -1111;
00678 fMaximum = -1111;
00679 fNdim = 1;
00680
00681
00682 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00683 gROOT->GetListOfFunctions()->Remove(f1old);
00684 SetName(name);
00685 gROOT->GetListOfFunctions()->Add(this);
00686
00687 if (!gStyle) return;
00688 SetLineColor(gStyle->GetFuncColor());
00689 SetLineWidth(gStyle->GetFuncWidth());
00690 SetLineStyle(gStyle->GetFuncStyle());
00691 SetFillStyle(0);
00692
00693 }
00694
00695
00696 TF1::TF1(const char *name,Double_t (*fcn)(const Double_t *, const Double_t *), Double_t xmin, Double_t xmax, Int_t npar)
00697 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00698 {
00699
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710
00711 fXmin = xmin;
00712 fXmax = xmax;
00713 fNpx = 100;
00714
00715 fType = 1;
00716 fMethodCall = 0;
00717 fCintFunc = 0;
00718 fFunctor = ROOT::Math::ParamFunctor(fcn);
00719
00720 if (npar > 0 ) fNpar = npar;
00721 if (fNpar) {
00722 fNames = new TString[fNpar];
00723 fParams = new Double_t[fNpar];
00724 fParErrors = new Double_t[fNpar];
00725 fParMin = new Double_t[fNpar];
00726 fParMax = new Double_t[fNpar];
00727 for (int i = 0; i < fNpar; i++) {
00728 fParams[i] = 0;
00729 fParErrors[i] = 0;
00730 fParMin[i] = 0;
00731 fParMax[i] = 0;
00732 }
00733 } else {
00734 fParErrors = 0;
00735 fParMin = 0;
00736 fParMax = 0;
00737 }
00738 fChisquare = 0;
00739 fIntegral = 0;
00740 fAlpha = 0;
00741 fBeta = 0;
00742 fGamma = 0;
00743 fNsave = 0;
00744 fSave = 0;
00745 fParent = 0;
00746 fNpfits = 0;
00747 fNDF = 0;
00748 fHistogram = 0;
00749 fMinimum = -1111;
00750 fMaximum = -1111;
00751 fNdim = 1;
00752
00753
00754 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00755 gROOT->GetListOfFunctions()->Remove(f1old);
00756 SetName(name);
00757 gROOT->GetListOfFunctions()->Add(this);
00758
00759 if (!gStyle) return;
00760 SetLineColor(gStyle->GetFuncColor());
00761 SetLineWidth(gStyle->GetFuncWidth());
00762 SetLineStyle(gStyle->GetFuncStyle());
00763 SetFillStyle(0);
00764
00765 }
00766
00767
00768
00769 TF1::TF1(const char *name, ROOT::Math::ParamFunctor f, Double_t xmin, Double_t xmax, Int_t npar ) :
00770 TFormula(),
00771 TAttLine(),
00772 TAttFill(),
00773 TAttMarker(),
00774 fXmin ( xmin ),
00775 fXmax ( xmax ),
00776 fNpx ( 100 ),
00777 fType ( 1 ),
00778 fNpfits ( 0 ),
00779 fNDF ( 0 ),
00780 fNsave ( 0 ),
00781 fChisquare ( 0 ),
00782 fIntegral ( 0 ),
00783 fParErrors ( 0 ),
00784 fParMin ( 0 ),
00785 fParMax ( 0 ),
00786 fSave ( 0 ),
00787 fAlpha ( 0 ),
00788 fBeta ( 0 ),
00789 fGamma ( 0 ),
00790 fParent ( 0 ),
00791 fHistogram ( 0 ),
00792 fMaximum ( -1111 ),
00793 fMinimum ( -1111 ),
00794 fMethodCall( 0 ),
00795 fCintFunc ( 0 ),
00796 fFunctor ( ROOT::Math::ParamFunctor(f) )
00797 {
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807 CreateFromFunctor(name, npar);
00808 }
00809
00810
00811
00812 void TF1::CreateFromFunctor(const char *name, Int_t npar)
00813 {
00814
00815
00816
00817
00818 fNdim = 1;
00819
00820 if (npar > 0 ) fNpar = npar;
00821 if (fNpar) {
00822 fNames = new TString[fNpar];
00823 fParams = new Double_t[fNpar];
00824 fParErrors = new Double_t[fNpar];
00825 fParMin = new Double_t[fNpar];
00826 fParMax = new Double_t[fNpar];
00827 for (int i = 0; i < fNpar; i++) {
00828 fParams[i] = 0;
00829 fParErrors[i] = 0;
00830 fParMin[i] = 0;
00831 fParMax[i] = 0;
00832 }
00833 } else {
00834 fParErrors = 0;
00835 fParMin = 0;
00836 fParMax = 0;
00837 }
00838
00839
00840 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00841 gROOT->GetListOfFunctions()->Remove(f1old);
00842 SetName(name);
00843 gROOT->GetListOfFunctions()->Add(this);
00844
00845 if (!gStyle) return;
00846 SetLineColor(gStyle->GetFuncColor());
00847 SetLineWidth(gStyle->GetFuncWidth());
00848 SetLineStyle(gStyle->GetFuncStyle());
00849 SetFillStyle(0);
00850
00851 }
00852
00853
00854 TF1::TF1(const char *name,void *ptr, Double_t xmin, Double_t xmax, Int_t npar, const char * className )
00855 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00856 {
00857
00858
00859
00860
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872 CreateFromCintClass(name, ptr, xmin, xmax, npar, className, 0 );
00873 }
00874
00875
00876 TF1::TF1(const char *name,void *ptr, void * , Double_t xmin, Double_t xmax, Int_t npar, const char * className, const char * methodName)
00877 :TFormula(), TAttLine(), TAttFill(), TAttMarker()
00878 {
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897 CreateFromCintClass(name, ptr, xmin, xmax, npar, className, methodName);
00898 }
00899
00900
00901 void TF1::CreateFromCintClass(const char *name,void *ptr, Double_t xmin, Double_t xmax, Int_t npar, const char * className, const char * methodName)
00902 {
00903
00904
00905
00906
00907
00908 fXmin = xmin;
00909 fXmax = xmax;
00910 fNpx = 100;
00911 fType = 3;
00912 if (npar > 0 ) fNpar = npar;
00913 if (fNpar) {
00914 fNames = new TString[fNpar];
00915 fParams = new Double_t[fNpar];
00916 fParErrors = new Double_t[fNpar];
00917 fParMin = new Double_t[fNpar];
00918 fParMax = new Double_t[fNpar];
00919 for (int i = 0; i < fNpar; i++) {
00920 fParams[i] = 0;
00921 fParErrors[i] = 0;
00922 fParMin[i] = 0;
00923 fParMax[i] = 0;
00924 }
00925 } else {
00926 fParErrors = 0;
00927 fParMin = 0;
00928 fParMax = 0;
00929 }
00930 fChisquare = 0;
00931 fIntegral = 0;
00932 fAlpha = 0;
00933 fBeta = 0;
00934 fGamma = 0;
00935 fParent = 0;
00936 fNpfits = 0;
00937 fNDF = 0;
00938 fNsave = 0;
00939 fSave = 0;
00940 fHistogram = 0;
00941 fMinimum = -1111;
00942 fMaximum = -1111;
00943 fMethodCall = 0;
00944 fNdim = 1;
00945
00946 TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name);
00947 gROOT->GetListOfFunctions()->Remove(f1old);
00948 SetName(name);
00949
00950 if (gStyle) {
00951 SetLineColor(gStyle->GetFuncColor());
00952 SetLineWidth(gStyle->GetFuncWidth());
00953 SetLineStyle(gStyle->GetFuncStyle());
00954 }
00955 SetFillStyle(0);
00956
00957 if (!ptr) return;
00958 fCintFunc = ptr;
00959
00960 if (!className) return;
00961
00962 TClass *cl = TClass::GetClass(className);
00963
00964 if (cl) {
00965 fMethodCall = new TMethodCall();
00966
00967
00968 if (methodName)
00969 fMethodCall->InitWithPrototype(cl,methodName,"Double_t*,Double_t*");
00970 else {
00971 fMethodCall->InitWithPrototype(cl,"operator()","Double_t*,Double_t*");
00972 if (! fMethodCall->IsValid() )
00973
00974 fMethodCall->InitWithPrototype(cl,"Eval","Double_t*,Double_t*");
00975 }
00976
00977 fNumber = -1;
00978 gROOT->GetListOfFunctions()->Add(this);
00979 if (! fMethodCall->IsValid() ) {
00980 if (methodName)
00981 Error("TF1","No function found in class %s with the signature %s(Double_t*,Double_t*)",className,methodName);
00982 else
00983 Error("TF1","No function found in class %s with the signature operator() (Double_t*,Double_t*) or Eval(Double_t*,Double_t*)",className);
00984 }
00985 } else {
00986 Error("TF1","can not find any class with name %s at the address 0x%lx",className,(Long_t)ptr);
00987 }
00988
00989
00990 }
00991
00992
00993
00994
00995 TF1& TF1::operator=(const TF1 &rhs)
00996 {
00997
00998
00999 if (this != &rhs) {
01000 rhs.Copy(*this);
01001 }
01002 return *this;
01003 }
01004
01005
01006
01007 TF1::~TF1()
01008 {
01009
01010
01011 if (fParMin) delete [] fParMin;
01012 if (fParMax) delete [] fParMax;
01013 if (fParErrors) delete [] fParErrors;
01014 if (fIntegral) delete [] fIntegral;
01015 if (fAlpha) delete [] fAlpha;
01016 if (fBeta) delete [] fBeta;
01017 if (fGamma) delete [] fGamma;
01018 if (fSave) delete [] fSave;
01019 delete fHistogram;
01020 delete fMethodCall;
01021
01022 if (fParent) fParent->RecursiveRemove(this);
01023 }
01024
01025
01026
01027 TF1::TF1(const TF1 &f1) : TFormula(), TAttLine(f1), TAttFill(f1), TAttMarker(f1)
01028 {
01029
01030
01031 fXmin = 0;
01032 fXmax = 0;
01033 fNpx = 100;
01034 fType = 0;
01035 fNpfits = 0;
01036 fNDF = 0;
01037 fNsave = 0;
01038 fChisquare = 0;
01039 fIntegral = 0;
01040 fParErrors = 0;
01041 fParMin = 0;
01042 fParMax = 0;
01043 fAlpha = 0;
01044 fBeta = 0;
01045 fGamma = 0;
01046 fParent = 0;
01047 fSave = 0;
01048 fHistogram = 0;
01049 fMinimum = -1111;
01050 fMaximum = -1111;
01051 fMethodCall = 0;
01052 fCintFunc = 0;
01053 SetFillStyle(0);
01054
01055 ((TF1&)f1).Copy(*this);
01056 }
01057
01058
01059
01060 void TF1::AbsValue(Bool_t flag)
01061 {
01062
01063
01064
01065
01066
01067 fgAbsValue = flag;
01068 }
01069
01070
01071
01072 void TF1::Browse(TBrowser *b)
01073 {
01074
01075
01076 Draw(b ? b->GetDrawOption() : "");
01077 gPad->Update();
01078 }
01079
01080
01081
01082 void TF1::Copy(TObject &obj) const
01083 {
01084
01085
01086 if (((TF1&)obj).fParMin) delete [] ((TF1&)obj).fParMin;
01087 if (((TF1&)obj).fParMax) delete [] ((TF1&)obj).fParMax;
01088 if (((TF1&)obj).fParErrors) delete [] ((TF1&)obj).fParErrors;
01089 if (((TF1&)obj).fIntegral) delete [] ((TF1&)obj).fIntegral;
01090 if (((TF1&)obj).fAlpha) delete [] ((TF1&)obj).fAlpha;
01091 if (((TF1&)obj).fBeta) delete [] ((TF1&)obj).fBeta;
01092 if (((TF1&)obj).fGamma) delete [] ((TF1&)obj).fGamma;
01093 if (((TF1&)obj).fSave) delete [] ((TF1&)obj).fSave;
01094 delete ((TF1&)obj).fHistogram;
01095 delete ((TF1&)obj).fMethodCall;
01096
01097 TFormula::Copy(obj);
01098 TAttLine::Copy((TF1&)obj);
01099 TAttFill::Copy((TF1&)obj);
01100 TAttMarker::Copy((TF1&)obj);
01101 ((TF1&)obj).fXmin = fXmin;
01102 ((TF1&)obj).fXmax = fXmax;
01103 ((TF1&)obj).fNpx = fNpx;
01104 ((TF1&)obj).fType = fType;
01105 ((TF1&)obj).fCintFunc = fCintFunc;
01106 ((TF1&)obj).fFunctor = fFunctor;
01107 ((TF1&)obj).fChisquare = fChisquare;
01108 ((TF1&)obj).fNpfits = fNpfits;
01109 ((TF1&)obj).fNDF = fNDF;
01110 ((TF1&)obj).fMinimum = fMinimum;
01111 ((TF1&)obj).fMaximum = fMaximum;
01112
01113 ((TF1&)obj).fParErrors = 0;
01114 ((TF1&)obj).fParMin = 0;
01115 ((TF1&)obj).fParMax = 0;
01116 ((TF1&)obj).fIntegral = 0;
01117 ((TF1&)obj).fAlpha = 0;
01118 ((TF1&)obj).fBeta = 0;
01119 ((TF1&)obj).fGamma = 0;
01120 ((TF1&)obj).fParent = fParent;
01121 ((TF1&)obj).fNsave = fNsave;
01122 ((TF1&)obj).fSave = 0;
01123 ((TF1&)obj).fHistogram = 0;
01124 ((TF1&)obj).fMethodCall = 0;
01125 if (fNsave) {
01126 ((TF1&)obj).fSave = new Double_t[fNsave];
01127 for (Int_t j=0;j<fNsave;j++) ((TF1&)obj).fSave[j] = fSave[j];
01128 }
01129 if (fNpar) {
01130 ((TF1&)obj).fParErrors = new Double_t[fNpar];
01131 ((TF1&)obj).fParMin = new Double_t[fNpar];
01132 ((TF1&)obj).fParMax = new Double_t[fNpar];
01133 Int_t i;
01134 for (i=0;i<fNpar;i++) ((TF1&)obj).fParErrors[i] = fParErrors[i];
01135 for (i=0;i<fNpar;i++) ((TF1&)obj).fParMin[i] = fParMin[i];
01136 for (i=0;i<fNpar;i++) ((TF1&)obj).fParMax[i] = fParMax[i];
01137 }
01138 if (fMethodCall) {
01139
01140 TMethodCall *m = new TMethodCall(*fMethodCall);
01141
01142 ((TF1&)obj).fMethodCall = m;
01143 }
01144 }
01145
01146
01147
01148 Double_t TF1::Derivative(Double_t x, Double_t *params, Double_t eps) const
01149 {
01150
01151
01152
01153
01154
01155
01156
01157
01158
01159
01160
01161
01162
01163
01164
01165
01166
01167
01168
01169
01170
01171
01172
01173
01174
01175
01176
01177
01178
01179
01180
01181 if (GetNdim() > 1) {
01182 Warning("Derivative","Function dimension is larger than one");
01183 }
01184
01185 ROOT::Math::RichardsonDerivator rd;
01186 double xmin, xmax;
01187 GetRange(xmin, xmax);
01188
01189 double h = eps* std::abs(xmax-xmin);
01190 if ( h <= 0 ) h = 0.001;
01191 double der = 0;
01192 if (params) {
01193 ROOT::Math::WrappedTF1 wtf(*( const_cast<TF1 *> (this) ));
01194 wtf.SetParameters(params);
01195 der = rd.Derivative1(wtf,x,h);
01196 }
01197 else {
01198
01199
01200 ROOT::Math::WrappedFunction<const TF1 & > wf( *this);
01201 der = rd.Derivative1(wf,x,h);
01202 }
01203
01204 gErrorTF1 = rd.Error();
01205 return der;
01206
01207 }
01208
01209
01210
01211 Double_t TF1::Derivative2(Double_t x, Double_t *params, Double_t eps) const
01212 {
01213
01214
01215
01216
01217
01218
01219
01220
01221
01222
01223
01224
01225
01226
01227
01228
01229
01230
01231
01232
01233
01234
01235
01236
01237
01238
01239
01240
01241
01242
01243
01244 if (GetNdim() > 1) {
01245 Warning("Derivative2","Function dimension is larger than one");
01246 }
01247
01248 ROOT::Math::RichardsonDerivator rd;
01249 double xmin, xmax;
01250 GetRange(xmin, xmax);
01251
01252 double h = eps* std::abs(xmax-xmin);
01253 if ( h <= 0 ) h = 0.001;
01254 double der = 0;
01255 if (params) {
01256 ROOT::Math::WrappedTF1 wtf(*( const_cast<TF1 *> (this) ));
01257 wtf.SetParameters(params);
01258 der = rd.Derivative2(wtf,x,h);
01259 }
01260 else {
01261
01262
01263 ROOT::Math::WrappedFunction<const TF1 & > wf( *this);
01264 der = rd.Derivative2(wf,x,h);
01265 }
01266
01267 gErrorTF1 = rd.Error();
01268
01269 return der;
01270 }
01271
01272
01273
01274 Double_t TF1::Derivative3(Double_t x, Double_t *params, Double_t eps) const
01275 {
01276
01277
01278
01279
01280
01281
01282
01283
01284
01285
01286
01287
01288
01289
01290
01291
01292
01293
01294
01295
01296
01297
01298
01299
01300
01301
01302
01303
01304
01305
01306
01307 if (GetNdim() > 1) {
01308 Warning("Derivative3","Function dimension is larger than one");
01309 }
01310
01311 ROOT::Math::RichardsonDerivator rd;
01312 double xmin, xmax;
01313 GetRange(xmin, xmax);
01314
01315 double h = eps* std::abs(xmax-xmin);
01316 if ( h <= 0 ) h = 0.001;
01317 double der = 0;
01318 if (params) {
01319 ROOT::Math::WrappedTF1 wtf(*( const_cast<TF1 *> (this) ));
01320 wtf.SetParameters(params);
01321 der = rd.Derivative3(wtf,x,h);
01322 }
01323 else {
01324
01325
01326 ROOT::Math::WrappedFunction<const TF1 & > wf( *this);
01327 der = rd.Derivative3(wf,x,h);
01328 }
01329
01330 gErrorTF1 = rd.Error();
01331 return der;
01332
01333 }
01334
01335
01336
01337 Double_t TF1::DerivativeError()
01338 {
01339
01340
01341
01342 return gErrorTF1;
01343 }
01344
01345
01346
01347 Int_t TF1::DistancetoPrimitive(Int_t px, Int_t py)
01348 {
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358 if (!fHistogram) return 9999;
01359 Int_t distance = 9999;
01360 if (px >= 0) {
01361 distance = fHistogram->DistancetoPrimitive(px,py);
01362 if (distance <= 1) return distance;
01363 } else {
01364 px = -px;
01365 }
01366
01367 Double_t xx[1];
01368 Double_t x = gPad->AbsPixeltoX(px);
01369 xx[0] = gPad->PadtoX(x);
01370 if (xx[0] < fXmin || xx[0] > fXmax) return distance;
01371 Double_t fval = Eval(xx[0]);
01372 Double_t y = gPad->YtoPad(fval);
01373 Int_t pybin = gPad->YtoAbsPixel(y);
01374 return TMath::Abs(py - pybin);
01375 }
01376
01377
01378
01379 void TF1::Draw(Option_t *option)
01380 {
01381
01382
01383
01384
01385
01386
01387
01388
01389
01390
01391
01392
01393
01394
01395 TString opt = option;
01396 opt.ToLower();
01397 if (gPad && !opt.Contains("same")) gPad->Clear();
01398
01399 AppendPad(option);
01400 }
01401
01402
01403
01404 TF1 *TF1::DrawCopy(Option_t *option) const
01405 {
01406
01407
01408
01409
01410
01411
01412
01413
01414
01415
01416
01417
01418
01419
01420 TF1 *newf1 = (TF1*)this->IsA()->New();
01421 Copy(*newf1);
01422 newf1->AppendPad(option);
01423 newf1->SetBit(kCanDelete);
01424 return newf1;
01425 }
01426
01427
01428
01429 TObject *TF1::DrawDerivative(Option_t *option)
01430 {
01431
01432
01433
01434
01435
01436
01437
01438
01439
01440
01441 TVirtualPad *pad = gROOT->GetSelectedPad();
01442 TVirtualPad *padsav = gPad;
01443 if (pad) pad->cd();
01444
01445 TGraph *gr = new TGraph(this,"d");
01446 gr->Draw(option);
01447 if (padsav) padsav->cd();
01448 return gr;
01449 }
01450
01451
01452
01453 TObject *TF1::DrawIntegral(Option_t *option)
01454 {
01455
01456
01457
01458
01459
01460
01461
01462
01463
01464
01465 TVirtualPad *pad = gROOT->GetSelectedPad();
01466 TVirtualPad *padsav = gPad;
01467 if (pad) pad->cd();
01468
01469 TGraph *gr = new TGraph(this,"i");
01470 gr->Draw(option);
01471 if (padsav) padsav->cd();
01472 return gr;
01473 }
01474
01475
01476
01477 void TF1::DrawF1(const char *formula, Double_t xmin, Double_t xmax, Option_t *option)
01478 {
01479
01480
01481 if (Compile(formula)) return;
01482
01483 SetRange(xmin, xmax);
01484
01485 Draw(option);
01486 }
01487
01488
01489
01490 Double_t TF1::Eval(Double_t x, Double_t y, Double_t z, Double_t t) const
01491 {
01492
01493
01494
01495
01496
01497
01498
01499
01500
01501 Double_t xx[4];
01502 xx[0] = x;
01503 xx[1] = y;
01504 xx[2] = z;
01505 xx[3] = t;
01506
01507 ((TF1*)this)->InitArgs(xx,fParams);
01508
01509 return ((TF1*)this)->EvalPar(xx,fParams);
01510 }
01511
01512
01513
01514 Double_t TF1::EvalPar(const Double_t *x, const Double_t *params)
01515 {
01516
01517
01518
01519
01520
01521
01522
01523
01524
01525
01526
01527
01528
01529
01530
01531
01532
01533 fgCurrent = this;
01534
01535 if (fType == 0) return TFormula::EvalPar(x,params);
01536 Double_t result = 0;
01537 if (fType == 1) {
01538
01539
01540
01541 if (!fFunctor.Empty()) {
01542 if (params) result = fFunctor((Double_t*)x,(Double_t*)params);
01543 else result = fFunctor((Double_t*)x,fParams);
01544
01545 }else result = GetSave(x);
01546 return result;
01547 }
01548 if (fType == 2) {
01549 if (fMethodCall) fMethodCall->Execute(result);
01550 else result = GetSave(x);
01551 return result;
01552 }
01553 if (fType == 3) {
01554
01555 if (fMethodCall) fMethodCall->Execute(fCintFunc,result);
01556 else result = GetSave(x);
01557 return result;
01558 }
01559 return result;
01560 }
01561
01562
01563
01564 void TF1::ExecuteEvent(Int_t event, Int_t px, Int_t py)
01565 {
01566
01567
01568
01569
01570 if (fHistogram) fHistogram->ExecuteEvent(event,px,py);
01571
01572 if (!gPad->GetView()) {
01573 if (event == kMouseMotion) gPad->SetCursor(kHand);
01574 }
01575 }
01576
01577
01578
01579 void TF1::FixParameter(Int_t ipar, Double_t value)
01580 {
01581
01582
01583
01584 if (ipar < 0 || ipar > fNpar-1) return;
01585 SetParameter(ipar,value);
01586 if (value != 0) SetParLimits(ipar,value,value);
01587 else SetParLimits(ipar,1,1);
01588 }
01589
01590
01591
01592 TF1 *TF1::GetCurrent()
01593 {
01594
01595
01596 return fgCurrent;
01597 }
01598
01599
01600
01601 TH1 *TF1::GetHistogram() const
01602 {
01603
01604
01605 if (fHistogram) return fHistogram;
01606
01607
01608 ((TF1*)this)->Paint();
01609 return fHistogram;
01610 }
01611
01612
01613
01614 Double_t TF1::GetMaximum(Double_t xmin, Double_t xmax, Double_t epsilon, Int_t maxiter,Bool_t logx) const
01615 {
01616
01617
01618
01619
01620
01621
01622
01623
01624
01625
01626
01627
01628
01629
01630
01631
01632 if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;}
01633
01634 if (!logx && gPad != 0) logx = gPad->GetLogx();
01635
01636 ROOT::Math::BrentMinimizer1D bm;
01637 GInverseFunc g(this);
01638 ROOT::Math::WrappedFunction<GInverseFunc> wf1(g);
01639 bm.SetFunction( wf1, xmin, xmax );
01640 bm.SetNpx(fNpx);
01641 bm.SetLogScan(logx);
01642 bm.Minimize(maxiter, epsilon, epsilon );
01643 Double_t x;
01644 x = - bm.FValMinimum();
01645
01646 return x;
01647 }
01648
01649
01650
01651 Double_t TF1::GetMaximumX(Double_t xmin, Double_t xmax, Double_t epsilon, Int_t maxiter,Bool_t logx) const
01652 {
01653
01654
01655
01656
01657
01658
01659
01660
01661
01662
01663
01664
01665
01666
01667
01668
01669 if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;}
01670
01671 if (!logx && gPad != 0) logx = gPad->GetLogx();
01672
01673 ROOT::Math::BrentMinimizer1D bm;
01674 GInverseFunc g(this);
01675 ROOT::Math::WrappedFunction<GInverseFunc> wf1(g);
01676 bm.SetFunction( wf1, xmin, xmax );
01677 bm.SetNpx(fNpx);
01678 bm.SetLogScan(logx);
01679 bm.Minimize(maxiter, epsilon, epsilon );
01680 Double_t x;
01681 x = bm.XMinimum();
01682
01683 return x;
01684 }
01685
01686
01687
01688 Double_t TF1::GetMinimum(Double_t xmin, Double_t xmax, Double_t epsilon, Int_t maxiter, Bool_t logx) const
01689 {
01690
01691
01692
01693
01694
01695
01696
01697
01698
01699
01700
01701
01702
01703
01704
01705
01706 if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;}
01707
01708 if (!logx && gPad != 0) logx = gPad->GetLogx();
01709
01710 ROOT::Math::BrentMinimizer1D bm;
01711 ROOT::Math::WrappedFunction<const TF1&> wf1(*this);
01712 bm.SetFunction( wf1, xmin, xmax );
01713 bm.SetNpx(fNpx);
01714 bm.SetLogScan(logx);
01715 bm.Minimize(maxiter, epsilon, epsilon );
01716 Double_t x;
01717 x = bm.FValMinimum();
01718
01719 return x;
01720 }
01721
01722
01723
01724 Double_t TF1::GetMinimumX(Double_t xmin, Double_t xmax, Double_t epsilon, Int_t maxiter, Bool_t logx) const
01725 {
01726
01727
01728
01729
01730
01731
01732
01733
01734
01735
01736
01737
01738
01739
01740
01741
01742
01743 if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;}
01744
01745 ROOT::Math::BrentMinimizer1D bm;
01746 ROOT::Math::WrappedFunction<const TF1&> wf1(*this);
01747 bm.SetFunction( wf1, xmin, xmax );
01748 bm.SetNpx(fNpx);
01749 bm.SetLogScan(logx);
01750 bm.Minimize(maxiter, epsilon, epsilon );
01751 Double_t x;
01752 x = bm.XMinimum();
01753
01754 return x;
01755 }
01756
01757
01758
01759 Double_t TF1::GetX(Double_t fy, Double_t xmin, Double_t xmax, Double_t epsilon, Int_t maxiter, Bool_t logx) const
01760 {
01761
01762
01763
01764
01765
01766
01767
01768
01769
01770
01771
01772
01773
01774
01775
01776
01777
01778
01779 if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;}
01780
01781 if (!logx && gPad != 0) logx = gPad->GetLogx();
01782
01783 GFunc g(this, fy);
01784 ROOT::Math::WrappedFunction<GFunc> wf1(g);
01785 ROOT::Math::BrentRootFinder brf;
01786 brf.SetFunction(wf1,xmin,xmax);
01787 brf.SetNpx(fNpx);
01788 brf.SetLogScan(logx);
01789 brf.Solve(maxiter, epsilon, epsilon);
01790 return brf.Root();
01791
01792 }
01793
01794
01795 Int_t TF1::GetNDF() const
01796 {
01797
01798
01799
01800
01801
01802 if (fNDF == 0 && (fNpfits > fNpar) ) return fNpfits-fNpar;
01803 return fNDF;
01804 }
01805
01806
01807
01808 Int_t TF1::GetNumberFreeParameters() const
01809 {
01810
01811
01812 Int_t nfree = fNpar;
01813 Double_t al,bl;
01814 for (Int_t i=0;i<fNpar;i++) {
01815 ((TF1*)this)->GetParLimits(i,al,bl);
01816 if (al*bl != 0 && al >= bl) nfree--;
01817 }
01818 return nfree;
01819 }
01820
01821
01822
01823 char *TF1::GetObjectInfo(Int_t px, Int_t ) const
01824 {
01825
01826
01827
01828
01829 static char info[64];
01830 Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
01831 snprintf(info,64,"(x=%g, f=%g)",x,((TF1*)this)->Eval(x));
01832 return info;
01833 }
01834
01835
01836
01837 Double_t TF1::GetParError(Int_t ipar) const
01838 {
01839
01840
01841 if (ipar < 0 || ipar > fNpar-1) return 0;
01842 return fParErrors[ipar];
01843 }
01844
01845
01846
01847 void TF1::GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
01848 {
01849
01850
01851 parmin = 0;
01852 parmax = 0;
01853 if (ipar < 0 || ipar > fNpar-1) return;
01854 if (fParMin) parmin = fParMin[ipar];
01855 if (fParMax) parmax = fParMax[ipar];
01856 }
01857
01858
01859
01860 Double_t TF1::GetProb() const
01861 {
01862
01863
01864 if (fNDF <= 0) return 0;
01865 return TMath::Prob(fChisquare,fNDF);
01866 }
01867
01868
01869
01870 Int_t TF1::GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum)
01871 {
01872
01873
01874
01875
01876
01877
01878
01879
01880
01881
01882
01883
01884
01885
01886
01887
01888
01889
01890
01891
01892
01893
01894
01895
01896
01897
01898
01899
01900
01901
01902
01903
01904 const Int_t npx = TMath::Max(fNpx, 2*nprobSum);
01905 const Double_t xMin = GetXmin();
01906 const Double_t xMax = GetXmax();
01907 const Double_t dx = (xMax-xMin)/npx;
01908
01909 TArrayD integral(npx+1);
01910 TArrayD alpha(npx);
01911 TArrayD beta(npx);
01912 TArrayD gamma(npx);
01913
01914 integral[0] = 0;
01915 Int_t intNegative = 0;
01916 Int_t i;
01917 for (i = 0; i < npx; i++) {
01918 const Double_t *params = 0;
01919 Double_t integ = Integral(Double_t(xMin+i*dx),Double_t(xMin+i*dx+dx),params);
01920 if (integ < 0) {intNegative++; integ = -integ;}
01921 integral[i+1] = integral[i] + integ;
01922 }
01923
01924 if (intNegative > 0)
01925 Warning("GetQuantiles","function:%s has %d negative values: abs assumed",
01926 GetName(),intNegative);
01927 if (integral[npx] == 0) {
01928 Error("GetQuantiles","Integral of function is zero");
01929 return 0;
01930 }
01931
01932 const Double_t total = integral[npx];
01933 for (i = 1; i <= npx; i++) integral[i] /= total;
01934
01935
01936
01937 for (i = 0; i < npx; i++) {
01938 const Double_t x0 = xMin+dx*i;
01939 const Double_t r2 = integral[i+1]-integral[i];
01940 const Double_t r1 = Integral(x0,x0+0.5*dx)/total;
01941 gamma[i] = (2*r2-4*r1)/(dx*dx);
01942 beta[i] = r2/dx-gamma[i]*dx;
01943 alpha[i] = x0;
01944 gamma[i] *= 2;
01945 }
01946
01947
01948
01949 for (i = 0; i < nprobSum; i++) {
01950 const Double_t r = probSum[i];
01951 Int_t bin = TMath::Max(TMath::BinarySearch(npx+1,integral.GetArray(),r),(Long64_t)0);
01952
01953 while (bin < npx-1 && TMath::AreEqualRel(integral[bin+1], r, 1E-12) ) {
01954 if (TMath::AreEqualRel(integral[bin+2], r, 1E-12) ) bin++;
01955 else break;
01956 }
01957
01958 const Double_t rr = r-integral[bin];
01959 if (rr != 0.0) {
01960 Double_t xx = 0.0;
01961 const Double_t fac = -2.*gamma[bin]*rr/beta[bin]/beta[bin];
01962 if (fac != 0 && fac <= 1)
01963 xx = (-beta[bin]+TMath::Sqrt(beta[bin]*beta[bin]+2*gamma[bin]*rr))/gamma[bin];
01964 else if (beta[bin] != 0.)
01965 xx = rr/beta[bin];
01966 q[i] = alpha[bin]+xx;
01967 } else {
01968 q[i] = alpha[bin];
01969 if (integral[bin+1] == r) q[i] += dx;
01970 }
01971 }
01972
01973 return nprobSum;
01974 }
01975
01976
01977
01978 Double_t TF1::GetRandom()
01979 {
01980
01981
01982
01983
01984
01985
01986
01987
01988
01989
01990
01991
01992
01993
01994
01995
01996
01997 if (fIntegral == 0) {
01998 fIntegral = new Double_t[fNpx+1];
01999 fAlpha = new Double_t[fNpx+1];
02000 fBeta = new Double_t[fNpx];
02001 fGamma = new Double_t[fNpx];
02002 fIntegral[0] = 0;
02003 fAlpha[fNpx] = 0;
02004 Double_t integ;
02005 Int_t intNegative = 0;
02006 Int_t i;
02007 Bool_t logbin = kFALSE;
02008 Double_t dx;
02009 Double_t xmin = fXmin;
02010 Double_t xmax = fXmax;
02011 if (xmin > 0 && xmax/xmin> fNpx) {
02012 logbin = kTRUE;
02013 fAlpha[fNpx] = 1;
02014 xmin = TMath::Log10(fXmin);
02015 xmax = TMath::Log10(fXmax);
02016 }
02017 dx = (xmax-xmin)/fNpx;
02018
02019 Double_t *xx = new Double_t[fNpx+1];
02020 for (i=0;i<fNpx;i++) {
02021 xx[i] = xmin +i*dx;
02022 }
02023 xx[fNpx] = xmax;
02024 for (i=0;i<fNpx;i++) {
02025 if (logbin) {
02026 integ = Integral(TMath::Power(10,xx[i]), TMath::Power(10,xx[i+1]));
02027 } else {
02028 integ = Integral(xx[i],xx[i+1]);
02029 }
02030 if (integ < 0) {intNegative++; integ = -integ;}
02031 fIntegral[i+1] = fIntegral[i] + integ;
02032 }
02033 if (intNegative > 0) {
02034 Warning("GetRandom","function:%s has %d negative values: abs assumed",GetName(),intNegative);
02035 }
02036 if (fIntegral[fNpx] == 0) {
02037 delete [] xx;
02038 Error("GetRandom","Integral of function is zero");
02039 return 0;
02040 }
02041 Double_t total = fIntegral[fNpx];
02042 for (i=1;i<=fNpx;i++) {
02043 fIntegral[i] /= total;
02044 }
02045
02046
02047
02048 Double_t x0,r1,r2,r3;
02049 for (i=0;i<fNpx;i++) {
02050 x0 = xx[i];
02051 r2 = fIntegral[i+1] - fIntegral[i];
02052 if (logbin) r1 = Integral(TMath::Power(10,x0),TMath::Power(10,x0+0.5*dx))/total;
02053 else r1 = Integral(x0,x0+0.5*dx)/total;
02054 r3 = 2*r2 - 4*r1;
02055 if (TMath::Abs(r3) > 1e-8) fGamma[i] = r3/(dx*dx);
02056 else fGamma[i] = 0;
02057 fBeta[i] = r2/dx - fGamma[i]*dx;
02058 fAlpha[i] = x0;
02059 fGamma[i] *= 2;
02060 }
02061 delete [] xx;
02062 }
02063
02064
02065 Double_t r = gRandom->Rndm();
02066 Int_t bin = TMath::BinarySearch(fNpx,fIntegral,r);
02067 Double_t rr = r - fIntegral[bin];
02068
02069 Double_t yy;
02070 if(fGamma[bin] != 0)
02071 yy = (-fBeta[bin] + TMath::Sqrt(fBeta[bin]*fBeta[bin]+2*fGamma[bin]*rr))/fGamma[bin];
02072 else
02073 yy = rr/fBeta[bin];
02074 Double_t x = fAlpha[bin] + yy;
02075 if (fAlpha[fNpx] > 0) return TMath::Power(10,x);
02076 return x;
02077 }
02078
02079
02080
02081 Double_t TF1::GetRandom(Double_t xmin, Double_t xmax)
02082 {
02083
02084
02085
02086
02087
02088
02089
02090
02091
02092
02093
02094
02095
02096
02097
02098
02099
02100
02101
02102
02103
02104 if (fIntegral == 0) {
02105 Double_t dx = (fXmax-fXmin)/fNpx;
02106 fIntegral = new Double_t[fNpx+1];
02107 fAlpha = new Double_t[fNpx];
02108 fBeta = new Double_t[fNpx];
02109 fGamma = new Double_t[fNpx];
02110 fIntegral[0] = 0;
02111 Double_t integ;
02112 Int_t intNegative = 0;
02113 Int_t i;
02114 for (i=0;i<fNpx;i++) {
02115 integ = Integral(Double_t(fXmin+i*dx), Double_t(fXmin+i*dx+dx));
02116 if (integ < 0) {intNegative++; integ = -integ;}
02117 fIntegral[i+1] = fIntegral[i] + integ;
02118 }
02119 if (intNegative > 0) {
02120 Warning("GetRandom","function:%s has %d negative values: abs assumed",GetName(),intNegative);
02121 }
02122 if (fIntegral[fNpx] == 0) {
02123 Error("GetRandom","Integral of function is zero");
02124 return 0;
02125 }
02126 Double_t total = fIntegral[fNpx];
02127 for (i=1;i<=fNpx;i++) {
02128 fIntegral[i] /= total;
02129 }
02130
02131
02132
02133 Double_t x0,r1,r2,r3;
02134 for (i=0;i<fNpx;i++) {
02135 x0 = fXmin+i*dx;
02136 r2 = fIntegral[i+1] - fIntegral[i];
02137 r1 = Integral(x0,x0+0.5*dx)/total;
02138 r3 = 2*r2 - 4*r1;
02139 if (TMath::Abs(r3) > 1e-8) fGamma[i] = r3/(dx*dx);
02140 else fGamma[i] = 0;
02141 fBeta[i] = r2/dx - fGamma[i]*dx;
02142 fAlpha[i] = x0;
02143 fGamma[i] *= 2;
02144 }
02145 }
02146
02147
02148 Double_t dx = (fXmax-fXmin)/fNpx;
02149 Int_t nbinmin = (Int_t)((xmin-fXmin)/dx);
02150 Int_t nbinmax = (Int_t)((xmax-fXmin)/dx)+2;
02151 if(nbinmax>fNpx) nbinmax=fNpx;
02152
02153 Double_t pmin=fIntegral[nbinmin];
02154 Double_t pmax=fIntegral[nbinmax];
02155
02156 Double_t r,x,xx,rr;
02157 do {
02158 r = gRandom->Uniform(pmin,pmax);
02159
02160 Int_t bin = TMath::BinarySearch(fNpx,fIntegral,r);
02161 rr = r - fIntegral[bin];
02162
02163 if(fGamma[bin] != 0)
02164 xx = (-fBeta[bin] + TMath::Sqrt(fBeta[bin]*fBeta[bin]+2*fGamma[bin]*rr))/fGamma[bin];
02165 else
02166 xx = rr/fBeta[bin];
02167 x = fAlpha[bin] + xx;
02168 } while(x<xmin || x>xmax);
02169 return x;
02170 }
02171
02172
02173
02174 void TF1::GetRange(Double_t &xmin, Double_t &xmax) const
02175 {
02176
02177
02178 xmin = fXmin;
02179 xmax = fXmax;
02180 }
02181
02182
02183
02184 void TF1::GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
02185 {
02186
02187
02188 xmin = fXmin;
02189 xmax = fXmax;
02190 ymin = 0;
02191 ymax = 0;
02192 }
02193
02194
02195
02196 void TF1::GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const
02197 {
02198
02199
02200 xmin = fXmin;
02201 xmax = fXmax;
02202 ymin = 0;
02203 ymax = 0;
02204 zmin = 0;
02205 zmax = 0;
02206 }
02207
02208
02209
02210 Double_t TF1::GetSave(const Double_t *xx)
02211 {
02212
02213
02214 if (fNsave <= 0) return 0;
02215 if (fSave == 0) return 0;
02216 Double_t x = Double_t(xx[0]);
02217 Double_t y,dx,xmin,xmax,xlow,xup,ylow,yup;
02218 if (fParent && fParent->InheritsFrom(TH1::Class())) {
02219
02220
02221 xmin = fSave[fNsave-3];
02222 xmax = fSave[fNsave-2];
02223 if (fSave[fNsave-1] == xmax) {
02224 TH1 *h = (TH1*)fParent;
02225 TAxis *xaxis = h->GetXaxis();
02226 Int_t bin1 = xaxis->FindBin(xmin);
02227 Int_t binup = xaxis->FindBin(xmax);
02228 Int_t bin = xaxis->FindBin(x);
02229 if (bin < binup) {
02230 xlow = xaxis->GetBinCenter(bin);
02231 xup = xaxis->GetBinCenter(bin+1);
02232 ylow = fSave[bin-bin1];
02233 yup = fSave[bin-bin1+1];
02234 } else {
02235 xlow = xaxis->GetBinCenter(bin-1);
02236 xup = xaxis->GetBinCenter(bin);
02237 ylow = fSave[bin-bin1-1];
02238 yup = fSave[bin-bin1];
02239 }
02240 dx = xup-xlow;
02241 y = ((xup*ylow-xlow*yup) + x*(yup-ylow))/dx;
02242 return y;
02243 }
02244 }
02245 Int_t np = fNsave - 3;
02246 xmin = Double_t(fSave[np+1]);
02247 xmax = Double_t(fSave[np+2]);
02248 dx = (xmax-xmin)/np;
02249 if (x < xmin || x > xmax) return 0;
02250
02251 if (TMath::IsNaN(x) ) return x;
02252 if (dx <= 0) return 0;
02253
02254 Int_t bin = Int_t((x-xmin)/dx);
02255 xlow = xmin + bin*dx;
02256 xup = xlow + dx;
02257 ylow = fSave[bin];
02258 yup = fSave[bin+1];
02259 y = ((xup*ylow-xlow*yup) + x*(yup-ylow))/dx;
02260 return y;
02261 }
02262
02263
02264
02265 TAxis *TF1::GetXaxis() const
02266 {
02267
02268
02269 TH1 *h = GetHistogram();
02270 if (!h) return 0;
02271 return h->GetXaxis();
02272 }
02273
02274
02275
02276 TAxis *TF1::GetYaxis() const
02277 {
02278
02279
02280 TH1 *h = GetHistogram();
02281 if (!h) return 0;
02282 return h->GetYaxis();
02283 }
02284
02285
02286
02287 TAxis *TF1::GetZaxis() const
02288 {
02289
02290
02291 TH1 *h = GetHistogram();
02292 if (!h) return 0;
02293 return h->GetZaxis();
02294 }
02295
02296
02297
02298
02299 Double_t TF1::GradientPar(Int_t ipar, const Double_t *x, Double_t eps)
02300 {
02301
02302
02303
02304
02305
02306
02307
02308
02309
02310
02311
02312
02313 if (fNpar == 0) return 0;
02314
02315 if(eps< 1e-10 || eps > 1) {
02316 Warning("Derivative","parameter esp=%g out of allowed range[1e-10,1], reset to 0.01",eps);
02317 eps = 0.01;
02318 }
02319 Double_t h;
02320 TF1 *func = (TF1*)this;
02321
02322 Double_t par0 = fParams[ipar];
02323
02324
02325 func->InitArgs(x, fParams);
02326
02327 Double_t al, bl;
02328 Double_t f1, f2, g1, g2, h2, d0, d2;
02329
02330 ((TF1*)this)->GetParLimits(ipar,al,bl);
02331 if (al*bl != 0 && al >= bl) {
02332
02333 return 0;
02334 }
02335
02336
02337 if (func->GetParError(ipar)!=0)
02338 h = eps*func->GetParError(ipar);
02339 else
02340 h=eps;
02341
02342
02343
02344 fParams[ipar] = par0 + h; f1 = func->EvalPar(x,fParams);
02345 fParams[ipar] = par0 - h; f2 = func->EvalPar(x,fParams);
02346 fParams[ipar] = par0 + h/2; g1 = func->EvalPar(x,fParams);
02347 fParams[ipar] = par0 - h/2; g2 = func->EvalPar(x,fParams);
02348
02349
02350 h2 = 1/(2.*h);
02351 d0 = f1 - f2;
02352 d2 = 2*(g1 - g2);
02353
02354 Double_t grad = h2*(4*d2 - d0)/3.;
02355
02356
02357 fParams[ipar] = par0;
02358
02359 return grad;
02360 }
02361
02362
02363 void TF1::GradientPar(const Double_t *x, Double_t *grad, Double_t eps)
02364 {
02365
02366
02367
02368
02369
02370
02371
02372
02373
02374
02375
02376
02377 if(eps< 1e-10 || eps > 1) {
02378 Warning("Derivative","parameter esp=%g out of allowed range[1e-10,1], reset to 0.01",eps);
02379 eps = 0.01;
02380 }
02381
02382 for (Int_t ipar=0; ipar<fNpar; ipar++){
02383 grad[ipar] = GradientPar(ipar,x,eps);
02384 }
02385 }
02386
02387
02388 void TF1::InitArgs(const Double_t *x, const Double_t *params)
02389 {
02390
02391
02392 if (fMethodCall) {
02393 Long_t args[2];
02394 args[0] = (Long_t)x;
02395 if (params) args[1] = (Long_t)params;
02396 else args[1] = (Long_t)fParams;
02397 fMethodCall->SetParamPtrs(args);
02398 }
02399 }
02400
02401
02402
02403 void TF1::InitStandardFunctions()
02404 {
02405
02406
02407 TF1 *f1;
02408 if (!gROOT->GetListOfFunctions()->FindObject("gaus")) {
02409 f1 = new TF1("gaus","gaus",-1,1); f1->SetParameters(1,0,1);
02410 f1 = new TF1("gausn","gausn",-1,1); f1->SetParameters(1,0,1);
02411 f1 = new TF1("landau","landau",-1,1); f1->SetParameters(1,0,1);
02412 f1 = new TF1("landaun","landaun",-1,1); f1->SetParameters(1,0,1);
02413 f1 = new TF1("expo","expo",-1,1); f1->SetParameters(1,1);
02414 for (Int_t i=0;i<10;i++) {
02415 f1 = new TF1(Form("pol%d",i),Form("pol%d",i),-1,1);
02416 f1->SetParameters(1,1,1,1,1,1,1,1,1,1);
02417 }
02418 }
02419 }
02420
02421
02422
02423 Double_t TF1::Integral(Double_t a, Double_t b, const Double_t *params, Double_t epsilon)
02424 {
02425
02426
02427
02428
02429
02430
02431
02432
02433
02434
02435
02436
02437
02438
02439
02440
02441
02442
02443
02444
02445
02446
02447
02448
02449
02450
02451
02452
02453
02454
02455
02456
02457
02458
02459
02460
02461
02462
02463
02464
02465
02466
02467
02468
02469
02470
02471
02472
02473
02474
02475
02476
02477
02478
02479
02480
02481
02482
02483
02484
02485
02486
02487
02488
02489
02490
02491
02492
02493
02494
02495
02496
02497
02498
02499
02500
02501
02502
02503
02504
02505
02506
02507
02508
02509
02510
02511
02512
02513
02514
02515
02516
02517
02518
02519
02520
02521
02522
02523
02524
02525
02526
02527
02528
02529
02530
02531
02532
02533
02534
02535
02536
02537
02538
02539
02540
02541
02542 TF1_EvalWrapper wf1( this, params, fgAbsValue );
02543
02544 ROOT::Math::GaussIntegrator giod;
02545
02546 giod.SetFunction(wf1);
02547 giod.SetRelTolerance(epsilon);
02548
02549
02550 return giod.Integral(a, b);
02551 }
02552
02553
02554
02555 Double_t TF1::Integral(Double_t, Double_t, Double_t, Double_t, Double_t)
02556 {
02557
02558
02559 Error("Integral","Must be called with a TF2 only");
02560 return 0;
02561 }
02562
02563
02564
02565 Double_t TF1::Integral(Double_t, Double_t, Double_t, Double_t, Double_t, Double_t, Double_t)
02566 {
02567
02568
02569 Error("Integral","Must be called with a TF3 only");
02570 return 0;
02571 }
02572
02573
02574 Double_t TF1::IntegralError(Double_t a, Double_t b, const Double_t * params, const Double_t * covmat, Double_t epsilon )
02575 {
02576
02577
02578
02579
02580
02581
02582
02583
02584
02585
02586
02587
02588
02589
02590
02591
02592
02593
02594
02595
02596
02597
02598
02599 Double_t x1[1];
02600 Double_t x2[1];
02601 x1[0] = a, x2[0] = b;
02602 return ROOT::TF1Helper::IntegralError(this,1,x1,x2,params,covmat,epsilon);
02603 }
02604
02605
02606 Double_t TF1::IntegralError(Int_t n, const Double_t * a, const Double_t * b, const Double_t * params, const Double_t * covmat, Double_t epsilon )
02607 {
02608
02609
02610
02611
02612
02613
02614
02615
02616
02617
02618
02619
02620
02621
02622
02623
02624
02625
02626
02627
02628
02629
02630
02631
02632
02633 return ROOT::TF1Helper::IntegralError(this,n,a,b,params,covmat,epsilon);
02634 }
02635
02636 #ifdef INTHEFUTURE
02637
02638 Double_t TF1::IntegralFast(const TGraph *g, Double_t a, Double_t b, Double_t *params)
02639 {
02640
02641
02642 if (!g) return 0;
02643 return IntegralFast(g->GetN(), g->GetX(), g->GetY(), a, b, params);
02644 }
02645 #endif
02646
02647
02648
02649 Double_t TF1::IntegralFast(Int_t num, Double_t * , Double_t * , Double_t a, Double_t b, Double_t *params, Double_t epsilon)
02650 {
02651
02652
02653
02654
02655 ROOT::Math::WrappedTF1 wf1(*this);
02656 if ( params )
02657 wf1.SetParameters( params );
02658 ROOT::Math::GaussLegendreIntegrator gli(num,epsilon);
02659 gli.SetFunction( wf1 );
02660 return gli.Integral(a, b);
02661
02662 }
02663
02664
02665
02666 Double_t TF1::IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Double_t eps, Double_t &relerr)
02667 {
02668
02669
02670
02671 Int_t nfnevl,ifail;
02672 Int_t minpts = 2+2*n*(n+1)+1;
02673 Int_t maxpts = 1000;
02674 Double_t result = IntegralMultiple(n,a,b,minpts, maxpts,eps,relerr,nfnevl,ifail);
02675 if (ifail > 0) {
02676 Warning("IntegralMultiple","failed code=%d, ",ifail);
02677 }
02678 return result;
02679 }
02680
02681
02682
02683 Double_t TF1::IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Int_t minpts, Int_t maxpts, Double_t eps, Double_t &relerr,Int_t &nfnevl, Int_t &ifail)
02684 {
02685
02686
02687
02688
02689
02690
02691
02692
02693
02694
02695
02696
02697
02698
02699
02700
02701
02702
02703
02704
02705
02706
02707
02708
02709
02710
02711
02712
02713
02714
02715
02716
02717
02718
02719
02720
02721
02722
02723
02724
02725
02726
02727
02728
02729
02730
02731
02732
02733
02734
02735
02736
02737
02738
02739
02740
02741
02742
02743
02744
02745
02746
02747
02748 ROOT::Math::WrappedMultiFunction<TF1&> wf1(*this, n);
02749
02750 ROOT::Math::AdaptiveIntegratorMultiDim aimd(wf1, eps, eps, maxpts);
02751 aimd.SetMinPts(minpts);
02752 double result = aimd.Integral(a,b);
02753 relerr = aimd.RelError();
02754 nfnevl = aimd.NEval();
02755 ifail = 0;
02756
02757 return result;
02758 }
02759
02760
02761
02762 Bool_t TF1::IsInside(const Double_t *x) const
02763 {
02764
02765
02766 if (x[0] < fXmin || x[0] > fXmax) return kFALSE;
02767 return kTRUE;
02768 }
02769
02770
02771
02772 void TF1::Paint(Option_t *option)
02773 {
02774
02775
02776 Int_t i;
02777 Double_t xv[1];
02778
02779 fgCurrent = this;
02780
02781 TString opt = option;
02782 opt.ToLower();
02783 Bool_t optSAME = kFALSE;
02784 if (opt.Contains("same")) optSAME = kTRUE;
02785
02786 Double_t xmin=fXmin, xmax=fXmax, pmin=fXmin, pmax=fXmax;
02787 if (gPad) {
02788 pmin = gPad->PadtoX(gPad->GetUxmin());
02789 pmax = gPad->PadtoX(gPad->GetUxmax());
02790 }
02791 if (optSAME) {
02792 if (xmax < pmin) return;
02793 if (xmin > pmax) return;
02794 if (xmin < pmin) xmin = pmin;
02795 if (xmax > pmax) xmax = pmax;
02796 }
02797
02798
02799
02800 TString xtitle = "";
02801 TString ytitle = "";
02802 char *semicol = (char*)strstr(GetTitle(),";");
02803 if (semicol) {
02804 Int_t nxt = strlen(semicol);
02805 char *ctemp = new char[nxt];
02806 strlcpy(ctemp,semicol+1,nxt);
02807 semicol = (char*)strstr(ctemp,";");
02808 if (semicol) {
02809 *semicol = 0;
02810 ytitle = semicol+1;
02811 }
02812 xtitle = ctemp;
02813 delete [] ctemp;
02814 }
02815 if (fHistogram) {
02816 xtitle = fHistogram->GetXaxis()->GetTitle();
02817 ytitle = fHistogram->GetYaxis()->GetTitle();
02818 if (!gPad->GetLogx() && fHistogram->TestBit(TH1::kLogX)) { delete fHistogram; fHistogram = 0;}
02819 if ( gPad->GetLogx() && !fHistogram->TestBit(TH1::kLogX)) { delete fHistogram; fHistogram = 0;}
02820 }
02821
02822 if (fHistogram) {
02823 fHistogram->GetXaxis()->SetLimits(xmin,xmax);
02824 } else {
02825
02826
02827 if (xmin > 0 && gPad && gPad->GetLogx()) {
02828 Double_t *xbins = new Double_t[fNpx+1];
02829 Double_t xlogmin = TMath::Log10(xmin);
02830 Double_t xlogmax = TMath::Log10(xmax);
02831 Double_t dlogx = (xlogmax-xlogmin)/((Double_t)fNpx);
02832 for (i=0;i<=fNpx;i++) {
02833 xbins[i] = gPad->PadtoX(xlogmin+ i*dlogx);
02834 }
02835 fHistogram = new TH1D("Func",GetTitle(),fNpx,xbins);
02836 fHistogram->SetBit(TH1::kLogX);
02837 delete [] xbins;
02838 } else {
02839 fHistogram = new TH1D("Func",GetTitle(),fNpx,xmin,xmax);
02840 }
02841 if (!fHistogram) return;
02842 if (fMinimum != -1111) fHistogram->SetMinimum(fMinimum);
02843 if (fMaximum != -1111) fHistogram->SetMaximum(fMaximum);
02844 fHistogram->SetDirectory(0);
02845 }
02846
02847 fHistogram->GetXaxis()->SetTitle(xtitle.Data());
02848 fHistogram->GetYaxis()->SetTitle(ytitle.Data());
02849
02850 InitArgs(xv,fParams);
02851 for (i=1;i<=fNpx;i++) {
02852 xv[0] = fHistogram->GetBinCenter(i);
02853 fHistogram->SetBinContent(i,EvalPar(xv,fParams));
02854 }
02855
02856
02857 Double_t minimum = fHistogram->GetMinimumStored();
02858 Double_t maximum = fHistogram->GetMaximumStored();
02859 if (minimum <= 0 && gPad && gPad->GetLogy()) minimum = -1111;
02860 if (gPad && gPad->GetUymin() < fHistogram->GetMinimum() &&
02861 !fHistogram->TestBit(TH1::kIsZoomed)) minimum = -1111;
02862 if (minimum == -1111) {
02863 if (fHistogram->TestBit(TH1::kIsZoomed)) {
02864 minimum = fHistogram->GetYaxis()->GetXmin();
02865 } else {
02866 minimum = fMinimum;
02867
02868
02869 if (minimum == -1111) {
02870 Double_t hmin;
02871 if (optSAME) hmin = gPad->GetUymin();
02872 else hmin = fHistogram->GetMinimum();
02873 if (hmin > 0) {
02874 Double_t hmax;
02875 Double_t hminpos = hmin;
02876 if (optSAME) hmax = gPad->GetUymax();
02877 else hmax = fHistogram->GetMaximum();
02878 hmin -= 0.05*(hmax-hmin);
02879 if (hmin < 0) hmin = 0;
02880 if (hmin <= 0 && gPad && gPad->GetLogy()) hmin = hminpos;
02881 minimum = hmin;
02882 }
02883 }
02884 }
02885 fHistogram->SetMinimum(minimum);
02886 }
02887 if (maximum == -1111) {
02888 if (fHistogram->TestBit(TH1::kIsZoomed)) {
02889 maximum = fHistogram->GetYaxis()->GetXmax();
02890 } else {
02891 maximum = fMaximum;
02892 }
02893 fHistogram->SetMaximum(maximum);
02894 }
02895 fHistogram->SetBit(TH1::kNoStats);
02896 fHistogram->SetLineColor(GetLineColor());
02897 fHistogram->SetLineStyle(GetLineStyle());
02898 fHistogram->SetLineWidth(GetLineWidth());
02899 fHistogram->SetFillColor(GetFillColor());
02900 fHistogram->SetFillStyle(GetFillStyle());
02901 fHistogram->SetMarkerColor(GetMarkerColor());
02902 fHistogram->SetMarkerStyle(GetMarkerStyle());
02903 fHistogram->SetMarkerSize(GetMarkerSize());
02904
02905
02906 if (!gPad) return;
02907 if (opt.Length() == 0) fHistogram->Paint("lf");
02908 else if (optSAME) fHistogram->Paint("lfsame");
02909 else fHistogram->Paint(option);
02910 }
02911
02912
02913
02914 void TF1::Print(Option_t *option) const
02915 {
02916
02917
02918 TFormula::Print(option);
02919 if (fHistogram) fHistogram->Print(option);
02920 }
02921
02922
02923
02924 void TF1::ReleaseParameter(Int_t ipar)
02925 {
02926
02927
02928
02929 if (ipar < 0 || ipar > fNpar-1) return;
02930 SetParLimits(ipar,0,0);
02931 }
02932
02933
02934
02935 void TF1::Save(Double_t xmin, Double_t xmax, Double_t, Double_t, Double_t, Double_t)
02936 {
02937
02938
02939 if (fSave != 0) {delete [] fSave; fSave = 0;}
02940 if (fParent && fParent->InheritsFrom(TH1::Class())) {
02941
02942 if ((xmin >0 && xmax > 0) && TMath::Abs(TMath::Log10(xmax/xmin) > TMath::Log10(fNpx))) {
02943 TH1 *h = (TH1*)fParent;
02944 Int_t bin1 = h->GetXaxis()->FindBin(xmin);
02945 Int_t bin2 = h->GetXaxis()->FindBin(xmax);
02946 fNsave = bin2-bin1+4;
02947 fSave = new Double_t[fNsave];
02948 Double_t xv[1];
02949 InitArgs(xv,fParams);
02950 for (Int_t i=bin1;i<=bin2;i++) {
02951 xv[0] = h->GetXaxis()->GetBinCenter(i);
02952 fSave[i-bin1] = EvalPar(xv,fParams);
02953 }
02954 fSave[fNsave-3] = xmin;
02955 fSave[fNsave-2] = xmax;
02956 fSave[fNsave-1] = xmax;
02957 return;
02958 }
02959 }
02960 fNsave = fNpx+3;
02961 if (fNsave <= 3) {fNsave=0; return;}
02962 fSave = new Double_t[fNsave];
02963 Double_t dx = (xmax-xmin)/fNpx;
02964 if (dx <= 0) {
02965 dx = (fXmax-fXmin)/fNpx;
02966 fNsave--;
02967 xmin = fXmin +0.5*dx;
02968 xmax = fXmax -0.5*dx;
02969 }
02970 Double_t xv[1];
02971 InitArgs(xv,fParams);
02972 for (Int_t i=0;i<=fNpx;i++) {
02973 xv[0] = xmin + dx*i;
02974 fSave[i] = EvalPar(xv,fParams);
02975 }
02976 fSave[fNpx+1] = xmin;
02977 fSave[fNpx+2] = xmax;
02978 }
02979
02980
02981
02982 void TF1::SavePrimitive(ostream &out, Option_t *option )
02983 {
02984
02985
02986 Int_t i;
02987 char quote = '"';
02988 out<<" "<<endl;
02989
02990 if (!fType) {
02991 out<<" TF1 *"<<GetName()<<" = new TF1("<<quote<<GetName()<<quote<<","<<quote<<GetTitle()<<quote<<","<<fXmin<<","<<fXmax<<");"<<endl;
02992 if (fNpx != 100) {
02993 out<<" "<<GetName()<<"->SetNpx("<<fNpx<<");"<<endl;
02994 }
02995 } else {
02996 out<<" TF1 *"<<GetName()<<" = new TF1("<<quote<<"*"<<GetName()<<quote<<","<<fXmin<<","<<fXmax<<","<<GetNpar()<<");"<<endl;
02997 out<<" //The original function : "<<GetTitle()<<" had originally been created by:" <<endl;
02998 out<<" //TF1 *"<<GetName()<<" = new TF1("<<quote<<GetName()<<quote<<","<<GetTitle()<<","<<fXmin<<","<<fXmax<<","<<GetNpar()<<");"<<endl;
02999 out<<" "<<GetName()<<"->SetRange("<<fXmin<<","<<fXmax<<");"<<endl;
03000 out<<" "<<GetName()<<"->SetName("<<quote<<GetName()<<quote<<");"<<endl;
03001 out<<" "<<GetName()<<"->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
03002 if (fNpx != 100) {
03003 out<<" "<<GetName()<<"->SetNpx("<<fNpx<<");"<<endl;
03004 }
03005 Double_t dx = (fXmax-fXmin)/fNpx;
03006 Double_t xv[1];
03007 InitArgs(xv,fParams);
03008 for (i=0;i<=fNpx;i++) {
03009 xv[0] = fXmin + dx*i;
03010 Double_t save = EvalPar(xv,fParams);
03011 out<<" "<<GetName()<<"->SetSavedPoint("<<i<<","<<save<<");"<<endl;
03012 }
03013 out<<" "<<GetName()<<"->SetSavedPoint("<<fNpx+1<<","<<fXmin<<");"<<endl;
03014 out<<" "<<GetName()<<"->SetSavedPoint("<<fNpx+2<<","<<fXmax<<");"<<endl;
03015 }
03016
03017 if (TestBit(kNotDraw)) {
03018 out<<" "<<GetName()<<"->SetBit(TF1::kNotDraw);"<<endl;
03019 }
03020 if (GetFillColor() != 0) {
03021 if (GetFillColor() > 228) {
03022 TColor::SaveColor(out, GetFillColor());
03023 out<<" "<<GetName()<<"->SetFillColor(ci);" << endl;
03024 } else
03025 out<<" "<<GetName()<<"->SetFillColor("<<GetFillColor()<<");"<<endl;
03026 }
03027 if (GetFillStyle() != 1001) {
03028 out<<" "<<GetName()<<"->SetFillStyle("<<GetFillStyle()<<");"<<endl;
03029 }
03030 if (GetMarkerColor() != 1) {
03031 if (GetMarkerColor() > 228) {
03032 TColor::SaveColor(out, GetMarkerColor());
03033 out<<" "<<GetName()<<"->SetMarkerColor(ci);" << endl;
03034 } else
03035 out<<" "<<GetName()<<"->SetMarkerColor("<<GetMarkerColor()<<");"<<endl;
03036 }
03037 if (GetMarkerStyle() != 1) {
03038 out<<" "<<GetName()<<"->SetMarkerStyle("<<GetMarkerStyle()<<");"<<endl;
03039 }
03040 if (GetMarkerSize() != 1) {
03041 out<<" "<<GetName()<<"->SetMarkerSize("<<GetMarkerSize()<<");"<<endl;
03042 }
03043 if (GetLineColor() != 1) {
03044 if (GetLineColor() > 228) {
03045 TColor::SaveColor(out, GetLineColor());
03046 out<<" "<<GetName()<<"->SetLineColor(ci);" << endl;
03047 } else
03048 out<<" "<<GetName()<<"->SetLineColor("<<GetLineColor()<<");"<<endl;
03049 }
03050 if (GetLineWidth() != 4) {
03051 out<<" "<<GetName()<<"->SetLineWidth("<<GetLineWidth()<<");"<<endl;
03052 }
03053 if (GetLineStyle() != 1) {
03054 out<<" "<<GetName()<<"->SetLineStyle("<<GetLineStyle()<<");"<<endl;
03055 }
03056 if (GetChisquare() != 0) {
03057 out<<" "<<GetName()<<"->SetChisquare("<<GetChisquare()<<");"<<endl;
03058 out<<" "<<GetName()<<"->SetNDF("<<GetNDF()<<");"<<endl;
03059 }
03060
03061 GetXaxis()->SaveAttributes(out,GetName(),"->GetXaxis()");
03062 GetYaxis()->SaveAttributes(out,GetName(),"->GetYaxis()");
03063
03064 Double_t parmin, parmax;
03065 for (i=0;i<fNpar;i++) {
03066 out<<" "<<GetName()<<"->SetParameter("<<i<<","<<GetParameter(i)<<");"<<endl;
03067 out<<" "<<GetName()<<"->SetParError("<<i<<","<<GetParError(i)<<");"<<endl;
03068 GetParLimits(i,parmin,parmax);
03069 out<<" "<<GetName()<<"->SetParLimits("<<i<<","<<parmin<<","<<parmax<<");"<<endl;
03070 }
03071 if (!strstr(option,"nodraw")) {
03072 out<<" "<<GetName()<<"->Draw("
03073 <<quote<<option<<quote<<");"<<endl;
03074 }
03075 }
03076
03077
03078
03079 void TF1::SetCurrent(TF1 *f1)
03080 {
03081
03082
03083
03084
03085 fgCurrent = f1;
03086 }
03087
03088
03089 void TF1::SetFitResult(const ROOT::Fit::FitResult & result, const Int_t* indpar )
03090 {
03091
03092
03093
03094
03095
03096
03097 if (result.IsEmpty()) {
03098 Warning("SetFitResult","Empty Fit result - nathing is set in TF1");
03099 return;
03100 }
03101 if (indpar == 0 && fNpar != (int) result.NPar() ) {
03102 Error("SetFitResult","Invalid Fit result passed - number of parameter is %d , different than TF1::GetNpar() = %d",fNpar,result.NPar());
03103 return;
03104 }
03105 if (result.Chi2() > 0)
03106 SetChisquare(result.Chi2() );
03107 else
03108 SetChisquare(result.MinFcnValue() );
03109
03110 SetNDF(result.Ndf() );
03111 SetNumberFitPoints(result.Ndf() + result.NFreeParameters() );
03112
03113
03114 for (Int_t i = 0; i < fNpar; ++i) {
03115 Int_t ipar = (indpar != 0) ? indpar[i] : i;
03116 if (ipar < 0) continue;
03117 fParams[i] = result.Parameter(ipar);
03118
03119 if (ipar < (int) result.Errors().size() )
03120 fParErrors[i] = result.Error(ipar);
03121 }
03122
03123 Update();
03124
03125 }
03126
03127
03128
03129 void TF1::SetMaximum(Double_t maximum)
03130 {
03131
03132
03133
03134
03135 fMaximum = maximum;
03136 if (fHistogram) fHistogram->SetMaximum(maximum);
03137 if (gPad) gPad->Modified();
03138 }
03139
03140
03141
03142 void TF1::SetMinimum(Double_t minimum)
03143 {
03144
03145
03146
03147
03148 fMinimum = minimum;
03149 if (fHistogram) fHistogram->SetMinimum(minimum);
03150 if (gPad) gPad->Modified();
03151 }
03152
03153
03154
03155 void TF1::SetNDF(Int_t ndf)
03156 {
03157
03158
03159
03160 fNDF = ndf;
03161 }
03162
03163
03164
03165 void TF1::SetNpx(Int_t npx)
03166 {
03167
03168
03169
03170
03171
03172
03173
03174 if (npx < 4) {
03175 Warning("SetNpx","Number of points must be >4 && < 100000, fNpx set to 4");
03176 fNpx = 4;
03177 } else if(npx > 100000) {
03178 Warning("SetNpx","Number of points must be >4 && < 100000, fNpx set to 100000");
03179 fNpx = 100000;
03180 } else {
03181 fNpx = npx;
03182 }
03183 Update();
03184 }
03185
03186
03187
03188 void TF1::SetParError(Int_t ipar, Double_t error)
03189 {
03190
03191
03192 if (ipar < 0 || ipar > fNpar-1) return;
03193 fParErrors[ipar] = error;
03194 }
03195
03196
03197
03198 void TF1::SetParErrors(const Double_t *errors)
03199 {
03200
03201
03202
03203 if (!errors) return;
03204 for (Int_t i=0;i<fNpar;i++) fParErrors[i] = errors[i];
03205 }
03206
03207
03208
03209 void TF1::SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
03210 {
03211
03212
03213
03214
03215
03216
03217 if (ipar < 0 || ipar > fNpar-1) return;
03218 Int_t i;
03219 if (!fParMin) {fParMin = new Double_t[fNpar]; for (i=0;i<fNpar;i++) fParMin[i]=0;}
03220 if (!fParMax) {fParMax = new Double_t[fNpar]; for (i=0;i<fNpar;i++) fParMax[i]=0;}
03221 fParMin[ipar] = parmin;
03222 fParMax[ipar] = parmax;
03223 }
03224
03225
03226
03227 void TF1::SetRange(Double_t xmin, Double_t xmax)
03228 {
03229
03230
03231
03232
03233
03234 fXmin = xmin;
03235 fXmax = xmax;
03236 Update();
03237 }
03238
03239
03240
03241 void TF1::SetSavedPoint(Int_t point, Double_t value)
03242 {
03243
03244
03245 if (!fSave) {
03246 fNsave = fNpx+3;
03247 fSave = new Double_t[fNsave];
03248 }
03249 if (point < 0 || point >= fNsave) return;
03250 fSave[point] = value;
03251 }
03252
03253
03254
03255 void TF1::SetTitle(const char *title)
03256 {
03257
03258
03259
03260
03261
03262 if (!title) return;
03263 fTitle = title;
03264 if (!fHistogram) return;
03265 fHistogram->SetTitle(title);
03266 if (gPad) gPad->Modified();
03267 }
03268
03269
03270
03271 void TF1::Streamer(TBuffer &b)
03272 {
03273
03274
03275 if (b.IsReading()) {
03276 UInt_t R__s, R__c;
03277 Version_t v = b.ReadVersion(&R__s, &R__c);
03278 if (v > 4) {
03279 b.ReadClassBuffer(TF1::Class(), this, v, R__s, R__c);
03280 if (v == 5 && fNsave > 0) {
03281
03282 Int_t np = fNsave - 3;
03283 fSave[np] = fSave[np-1];
03284 fSave[np+1] = fXmin;
03285 fSave[np+2] = fXmax;
03286 }
03287 return;
03288 }
03289
03290 TFormula::Streamer(b);
03291 TAttLine::Streamer(b);
03292 TAttFill::Streamer(b);
03293 TAttMarker::Streamer(b);
03294 if (v < 4) {
03295 Float_t xmin,xmax;
03296 b >> xmin; fXmin = xmin;
03297 b >> xmax; fXmax = xmax;
03298 } else {
03299 b >> fXmin;
03300 b >> fXmax;
03301 }
03302 b >> fNpx;
03303 b >> fType;
03304 b >> fChisquare;
03305 b.ReadArray(fParErrors);
03306 if (v > 1) {
03307 b.ReadArray(fParMin);
03308 b.ReadArray(fParMax);
03309 } else {
03310 fParMin = new Double_t[fNpar+1];
03311 fParMax = new Double_t[fNpar+1];
03312 }
03313 b >> fNpfits;
03314 if (v == 1) {
03315 b >> fHistogram;
03316 delete fHistogram; fHistogram = 0;
03317 }
03318 if (v > 1) {
03319 if (v < 4) {
03320 Float_t minimum,maximum;
03321 b >> minimum; fMinimum =minimum;
03322 b >> maximum; fMaximum =maximum;
03323 } else {
03324 b >> fMinimum;
03325 b >> fMaximum;
03326 }
03327 }
03328 if (v > 2) {
03329 b >> fNsave;
03330 if (fNsave > 0) {
03331 fSave = new Double_t[fNsave+10];
03332 b.ReadArray(fSave);
03333
03334 fSave[fNsave] = fSave[fNsave-1];
03335 fSave[fNsave+1] = fSave[fNsave+2];
03336 fSave[fNsave+2] = fSave[fNsave+3];
03337 fNsave += 3;
03338 } else fSave = 0;
03339 }
03340 b.CheckByteCount(R__s, R__c, TF1::IsA());
03341
03342
03343 } else {
03344 Int_t saved = 0;
03345 if (fType > 0 && fNsave <= 0) { saved = 1; Save(fXmin,fXmax,0,0,0,0);}
03346
03347 b.WriteClassBuffer(TF1::Class(),this);
03348
03349 if (saved) {delete [] fSave; fSave = 0; fNsave = 0;}
03350 }
03351 }
03352
03353
03354
03355 void TF1::Update()
03356 {
03357
03358
03359
03360 delete fHistogram;
03361 fHistogram = 0;
03362 if (fIntegral) {
03363 delete [] fIntegral; fIntegral = 0;
03364 delete [] fAlpha; fAlpha = 0;
03365 delete [] fBeta; fBeta = 0;
03366 delete [] fGamma; fGamma = 0;
03367 }
03368 }
03369
03370
03371
03372 void TF1::RejectPoint(Bool_t reject)
03373 {
03374
03375
03376
03377
03378
03379
03380 fgRejectPoint = reject;
03381 }
03382
03383
03384
03385 Bool_t TF1::RejectedPoint()
03386 {
03387
03388
03389 return fgRejectPoint;
03390 }
03391
03392
03393 Double_t TF1::Moment(Double_t n, Double_t a, Double_t b, const Double_t *params, Double_t epsilon)
03394 {
03395
03396
03397
03398
03399
03400
03401
03402 TF1_EvalWrapper func(this, params, kTRUE, n);
03403
03404 ROOT::Math::GaussIntegrator giod;
03405
03406 giod.SetFunction(func);
03407 giod.SetRelTolerance(epsilon);
03408
03409 Double_t norm = giod.Integral(a, b);
03410 if (norm == 0) {
03411 Error("Moment", "Integral zero over range");
03412 return 0;
03413 }
03414
03415
03416
03417 ROOT::Math::Functor1D xnfunc( &func, &TF1_EvalWrapper::EvalNMom);
03418 giod.SetFunction(xnfunc);
03419
03420 Double_t res = giod.Integral(a,b)/norm;
03421
03422 return res;
03423 }
03424
03425
03426
03427 Double_t TF1::CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params, Double_t epsilon)
03428 {
03429
03430
03431
03432
03433
03434
03435 TF1_EvalWrapper func(this, params, kTRUE, n);
03436
03437 ROOT::Math::GaussIntegrator giod;
03438
03439 giod.SetFunction(func);
03440 giod.SetRelTolerance(epsilon);
03441
03442 Double_t norm = giod.Integral(a, b);
03443 if (norm == 0) {
03444 Error("Moment", "Integral zero over range");
03445 return 0;
03446 }
03447
03448
03449
03450 ROOT::Math::Functor1D xfunc( &func, &TF1_EvalWrapper::EvalFirstMom);
03451 giod.SetFunction(xfunc);
03452
03453
03454 Double_t xbar = giod.Integral(a,b)/norm;
03455
03456
03457 func.fX0 = xbar;
03458 ROOT::Math::Functor1D xnfunc( &func, &TF1_EvalWrapper::EvalNMom);
03459 giod.SetFunction(xnfunc);
03460
03461 Double_t res = giod.Integral(a,b)/norm;
03462 return res;
03463 }
03464
03465
03466
03467
03468
03469 #ifdef INTHEFUTURE
03470 void TF1::CalcGaussLegendreSamplingPoints(TGraph *g, Double_t eps)
03471 {
03472
03473
03474
03475 if (!g) return;
03476 CalcGaussLegendreSamplingPoints(g->GetN(), g->GetX(), g->GetY(), eps);
03477 }
03478
03479
03480
03481 TGraph *TF1::CalcGaussLegendreSamplingPoints(Int_t num, Double_t eps)
03482 {
03483
03484
03485
03486
03487
03488
03489 if (num<=0)
03490 return 0;
03491
03492 TGraph *g = new TGraph(num);
03493 CalcGaussLegendreSamplingPoints(g->GetN(), g->GetX(), g->GetY(), eps);
03494 return g;
03495 }
03496 #endif
03497
03498
03499
03500 void TF1::CalcGaussLegendreSamplingPoints(Int_t num, Double_t *x, Double_t *w, Double_t eps)
03501 {
03502
03503
03504
03505
03506
03507
03508
03509
03510
03511
03512
03513
03514
03515
03516
03517
03518
03519
03520
03521 ROOT::Math::GaussLegendreIntegrator gli(num,eps);
03522 gli.GetWeightVectors(x, w);
03523
03524
03525 }