00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 #include "TROOT.h"
00102 #include "TFumili.h"
00103 #include "TMath.h"
00104 #include "TF1.h"
00105 #include "TH1.h"
00106 #include "TGraphAsymmErrors.h"
00107
00108 #include "Riostream.h"
00109
00110
00111 extern void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00112 extern void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00113 extern void GraphFitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00114
00115
00116 ClassImp(TFumili);
00117
00118 TFumili *gFumili=0;
00119
00120
00121 static const Double_t gMAXDOUBLE=1e300;
00122 static const Double_t gMINDOUBLE=-1e300;
00123
00124
00125
00126 TFumili::TFumili(Int_t maxpar)
00127 {
00128
00129
00130 fMaxParam = TMath::Max(maxpar,25);
00131 if (fMaxParam>200) fMaxParam=200;
00132 fMaxParam2 *= fMaxParam;
00133 BuildArrays();
00134
00135 fNumericDerivatives = true;
00136 fLogLike = false;
00137 fNpar = fMaxParam;
00138 fGRAD = false;
00139 fWARN = true;
00140 fDEBUG = false;
00141 fNlog = 0;
00142 fSumLog = 0;
00143 fNED1 = 0;
00144 fNED2 = 0;
00145 fNED12 = fNED1+fNED2;
00146 fEXDA = 0;
00147 fFCN = 0;
00148 fNfcn = 0;
00149 fRP = 1.e-15;
00150 fS = 1e10;
00151 fEPS =0.01;
00152 fENDFLG = 0;
00153 fNlimMul = 2;
00154 fNmaxIter= 150;
00155 fNstepDec= 3;
00156 fLastFixed = -1;
00157
00158 SetName("Fumili");
00159 gFumili = this;
00160 gROOT->GetListOfSpecials()->Add(gFumili);
00161 }
00162
00163
00164 void TFumili::BuildArrays(){
00165
00166
00167
00168 fCmPar = new Double_t[fMaxParam];
00169 fA = new Double_t[fMaxParam];
00170 fPL0 = new Double_t[fMaxParam];
00171 fPL = new Double_t[fMaxParam];
00172 fParamError = new Double_t[fMaxParam];
00173 fDA = new Double_t[fMaxParam];
00174 fAMX = new Double_t[fMaxParam];
00175 fAMN = new Double_t[fMaxParam];
00176 fR = new Double_t[fMaxParam];
00177 fDF = new Double_t[fMaxParam];
00178 fGr = new Double_t[fMaxParam];
00179 fANames = new TString[fMaxParam];
00180
00181
00182
00183 Int_t zSize = fMaxParam*(fMaxParam+1)/2;
00184 fZ0 = new Double_t[zSize];
00185 fZ = new Double_t[zSize];
00186
00187 for (Int_t i=0;i<fMaxParam;i++) {
00188 fA[i] =0.;
00189 fDF[i]=0.;
00190 fAMN[i]=gMINDOUBLE;
00191 fAMX[i]=gMAXDOUBLE;
00192 fPL0[i]=.1;
00193 fPL[i] =.1;
00194 fParamError[i]=0.;
00195 fANames[i]=Form("%d",i);
00196 }
00197 }
00198
00199
00200
00201 TFumili::~TFumili() {
00202
00203
00204
00205 DeleteArrays();
00206 gROOT->GetListOfSpecials()->Remove(this);
00207 if (gFumili == this) gFumili = 0;
00208 }
00209
00210
00211 Double_t TFumili::Chisquare(Int_t npar, Double_t *params) const
00212 {
00213
00214
00215 Double_t amin = 0;
00216 H1FitChisquareFumili(npar,params,amin,params,1);
00217 return 2*amin;
00218 }
00219
00220
00221
00222 void TFumili::Clear(Option_t *)
00223 {
00224
00225
00226
00227
00228
00229
00230
00231 fNpar = fMaxParam;
00232 fNfcn = 0;
00233 for (Int_t i=0;i<fNpar;i++) {
00234 fA[i] =0.;
00235 fDF[i] =0.;
00236 fPL0[i] =.1;
00237 fPL[i] =.1;
00238 fAMN[i] = gMINDOUBLE;
00239 fAMX[i] = gMAXDOUBLE;
00240 fParamError[i]=0.;
00241 fANames[i]=Form("%d",i);
00242 }
00243 }
00244
00245
00246
00247 void TFumili::DeleteArrays(){
00248
00249
00250
00251 delete[] fCmPar;
00252 delete[] fANames;
00253 delete[] fDF;
00254
00255 delete[] fZ0;
00256 delete[] fZ;
00257 delete[] fGr;
00258 delete[] fA;
00259 delete[] fPL0;
00260 delete[] fPL;
00261 delete[] fDA;
00262 delete[] fAMN;
00263 delete[] fAMX;
00264 delete[] fParamError;
00265 delete[] fR;
00266 }
00267
00268
00269
00270 void TFumili::Derivatives(Double_t *df,Double_t *fX){
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282 Double_t ff,ai,hi,y,pi;
00283 y = EvalTFN(df,fX);
00284 for (Int_t i=0;i<fNpar;i++) {
00285 df[i]=0;
00286 if(fPL0[i]>0.) {
00287 ai = fA[i];
00288 hi = 0.01*fPL0[i];
00289 pi = fRP*TMath::Abs(ai);
00290 if (hi<pi) hi = pi;
00291 fA[i] = ai+hi;
00292
00293 if (fA[i]>fAMX[i]) {
00294 fA[i] = ai-hi;
00295 hi = -hi;
00296 if (fA[i]<fAMN[i]) {
00297 fA[i] = fAMX[i];
00298 hi = fAMX[i]-ai;
00299 if (fAMN[i]-ai+hi<0) {
00300 fA[i]=fAMN[i];
00301 hi=fAMN[i]-ai;
00302 }
00303 }
00304 }
00305 ff = EvalTFN(df,fX);
00306 df[i] = (ff-y)/hi;
00307 fA[i] = ai;
00308 }
00309 }
00310 }
00311
00312
00313
00314
00315 Int_t TFumili::Eval(Int_t& npar, Double_t *grad, Double_t &fval, Double_t *par, Int_t flag)
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 if (fFCN) (*fFCN)(npar,grad,fval,par,flag);
00343 return npar;
00344 }
00345
00346
00347
00348 Double_t TFumili::EvalTFN(Double_t * , Double_t *X)
00349 {
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359 TF1 *f1 = (TF1*)fUserFunc;
00360 return f1->EvalPar(X,fA);
00361
00362
00363 }
00364
00365
00366 Int_t TFumili::ExecuteCommand(const char *command, Double_t *args, Int_t nargs){
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382 TString comand = command;
00383 static TString clower = "abcdefghijklmnopqrstuvwxyz";
00384 static TString cupper = "ABCDEFGHIJKLMNOPQRSTUVWXYZ";
00385 const Int_t nntot = 40;
00386 const char *cname[nntot] = {
00387 "MINImize ",
00388 "SEEk ",
00389 "SIMplex ",
00390 "MIGrad ",
00391 "MINOs ",
00392 "SET xxx ",
00393 "SHOw xxx ",
00394 "TOP of pag",
00395 "fiX ",
00396 "REStore ",
00397 "RELease ",
00398 "SCAn ",
00399 "CONtour ",
00400 "HESse ",
00401 "SAVe ",
00402 "IMProve ",
00403 "CALl fcn ",
00404 "STAndard ",
00405 "END ",
00406 "EXIt ",
00407 "RETurn ",
00408 "CLEar ",
00409 "HELP ",
00410 "MNContour ",
00411 "STOp ",
00412 "JUMp ",
00413 " ",
00414 " ",
00415 "FUMili ",
00416 " ",
00417 " ",
00418 " ",
00419 " ",
00420 "COVARIANCE",
00421 "PRINTOUT ",
00422 "GRADIENT ",
00423 "MATOUT ",
00424 "ERROR DEF ",
00425 "LIMITS ",
00426 "PUNCH "};
00427
00428
00429 fCword = comand;
00430 fCword.ToUpper();
00431 if (nargs<=0) fCmPar[0] = 0;
00432 Int_t i;
00433 for(i=0;i<fMaxParam;i++) {
00434 if(i<=nargs) fCmPar[i] = args[i];
00435 }
00436
00437
00438
00439
00440
00441
00442
00443
00444 TString ctemp = fCword(0,3);
00445 Int_t ind;
00446 for (ind = 0; ind < nntot; ++ind) {
00447 if (strncmp(ctemp.Data(),cname[ind],3) == 0) break;
00448 }
00449 if (ind==nntot) return -3;
00450 if (fCword(0,4) == "MINO") ind=3;
00451 switch (ind) {
00452 case 0: case 3: case 2: case 28:
00453
00454
00455 if(nargs>=1)
00456 fNmaxIter=TMath::Max(Int_t(fCmPar[0]),fNmaxIter);
00457 if(nargs==2)
00458 fEPS=fCmPar[1];
00459 return Minimize();
00460 case 1:
00461
00462 return -10;
00463
00464 case 4:
00465 return -10;
00466
00467 case 5: case 6:
00468 return ExecuteSetCommand(nargs);
00469
00470 case 7:
00471 Printf("1");
00472 return 0;
00473 case 8:
00474 if (nargs<1) return -1;
00475 for (i=0;i<nargs;i++) {
00476 Int_t parnum = Int_t(fCmPar[i])-1;
00477 FixParameter(parnum);
00478 }
00479 return 0;
00480 case 9:
00481 if (nargs<1) return 0;
00482 if(fCmPar[0]==0.)
00483 for (i=0;i<fNpar;i++)
00484 ReleaseParameter(i);
00485 else
00486 if(fCmPar[0]==1.) {
00487 ReleaseParameter(fLastFixed);
00488 cout <<fLastFixed<<endl;
00489 }
00490 return 0;
00491 case 10:
00492 if (nargs<1) return -1;
00493 for (i=0;i<nargs;i++) {
00494 Int_t parnum = Int_t(fCmPar[i])-1;
00495 ReleaseParameter(parnum);
00496 }
00497 return 0;
00498 case 11:
00499 return -10;
00500 case 12:
00501 return -10;
00502
00503 case 13:
00504 return -10;
00505 case 14:
00506 Printf("SAVe command is obsolete");
00507 return -10;
00508 case 15:
00509 return -10;
00510 case 16:
00511 {if(nargs<1) return -1;
00512 Int_t flag = Int_t(fCmPar[0]);
00513 Double_t fval;
00514 Eval(fNpar,fGr,fval,fA,flag);
00515 return 0;}
00516 case 17:
00517 return 0;
00518 case 18: case 19:
00519 case 20: case 24: {
00520 Double_t fval;
00521 Int_t flag = 3;
00522 Eval(fNpar,fGr,fval,fA,flag);
00523 return 0;
00524 }
00525 case 21:
00526 Clear();
00527 return 0;
00528 case 22:
00529 case 23:
00530 case 25:
00531 return -10;
00532 case 26: case 27: case 29: case 30: case 31: case 32:
00533 return 0;
00534 case 33: case 34: case 35: case 36: case 37: case 38:
00535 case 39:
00536 Printf("Obsolete command. Use corresponding SET command instead");
00537 return -10;
00538 default:
00539 break;
00540 }
00541 return 0;
00542 }
00543
00544
00545
00546
00547 Int_t TFumili::ExecuteSetCommand(Int_t nargs){
00548
00549
00550
00551
00552 static Int_t nntot = 30;
00553 static const char *cname[30] = {
00554 "FCN value ",
00555 "PARameters",
00556 "LIMits ",
00557 "COVariance",
00558 "CORrelatio",
00559 "PRInt levl",
00560 "NOGradient",
00561 "GRAdient ",
00562 "ERRor def ",
00563 "INPut file",
00564 "WIDth page",
00565 "LINes page",
00566 "NOWarnings",
00567 "WARnings ",
00568 "RANdom gen",
00569 "TITle ",
00570 "STRategy ",
00571 "EIGenvalue",
00572 "PAGe throw",
00573 "MINos errs",
00574 "EPSmachine",
00575 "OUTputfile",
00576 "BATch ",
00577 "INTeractiv",
00578 "VERsion ",
00579 "reserve ",
00580 "NODebug ",
00581 "DEBug ",
00582 "SHOw ",
00583 "SET "};
00584
00585 TString cfname, cmode, ckind, cwarn, copt, ctemp, ctemp2;
00586 Int_t i, ind;
00587 Bool_t setCommand=kFALSE;
00588 for (ind = 0; ind < nntot; ++ind) {
00589 ctemp = cname[ind];
00590 ckind = ctemp(0,3);
00591 ctemp2 = fCword(4,6);
00592 if (strstr(ctemp2.Data(),ckind.Data())) break;
00593 }
00594 ctemp2 = fCword(0,3);
00595 if(ctemp2.Contains("SET")) setCommand=true;
00596 if(ctemp2.Contains("HEL") || ctemp2.Contains("SHO")) setCommand=false;
00597
00598 if (ind>=nntot) return -3;
00599
00600 switch(ind) {
00601 case 0:
00602 if(!setCommand) Printf("FCN=%f",fS);
00603 return 0;
00604 case 1:
00605 {
00606 if (nargs<2 && setCommand) return -1;
00607 Int_t parnum;
00608 Double_t val;
00609 if(setCommand) {
00610 parnum = Int_t(fCmPar[0])-1;
00611 val= fCmPar[1];
00612 if(parnum<0 || parnum>=fNpar) return -2;
00613 fA[parnum] = val;
00614 } else {
00615 if (nargs>0) {
00616 parnum = Int_t(fCmPar[0])-1;
00617 if(parnum<0 || parnum>=fNpar) return -2;
00618 Printf("Parameter %s = %E",fANames[parnum].Data(),fA[parnum]);
00619 } else
00620 for (i=0;i<fNpar;i++)
00621 Printf("Parameter %s = %E",fANames[i].Data(),fA[i]);
00622
00623 }
00624 return 0;
00625 }
00626 case 2:
00627 {
00628 Int_t parnum;
00629 Double_t lolim,uplim;
00630 if (nargs<1) {
00631 for(i=0;i<fNpar;i++)
00632 if(setCommand) {
00633 fAMN[i] = gMINDOUBLE;
00634 fAMX[i] = gMAXDOUBLE;
00635 } else
00636 Printf("Limits for param %s: Low=%E, High=%E",
00637 fANames[i].Data(),fAMN[i],fAMX[i]);
00638 } else {
00639 parnum = Int_t(fCmPar[0])-1;
00640 if(parnum<0 || parnum>=fNpar)return -1;
00641 if(setCommand) {
00642 if(nargs>2) {
00643 lolim = fCmPar[1];
00644 uplim = fCmPar[2];
00645 if(uplim==lolim) return -1;
00646 if(lolim>uplim) {
00647 Double_t tmp = lolim;
00648 lolim = uplim;
00649 uplim = tmp;
00650 }
00651 } else {
00652 lolim = gMINDOUBLE;
00653 uplim = gMAXDOUBLE;
00654 }
00655 fAMN[parnum] = lolim;
00656 fAMX[parnum] = uplim;
00657 } else
00658 Printf("Limits for param %s Low=%E, High=%E",
00659 fANames[parnum].Data(),fAMN[parnum],fAMX[parnum]);
00660 }
00661 return 0;
00662 }
00663 case 3:
00664 {
00665 if(setCommand) return 0;
00666 Printf("\nCovariant matrix ");
00667 Int_t l = 0,nn=0,nnn=0;
00668 for (i=0;i<fNpar;i++) if(fPL0[i]>0.) nn++;
00669 for (i=0;i<nn;i++) {
00670 for(;fPL0[nnn]<=0.;nnn++) { }
00671 printf("%5s: ",fANames[nnn++].Data());
00672 for (Int_t j=0;j<=i;j++)
00673 printf("%11.2E",fZ[l++]);
00674 cout<<endl;
00675 }
00676 cout<<endl;
00677 return 0;
00678 }
00679 case 4:
00680 if(setCommand) return 0;
00681 Printf("\nGlobal correlation factors (maximum correlation of the parameter\n with arbitrary linear combination of other parameters)");
00682 for(i=0;i<fNpar;i++) {
00683 printf("%5s: ",fANames[i].Data());
00684 printf("%11.3E\n",TMath::Sqrt(1-1/((fR[i]!=0.)?fR[i]:1.)) );
00685 }
00686 cout<<endl;
00687 return 0;
00688 case 5:
00689 return -10;
00690 case 6:
00691 if(!setCommand) return 0;
00692 fGRAD = false;
00693 return 0;
00694 case 7:
00695 if(!setCommand) return 0;
00696 fGRAD = true;
00697 return 0;
00698 case 8:
00699 return 0;
00700 case 9:
00701 return -10;
00702 case 10:
00703 return -10;
00704 case 11:
00705 return -10;
00706 case 12:
00707 if(!setCommand) return 0;
00708 fWARN = false;
00709 return 0;
00710 case 13:
00711 if(!setCommand) return 0;
00712 fWARN = true;
00713 return 0;
00714 case 14:
00715 return -10;
00716 case 15:
00717 return 0;
00718 case 16:
00719 return 0;
00720 case 17:
00721 return -10;
00722 case 18:
00723 return 0;
00724 case 19:
00725 return -10;
00726 case 20:
00727 if(!setCommand) {
00728 Printf("Relative floating point presicion RP=%E",fRP);
00729 } else
00730 if (nargs>0) {
00731 Double_t pres=fCmPar[0];
00732 if (pres<1e-5 && pres>1e-34) fRP=pres;
00733 }
00734 return 0;
00735 case 21:
00736 return -10;
00737 case 22:
00738 return 0;
00739 case 23:
00740 return 0;
00741 case 24:
00742 if(setCommand) return 0;
00743 Printf("FUMILI-ROOT version 0.1");
00744 return 0;
00745 case 25:
00746 return 0;
00747 case 26:
00748 if(!setCommand) return 0;
00749 fDEBUG = false;
00750 return 0;
00751 case 27:
00752 if(!setCommand) return 0;
00753 fDEBUG = true;
00754 return 0;
00755 case 28:
00756 case 29:
00757 return -3;
00758 default:
00759 break;
00760 }
00761 return -3;
00762 }
00763
00764
00765 void TFumili::FixParameter(Int_t ipar) {
00766
00767
00768 if(ipar>=0 && ipar<fNpar && fPL0[ipar]>0.) {
00769 fPL0[ipar] = -fPL0[ipar];
00770 fLastFixed = ipar;
00771 }
00772 }
00773
00774
00775 Double_t *TFumili::GetCovarianceMatrix() const
00776 {
00777
00778
00779 return fZ;
00780
00781 }
00782
00783
00784 Double_t TFumili::GetCovarianceMatrixElement(Int_t i, Int_t j) const
00785 {
00786
00787
00788 if (!fZ) return 0;
00789 if (i < 0 || i >= fNpar || j < 0 || j >= fNpar) {
00790 Error("GetCovarianceMatrixElement","Illegal arguments i=%d, j=%d",i,j);
00791 return 0;
00792 }
00793 return fZ[j+fNpar*i];
00794 }
00795
00796
00797
00798 Int_t TFumili::GetNumberTotalParameters() const
00799 {
00800
00801
00802 return fNpar;
00803 }
00804
00805
00806 Int_t TFumili::GetNumberFreeParameters() const
00807 {
00808
00809
00810 Int_t nfree = fNpar;
00811 for (Int_t i=0;i<fNpar;i++) {
00812 if (IsFixed(i)) nfree--;
00813 }
00814 return nfree;
00815 }
00816
00817
00818 Double_t TFumili::GetParError(Int_t ipar) const
00819 {
00820
00821
00822 if (ipar<0 || ipar>=fNpar) return 0;
00823 else return fParamError[ipar];
00824 }
00825
00826
00827 Double_t TFumili::GetParameter(Int_t ipar) const
00828 {
00829
00830
00831 if (ipar<0 || ipar>=fNpar) return 0;
00832 else return fA[ipar];
00833 }
00834
00835
00836
00837 Int_t TFumili::GetParameter(Int_t ipar,char *cname,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const
00838 {
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848 if (ipar<0 || ipar>=fNpar) {
00849 value = 0;
00850 verr = 0;
00851 vlow = 0;
00852 vhigh = 0;
00853 return -1;
00854 }
00855 strcpy(cname,fANames[ipar].Data());
00856 value = fA[ipar];
00857 verr = fParamError[ipar];
00858 vlow = fAMN[ipar];
00859 vhigh = fAMX[ipar];
00860 return 0;
00861 }
00862
00863
00864 const char *TFumili::GetParName(Int_t ipar) const
00865 {
00866
00867
00868 if (ipar < 0 || ipar > fNpar) return "";
00869 return fANames[ipar].Data();
00870 }
00871
00872
00873 Int_t TFumili::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
00874 {
00875
00876
00877 eparab = 0;
00878 globcc = 0;
00879 if (ipar<0 || ipar>=fNpar) {
00880 eplus = 0;
00881 eminus = 0;
00882 return -1;
00883 }
00884 eplus=fParamError[ipar];
00885 eminus=-eplus;
00886 return 0;
00887 }
00888
00889
00890 Int_t TFumili::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
00891 {
00892
00893
00894
00895
00896
00897
00898 amin = 2*fS;
00899 edm = fGT;
00900 errdef = 0;
00901 nparx = fNpar;
00902 nvpar = 0;
00903 for(Int_t ii=0; ii<fNpar; ii++) {
00904 if(fPL0[ii]>0.) nvpar++;
00905 }
00906 return 0;
00907 }
00908
00909
00910
00911
00912 Double_t TFumili::GetSumLog(Int_t n)
00913 {
00914
00915
00916
00917 if (n < 0) return 0;
00918 if (n > fNlog) {
00919 if (fSumLog) delete [] fSumLog;
00920 fNlog = 2*n+1000;
00921 fSumLog = new Double_t[fNlog+1];
00922 Double_t fobs = 0;
00923 for (Int_t j=0;j<=fNlog;j++) {
00924 if (j > 1) fobs += TMath::Log(j);
00925 fSumLog[j] = fobs;
00926 }
00927 }
00928 if (fSumLog) return fSumLog[n];
00929 return 0;
00930 }
00931
00932
00933
00934
00935 void TFumili::InvertZ(Int_t n)
00936 {
00937
00938
00939
00940
00941
00942
00943 static Double_t am = 3.4e138;
00944 static Double_t rp = 5.0e-14;
00945 Double_t ap, aps, c, d;
00946 Double_t *r_1=fR;
00947 Double_t *pl_1=fPL;
00948 Double_t *z_1=fZ;
00949 Int_t i, k, l, ii, ki, li, kk, ni, ll, nk, nl, ir, lk;
00950 if (n < 1) {
00951 return;
00952 }
00953 --pl_1;
00954 --r_1;
00955 --z_1;
00956 aps = am / n;
00957 aps = sqrt(aps);
00958 ap = 1.0e0 / (aps * aps);
00959 ir = 0;
00960 for (i = 1; i <= n; ++i) {
00961 L1:
00962 ++ir;
00963 if (pl_1[ir] <= 0.0e0) goto L1;
00964 else goto L2;
00965 L2:
00966 ni = i * (i - 1) / 2;
00967 ii = ni + i;
00968 k = n + 1;
00969 if (z_1[ii] <= rp * TMath::Abs(r_1[ir]) || z_1[ii] <= ap) {
00970 goto L19;
00971 }
00972 z_1[ii] = 1.0e0 / sqrt(z_1[ii]);
00973 nl = ii - 1;
00974 L3:
00975 if (nl - ni <= 0) goto L5;
00976 else goto L4;
00977 L4:
00978 z_1[nl] *= z_1[ii];
00979 if (TMath::Abs(z_1[nl]) >= aps) {
00980 goto L16;
00981 }
00982 --nl;
00983 goto L3;
00984 L5:
00985 if (i - n >= 0) goto L12;
00986 else goto L6;
00987 L6:
00988 --k;
00989 nk = k * (k - 1) / 2;
00990 nl = nk;
00991 kk = nk + i;
00992 d = z_1[kk] * z_1[ii];
00993 c = d * z_1[ii];
00994 l = k;
00995 L7:
00996 ll = nk + l;
00997 li = nl + i;
00998 z_1[ll] -= z_1[li] * c;
00999 --l;
01000 nl -= l;
01001 if (l - i <= 0) goto L9;
01002 else goto L7;
01003 L8:
01004 ll = nk + l;
01005 li = ni + l;
01006 z_1[ll] -= z_1[li] * d;
01007 L9:
01008 --l;
01009 if (l <= 0) goto L10;
01010 else goto L8;
01011 L10:
01012 z_1[kk] = -c;
01013 if (k - i - 1 <= 0) goto L11;
01014 else goto L6;
01015 L11:
01016 ;
01017 }
01018 L12:
01019 for (i = 1; i <= n; ++i) {
01020 for (k = i; k <= n; ++k) {
01021 nl = k * (k - 1) / 2;
01022 ki = nl + i;
01023 d = 0.0e0;
01024 for (l = k; l <= n; ++l) {
01025 li = nl + i;
01026 lk = nl + k;
01027 d += z_1[li] * z_1[lk];
01028 nl += l;
01029 }
01030 ki = k * (k - 1) / 2 + i;
01031 z_1[ki] = d;
01032 }
01033 }
01034 L15:
01035 return;
01036 L16:
01037 k = i + nl - ii;
01038 ir = 0;
01039 for (i = 1; i <= k; ++i) {
01040 L17:
01041 ++ir;
01042 if (pl_1[ir] <= 0.0e0) {
01043 goto L17;
01044 }
01045 }
01046 L19:
01047 pl_1[ir] = -2.0e0;
01048 r_1[ir] = 0.0e0;
01049 fINDFLG[0] = ir - 1;
01050 goto L15;
01051 }
01052
01053
01054
01055
01056
01057 Bool_t TFumili::IsFixed(Int_t ipar) const
01058 {
01059
01060
01061 if(ipar < 0 || ipar >= fNpar) {
01062 Warning("IsFixed","Illegal parameter number :%d",ipar);
01063 return kFALSE;
01064 }
01065 if (fPL0[ipar] < 0) return kTRUE;
01066 else return kFALSE;
01067 }
01068
01069
01070
01071 Int_t TFumili::Minimize()
01072 {
01073
01074
01075
01076
01077
01078
01079
01080
01081
01082
01083
01084
01085
01086
01087
01088
01089
01090
01091
01092
01093
01094 Int_t i;
01095
01096 fINDFLG[2]=0;
01097
01098
01099
01100 Int_t parn;
01101
01102 if(fFCN) {
01103 Eval(parn,fGr,fS,fA,9); fNfcn++;
01104 }
01105 for( i = 0; i < fNpar; i++) {
01106 if(fA[i] > fAMX[i]) fA[i] = fAMX[i];
01107 if(fA[i] < fAMN[i]) fA[i] = fAMN[i];
01108 }
01109
01110 Int_t nn2, n, fixFLG, ifix1, fi, nn3, nn1, n0;
01111 Double_t t1;
01112 Double_t sp, t, olds=0;
01113 Double_t bi, aiMAX=0, amb;
01114 Double_t afix, sigi, akap;
01115 Double_t alambd, al, bm, abi, abm;
01116 Int_t l1, k, ifix;
01117
01118 nn2=0;
01119
01120
01121 n=fNpar;
01122 fixFLG=0;
01123
01124
01125 fENDFLG=0;
01126
01127
01128 fINDFLG[1] = 0;
01129 ifix1=-1;
01130 fi=0;
01131 nn3=0;
01132
01133
01134 for( i=0; i < n; i++) {
01135 fR[i]=0.;
01136 if ( fEPS > 0.) fParamError[i] = 0.;
01137 fPL[i] = fPL0[i];
01138 }
01139
01140 L3:
01141
01142 nn1 = 1;
01143 t1 = 1.;
01144
01145 L4:
01146
01147
01148 fS = 0.;
01149
01150 n0 = 0;
01151 for( i = 0; i < n; i++) {
01152 fGr[i]=0.;
01153 if (fPL0[i] > .0) {
01154 n0=n0+1;
01155
01156 if (fPL[i] > .0) fPL0[i]=fPL[i];
01157 }
01158 }
01159 Int_t nn0;
01160
01161 nn0 = n0*(n0+1)/2;
01162
01163
01164 for( i=0; i < nn0; i++) fZ[i]=0.;
01165
01166
01167 fINDFLG[0] = 0;
01168 Int_t ijkl=1;
01169
01170
01171 if(fFCN) {
01172 Eval(parn,fGr,fS,fA,2);
01173 fNfcn++;
01174 } else
01175 ijkl = SGZ();
01176 if(!ijkl) return 10;
01177 if (ijkl == -1) fINDFLG[0]=1;
01178
01179
01180 sp=fRP*TMath::Abs(fS);
01181
01182
01183 for( i=0; i < nn0; i++) fZ0[i] = fZ[i];
01184 if (nn3 > 0) {
01185 if (nn1 <= fNstepDec) {
01186 t=2.*(fS-olds-fGT);
01187 if (fINDFLG[0] == 0) {
01188 if (TMath::Abs(fS-olds) <= sp && -fGT <= sp) goto L19;
01189 if( 0.59*t < -fGT) goto L19;
01190 t = -fGT/t;
01191 if (t < 0.25 ) t = 0.25;
01192 }
01193 else t = 0.25;
01194 fGT = fGT*t;
01195 t1 = t1*t;
01196 nn2=0;
01197 for( i = 0; i < n; i++) {
01198 if (fPL[i] > 0.) {
01199 fA[i]=fA[i]-fDA[i];
01200 fPL[i]=fPL[i]*t;
01201 fDA[i]=fDA[i]*t;
01202 fA[i]=fA[i]+fDA[i];
01203 }
01204 }
01205 nn1=nn1+1;
01206 goto L4;
01207 }
01208 }
01209
01210 L19:
01211
01212 if(fINDFLG[0] != 0) {
01213 fENDFLG=-4;
01214 printf("trying to execute an illegal junp at L85\n");
01215
01216 }
01217
01218
01219 Int_t k1, k2, i1, j, l;
01220 k1 = 1;
01221 k2 = 1;
01222 i1 = 1;
01223
01224
01225 for( i = 0; i < n; i++) {
01226 if (fPL0[i] > .0) {
01227
01228 if (fPL[i] == 0.) fPL[i]=fPL0[i];
01229 if (fPL[i] > .0) {
01230
01231
01232 if ((fA[i] >= fAMX[i] && fGr[i] < 0.) ||
01233 (fA[i] <= fAMN[i] && fGr[i] > 0.)) {
01234 fPL[i] = 0.;
01235 k1 = k1 + i1;
01236
01237
01238 } else {
01239 for( j=0; j <= i; j++) {
01240 if (fPL0[j] > .0) {
01241
01242
01243 if (fPL[j] > .0) {
01244 fZ[k2 -1] = fZ0[k1 -1];
01245 k2=k2+1;
01246 }
01247 k1=k1+1;
01248 }
01249 }
01250 }
01251 }
01252 else k1 = k1 + i1;
01253 i1=i1+1;
01254 }
01255 }
01256
01257
01258 i1 = 1;
01259 l = 1;
01260 for( i = 0; i < n; i++) {
01261 if (fPL[i] > .0) {
01262 fR[i] = fZ[l - 1];
01263 i1 = i1+1;
01264 l = l + i1;
01265 }
01266 }
01267
01268 n0 = i1 - 1;
01269 InvertZ(n0);
01270
01271
01272 if (fINDFLG[0] != 0) {
01273
01274 fINDFLG[0] = 0;
01275 fINDFLG[1] = 1;
01276 fixFLG = fixFLG + 1;
01277 fi = 0;
01278 goto L19;
01279 }
01280
01281
01282 i1 = 1;
01283 for( i = 0; i < n; i++) {
01284 fDA[i]=0.;
01285 if (fPL[i] > .0) {
01286 l1=1;
01287 for( l = 0; l < n; l++) {
01288 if (fPL[l] > .0) {
01289
01290
01291 if (i1 <= l1 ) k=l1*(l1-1)/2+i1;
01292 else k=i1*(i1-1)/2+l1;
01293
01294 fDA[i]=fDA[i]-fGr[l]*fZ[k - 1];
01295 l1=l1+1;
01296 }
01297 }
01298 i1=i1+1;
01299 }
01300 }
01301
01302
01303 afix=0.;
01304 ifix = -1;
01305 i1 = 1;
01306 l = i1;
01307 for( i = 0; i < n; i++)
01308 if (fPL[i] > .0) {
01309 sigi = TMath::Sqrt(TMath::Abs(fZ[l - 1]));
01310 fR[i] = fR[i]*fZ[l - 1];
01311 if (fEPS > .0) fParamError[i]=sigi;
01312 if ((fA[i] >= fAMX[i] && fDA[i] > 0.) || (fA[i] <= fAMN[i]
01313 && fDA[i] < .0)) {
01314
01315
01316 akap = TMath::Abs(fDA[i]/sigi);
01317
01318 if (akap > afix) {
01319 afix=akap;
01320 ifix=i;
01321 ifix1=i;
01322 }
01323 }
01324 i1=i1+1;
01325 l=l+i1;
01326 }
01327 if (ifix != -1) {
01328
01329
01330 fPL[ifix] = -1.;
01331 fixFLG = fixFLG + 1;
01332 fi = 0;
01333
01334 goto L19;
01335 }
01336
01337
01338
01339 alambd = 1.;
01340 fAKAPPA = 0.;
01341 Int_t imax;
01342 imax = -1;
01343
01344
01345 for( i = 0; i < n; i++) {
01346 if (fPL[i] > .0) {
01347 bm = fAMX[i] - fA[i];
01348 abi = fA[i] + fPL[i];
01349 abm = fAMX[i];
01350 if (fDA[i] <= .0) {
01351 bm = fA[i] - fAMN[i];
01352 abi = fA[i] - fPL[i];
01353 abm = fAMN[i];
01354 }
01355 bi = fPL[i];
01356
01357
01358 if ( bi > bm) {
01359 bi = bm;
01360 abi = abm;
01361 }
01362
01363 if ( TMath::Abs(fDA[i]) > bi) {
01364
01365 al = TMath::Abs(bi/fDA[i]);
01366 if (alambd > al) {
01367 imax=i;
01368 aiMAX=abi;
01369 alambd=al;
01370 }
01371 }
01372
01373 akap = TMath::Abs(fDA[i]/fParamError[i]);
01374 if (akap > fAKAPPA) fAKAPPA=akap;
01375 }
01376 }
01377
01378 fGT = 0.;
01379 amb = 1.e18;
01380
01381 if (alambd > .0) amb = 0.25/alambd;
01382 for( i = 0; i < n; i++) {
01383 if (fPL[i] > .0) {
01384 if (nn2 > fNlimMul ) {
01385 if (TMath::Abs(fDA[i]/fPL[i]) > amb ) {
01386 fPL[i] = 4.*fPL[i];
01387 t1=4.;
01388 }
01389 }
01390
01391 fDA[i] = fDA[i]*alambd;
01392
01393 fGT = fGT + fDA[i]*fGr[i];
01394 }
01395 }
01396
01397
01398
01399
01400 if (-fGT <= sp && t1 < 1. && alambd < 1.)fENDFLG = -1;
01401
01402 if (fENDFLG >= 0) {
01403 if (fAKAPPA < TMath::Abs(fEPS)) {
01404 if (fixFLG == 0)
01405 fENDFLG=1;
01406 else {
01407 if (fENDFLG == 0) {
01408
01409 fENDFLG = 1;
01410 fixFLG = 0;
01411 ifix1=-1;
01412
01413 for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
01414 fINDFLG[1] = 0;
01415
01416 goto L19;
01417 } else {
01418 if( ifix1 >= 0) {
01419 fi = fi + 1;
01420 fENDFLG = 0;
01421 }
01422 }
01423 }
01424 } else {
01425 if( fixFLG != 0) {
01426 if( fi > fixFLG ) {
01427
01428 fENDFLG = 1;
01429 fixFLG = 0;
01430 ifix1=-1;
01431 for( i = 0; i < fNpar; i++) fPL[i] = fPL0[i];
01432 fINDFLG[1] = 0;
01433 goto L19;
01434 } else {
01435 fi = fi + 1;
01436 fENDFLG = 0;
01437 }
01438 } else {
01439 fi = fi + 1;
01440 fENDFLG = 0;
01441 }
01442 }
01443 }
01444
01445
01446
01447 if(fENDFLG == 0 && nn3 >= fNmaxIter) fENDFLG=-3;
01448
01449
01450 if(fENDFLG > 0 && fINDFLG[1] > 0) fENDFLG=-2;
01451
01452
01453 if (fENDFLG == 0) {
01454 for ( i = 0; i < n; i++) fA[i] = fA[i] + fDA[i];
01455 if (imax >= 0) fA[imax] = aiMAX;
01456 olds=fS;
01457 nn2=nn2+1;
01458 nn3=nn3+1;
01459 } else {
01460
01461
01462 Int_t il;
01463 il = 0;
01464 for( Int_t ip = 0; ip < fNpar; ip++) {
01465 if( fPL0[ip] > .0) {
01466 for( Int_t jp = 0; jp <= ip; jp++) {
01467 if(fPL0[jp] > .0) {
01468
01469 il = il + 1;
01470 }
01471 }
01472 }
01473 }
01474 return fENDFLG - 1;
01475 }
01476 goto L3;
01477 }
01478
01479
01480 void TFumili::PrintResults(Int_t ikode,Double_t p) const
01481 {
01482
01483
01484
01485
01486
01487
01488
01489
01490
01491
01492 TString exitStatus="";
01493 TString xsexpl="";
01494 TString colhdu[3],colhdl[3],cx2,cx3;
01495 switch (fENDFLG) {
01496 case 1:
01497 exitStatus="CONVERGED";
01498 break;
01499 case -1:
01500 exitStatus="CONST FCN";
01501 xsexpl="****\n* FUNCTION IS NOT DECREASING OR BAD DERIVATIVES\n****";
01502 break;
01503 case -2:
01504 exitStatus="ERRORS INF";
01505 xsexpl="****\n* ESTIMATED ERRORS ARE INfiNITE\n****";
01506 break;
01507 case -3:
01508 exitStatus="MAX ITER.";
01509 xsexpl="****\n* MAXIMUM NUMBER OF ITERATIONS IS EXCEEDED\n****";
01510 break;
01511 case -4:
01512 exitStatus="ZERO PROBAB";
01513 xsexpl="****\n* PROBABILITY OF LIKLIHOOD FUNCTION IS NEGATIVE OR ZERO\n****";
01514 break;
01515 default:
01516 exitStatus="UNDEfiNED";
01517 xsexpl="****\n* fiT IS IN PROGRESS\n****";
01518 break;
01519 }
01520 if (ikode == 1) {
01521 colhdu[0] = " ";
01522 colhdl[0] = " ERROR ";
01523 colhdu[1] = " PHYSICAL";
01524 colhdu[2] = " LIMITS ";
01525 colhdl[1] = " NEGATIVE ";
01526 colhdl[2] = " POSITIVE ";
01527 }
01528 if (ikode == 2) {
01529 colhdu[0] = " ";
01530 colhdl[0] = " ERROR ";
01531 colhdu[1] = " INTERNAL ";
01532 colhdl[1] = " STEP SIZE ";
01533 colhdu[2] = " INTERNAL ";
01534 colhdl[2] = " VALUE ";
01535 }
01536 if (ikode == 3) {
01537 colhdu[0] = " ";
01538 colhdl[0] = " ERROR ";
01539 colhdu[1] = " STEP ";
01540 colhdl[1] = " SIZE ";
01541 colhdu[2] = " fiRST ";
01542 colhdl[2] = " DERIVATIVE";
01543 }
01544 if (ikode == 4) {
01545 colhdu[0] = " PARABOLIC ";
01546 colhdl[0] = " ERROR ";
01547 colhdu[1] = " MINOS ";
01548 colhdu[2] = "ERRORS ";
01549 colhdl[1] = " NEGATIVE ";
01550 colhdl[2] = " POSITIVE ";
01551 }
01552 if(fENDFLG<1)Printf("%s", xsexpl.Data());
01553 Printf(" FCN=%g FROM FUMILI STATUS=%-10s %9d CALLS OF FCN",
01554 p,exitStatus.Data(),fNfcn);
01555 Printf(" EDM=%g ",-fGT);
01556 Printf(" EXT PARAMETER %-14s%-14s%-14s",
01557 (const char*)colhdu[0].Data()
01558 ,(const char*)colhdu[1].Data()
01559 ,(const char*)colhdu[2].Data());
01560 Printf(" NO. NAME VALUE %-14s%-14s%-14s",
01561 (const char*)colhdl[0].Data()
01562 ,(const char*)colhdl[1].Data()
01563 ,(const char*)colhdl[2].Data());
01564
01565 for (Int_t i=0;i<fNpar;i++) {
01566 if (ikode==3) {
01567 cx2 = Form("%14.5e",fDA[i]);
01568 cx3 = Form("%14.5e",fGr[i]);
01569
01570 }
01571 if (ikode==1) {
01572 cx2 = Form("%14.5e",fAMN[i]);
01573 cx3 = Form("%14.5e",fAMX[i]);
01574 }
01575 if (ikode==2) {
01576 cx2 = Form("%14.5e",fDA[i]);
01577 cx3 = Form("%14.5e",fA[i]);
01578 }
01579 if(ikode==4){
01580 cx2 = " *undefined* ";
01581 cx3 = " *undefined* ";
01582 }
01583 if(fPL0[i]<=0.) { cx2=" *fixed* ";cx3=""; }
01584 Printf("%4d %-11s%14.5e%14.5e%-14s%-14s",i+1
01585 ,fANames[i].Data(),fA[i],fParamError[i]
01586 ,cx2.Data(),cx3.Data());
01587 }
01588 }
01589
01590
01591
01592 void TFumili::ReleaseParameter(Int_t ipar) {
01593
01594
01595 if(ipar>=0 && ipar<fNpar && fPL0[ipar]<=0.) {
01596 fPL0[ipar] = -fPL0[ipar];
01597 if (fPL0[ipar] == 0. || fPL0[ipar]>=1.) fPL0[ipar]=.1;
01598 }
01599 }
01600
01601
01602
01603 void TFumili::SetData(Double_t *exdata,Int_t numpoints,Int_t vecsize){
01604
01605
01606
01607
01608
01609
01610
01611
01612
01613
01614
01615
01616
01617
01618
01619
01620
01621
01622
01623
01624
01625
01626
01627 if(exdata){
01628 fNED1 = numpoints;
01629 fNED2 = vecsize;
01630 fEXDA = exdata;
01631 }
01632 }
01633
01634
01635
01636 void TFumili::SetFitMethod(const char *name)
01637 {
01638
01639
01640 if (!strcmp(name,"H1FitChisquare")) SetFCN(H1FitChisquareFumili);
01641 if (!strcmp(name,"H1FitLikelihood")) SetFCN(H1FitLikelihoodFumili);
01642 if (!strcmp(name,"GraphFitChisquare")) SetFCN(GraphFitChisquareFumili);
01643 }
01644
01645
01646
01647 Int_t TFumili::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh) {
01648
01649
01650
01651
01652
01653 if (ipar<0 || ipar>=fNpar) return -1;
01654 fANames[ipar] = parname;
01655 fA[ipar] = value;
01656 fParamError[ipar] = verr;
01657 if(vlow<vhigh) {
01658 fAMN[ipar] = vlow;
01659 fAMX[ipar] = vhigh;
01660 } else {
01661 if(vhigh<vlow) {
01662 fAMN[ipar] = vhigh;
01663 fAMX[ipar] = vlow;
01664 }
01665 if(vhigh==vlow) {
01666 if(vhigh==0.) {
01667 ReleaseParameter(ipar);
01668 fAMN[ipar] = gMINDOUBLE;
01669 fAMX[ipar] = gMAXDOUBLE;
01670 }
01671 if(vlow!=0) FixParameter(ipar);
01672 }
01673 }
01674 return 0;
01675 }
01676
01677
01678 Int_t TFumili::SGZ()
01679 {
01680
01681
01682
01683 fS = 0.;
01684 Int_t i,j,l,k2=1,k1,ki=0;
01685 Double_t *x = new Double_t[fNED2];
01686 Double_t *df = new Double_t[fNpar];
01687 Int_t nx = fNED2-2;
01688 for (l=0;l<fNED1;l++) {
01689 k1 = k2;
01690 if (fLogLike) {
01691 fNumericDerivatives = kTRUE;
01692 nx = fNED2;
01693 k1 -= 2;
01694 }
01695
01696 for (i=0;i<nx;i++){
01697 ki += 1+i;
01698 x[i] = fEXDA[ki];
01699 }
01700
01701 Double_t y = EvalTFN(df,x);
01702 if(fNumericDerivatives) Derivatives(df,x);
01703 Double_t sig=1.;
01704 if(fLogLike) {
01705 if(y>0.) {
01706 fS = fS - log(y);
01707 y = -y;
01708 sig= y;
01709 } else {
01710 delete [] x;
01711 delete [] df;
01712 fS = 1e10;
01713 return -1;
01714 }
01715 } else {
01716 sig = fEXDA[k2];
01717 y = y - fEXDA[k1-1];
01718 fS = fS + (y*y/(sig*sig))*.5;
01719 }
01720 Int_t n = 0;
01721 for (i=0;i<fNpar;i++) {
01722 if (fPL0[i]>0){
01723 df[n] = df[i]/sig;
01724 fGr[i] += df[n]*(y/sig);
01725 n++;
01726 }
01727 }
01728 l = 0;
01729 for (i=0;i<n;i++)
01730 for (j=0;j<=i;j++)
01731 fZ[l++] += df[i]*df[j];
01732 k2 += fNED2;
01733 }
01734
01735 delete[] df;
01736 delete[] x;
01737 return 1;
01738 }
01739
01740
01741
01742 void TFumili::FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01743 {
01744
01745
01746
01747
01748
01749
01750
01751 Foption_t fitOption = GetFitOption();
01752 if (fitOption.Integral) {
01753 FitChisquareI(npar,gin,f,u,flag);
01754 return;
01755 }
01756 Double_t cu,eu,fu,fsum;
01757 Double_t x[3];
01758 Double_t *zik=0;
01759 Double_t *pl0=0;
01760
01761 TH1 *hfit = (TH1*)GetObjectFit();
01762 TF1 *f1 = (TF1*)GetUserFunc();
01763 Int_t nd = hfit->GetDimension();
01764 Int_t j;
01765
01766 npar = f1->GetNpar();
01767 SetParNumber(npar);
01768 if(flag == 9) return;
01769 zik = GetZ();
01770 pl0 = GetPL0();
01771
01772 Double_t *df = new Double_t[npar];
01773 f1->InitArgs(x,u);
01774 f = 0;
01775
01776 Int_t npfit = 0;
01777 Double_t *cache = fCache;
01778 for (Int_t i=0;i<fNpoints;i++) {
01779 if (nd > 2) x[2] = cache[4];
01780 if (nd > 1) x[1] = cache[3];
01781 x[0] = cache[2];
01782 cu = cache[0];
01783 TF1::RejectPoint(kFALSE);
01784 fu = f1->EvalPar(x,u);
01785 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01786 eu = cache[1];
01787 Derivatives(df,x);
01788 Int_t n = 0;
01789 fsum = (fu-cu)/eu;
01790 if (flag!=1) {
01791 for (j=0;j<npar;j++) {
01792 if (pl0[j]>0) {
01793 df[n] = df[j]/eu;
01794
01795 gin[j] += df[n]*fsum;
01796 n++;
01797 }
01798 }
01799 Int_t l = 0;
01800 for (j=0;j<n;j++)
01801 for (Int_t k=0;k<=j;k++)
01802 zik[l++] += df[j]*df[k];
01803 }
01804 f += .5*fsum*fsum;
01805 npfit++;
01806 cache += fPointSize;
01807 }
01808 f1->SetNumberFitPoints(npfit);
01809 delete [] df;
01810 }
01811
01812
01813
01814 void TFumili::FitChisquareI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01815 {
01816
01817
01818
01819
01820
01821
01822
01823 Double_t cu,eu,fu,fsum;
01824 Double_t x[3];
01825 Double_t *zik=0;
01826 Double_t *pl0=0;
01827
01828 TH1 *hfit = (TH1*)GetObjectFit();
01829 TF1 *f1 = (TF1*)GetUserFunc();
01830 Int_t nd = hfit->GetDimension();
01831 Int_t j;
01832
01833 f1->InitArgs(x,u);
01834 npar = f1->GetNpar();
01835 SetParNumber(npar);
01836 if(flag == 9) return;
01837 zik = GetZ();
01838 pl0 = GetPL0();
01839
01840 Double_t *df=new Double_t[npar];
01841 f = 0;
01842
01843 Int_t npfit = 0;
01844 Double_t *cache = fCache;
01845 for (Int_t i=0;i<fNpoints;i++) {
01846 cu = cache[0];
01847 TF1::RejectPoint(kFALSE);
01848 f1->SetParameters(u);
01849 if (nd < 2) {
01850 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
01851 } else if (nd < 3) {
01852 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
01853 } else {
01854 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
01855 }
01856 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01857 eu = cache[1];
01858 Derivatives(df,x);
01859 Int_t n = 0;
01860 fsum = (fu-cu)/eu;
01861 if (flag!=1) {
01862 for (j=0;j<npar;j++) {
01863 if (pl0[j]>0){
01864 df[n] = df[j]/eu;
01865
01866 gin[j] += df[n]*fsum;
01867 n++;
01868 }
01869 }
01870 Int_t l = 0;
01871 for (j=0;j<n;j++)
01872 for (Int_t k=0;k<=j;k++)
01873 zik[l++] += df[j]*df[k];
01874 }
01875 f += .5*fsum*fsum;
01876 npfit++;
01877 cache += fPointSize;
01878 }
01879 f1->SetNumberFitPoints(npfit);
01880 delete[] df;
01881 }
01882
01883
01884
01885 void TFumili::FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01886 {
01887
01888
01889
01890
01891
01892
01893
01894
01895
01896
01897
01898 Foption_t fitOption = GetFitOption();
01899 if (fitOption.Integral) {
01900 FitLikelihoodI(npar,gin,f,u,flag);
01901 return;
01902 }
01903 Double_t cu,fu,fobs,fsub;
01904 Double_t dersum[100];
01905 Double_t x[3];
01906 Int_t icu;
01907
01908 TH1 *hfit = (TH1*)GetObjectFit();
01909 TF1 *f1 = (TF1*)GetUserFunc();
01910 Int_t nd = hfit->GetDimension();
01911 Int_t j;
01912 Double_t *zik = GetZ();
01913 Double_t *pl0 = GetPL0();
01914
01915 npar = f1->GetNpar();
01916 SetParNumber(npar);
01917 if(flag == 9) return;
01918 Double_t *df=new Double_t[npar];
01919 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
01920 f1->InitArgs(x,u);
01921 f = 0;
01922
01923 Int_t npfit = 0;
01924 Double_t *cache = fCache;
01925 for (Int_t i=0;i<fNpoints;i++) {
01926 if (nd > 2) x[2] = cache[4];
01927 if (nd > 1) x[1] = cache[3];
01928 x[0] = cache[2];
01929 cu = cache[0];
01930 TF1::RejectPoint(kFALSE);
01931 fu = f1->EvalPar(x,u);
01932 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
01933 if (flag == 2) {
01934 for (j=0;j<npar;j++) {
01935 dersum[j] += 1;
01936
01937 }
01938 }
01939 if (fu < 1.e-9) fu = 1.e-9;
01940 icu = Int_t(cu);
01941 fsub = -fu +icu*TMath::Log(fu);
01942 fobs = GetSumLog(icu);
01943 fsub -= fobs;
01944 Derivatives(df,x);
01945 int n=0;
01946
01947
01948 for (j=0;j<npar;j++) {
01949 if (pl0[j]>0){
01950 df[n] = df[j]*(icu/fu-1);
01951 gin[j] -= df[n];
01952 n++;
01953 }
01954 }
01955 Int_t l = 0;
01956
01957
01958 for (j=0;j<n;j++)
01959 for (Int_t k=0;k<=j;k++)
01960 zik[l++] += df[j]*df[k];
01961
01962 f -= fsub;
01963 npfit++;
01964 cache += fPointSize;
01965 }
01966 f *= 2;
01967 f1->SetNumberFitPoints(npfit);
01968 delete[] df;
01969 }
01970
01971
01972
01973 void TFumili::FitLikelihoodI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
01974 {
01975
01976
01977
01978
01979
01980
01981
01982
01983
01984
01985
01986 Double_t cu,fu,fobs,fsub;
01987 Double_t dersum[100];
01988 Double_t x[3];
01989 Int_t icu;
01990
01991 TH1 *hfit = (TH1*)GetObjectFit();
01992 TF1 *f1 = (TF1*)GetUserFunc();
01993 Int_t nd = hfit->GetDimension();
01994 Int_t j;
01995 Double_t *zik = GetZ();
01996 Double_t *pl0 = GetPL0();
01997
01998 Double_t *df=new Double_t[npar];
01999
02000 npar = f1->GetNpar();
02001 SetParNumber(npar);
02002 if(flag == 9) {delete [] df; return;}
02003 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
02004 f1->InitArgs(x,u);
02005 f = 0;
02006
02007 Int_t npfit = 0;
02008 Double_t *cache = fCache;
02009 for (Int_t i=0;i<fNpoints;i++) {
02010 if (nd > 2) x[2] = cache[4];
02011 if (nd > 1) x[1] = cache[3];
02012 x[0] = cache[2];
02013 cu = cache[0];
02014 TF1::RejectPoint(kFALSE);
02015 if (nd < 2) {
02016 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
02017 } else if (nd < 3) {
02018 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5])/(cache[3]*cache[5]);
02019 } else {
02020 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],cache[4] - 0.5*cache[5],cache[4] + 0.5*cache[5],cache[6] - 0.5*cache[7],cache[6] + 0.5*cache[7])/(cache[3]*cache[5]*cache[7]);
02021 }
02022 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
02023 if (flag == 2) {
02024 for (j=0;j<npar;j++) {
02025 dersum[j] += 1;
02026
02027 }
02028 }
02029 if (fu < 1.e-9) fu = 1.e-9;
02030 icu = Int_t(cu);
02031 fsub = -fu +icu*TMath::Log(fu);
02032 fobs = GetSumLog(icu);
02033 fsub -= fobs;
02034 Derivatives(df,x);
02035 int n=0;
02036
02037
02038 for (j=0;j<npar;j++) {
02039 if (pl0[j]>0){
02040 df[n] = df[j]*(icu/fu-1);
02041 gin[j] -= df[n];
02042 n++;
02043 }
02044 }
02045 Int_t l = 0;
02046
02047
02048 for (j=0;j<n;j++)
02049 for (Int_t k=0;k<=j;k++)
02050 zik[l++] += df[j]*df[k];
02051
02052 f -= fsub;
02053 npfit++;
02054 cache += fPointSize;
02055 }
02056 f *= 2;
02057 f1->SetNumberFitPoints(npfit);
02058 delete[] df;
02059 }
02060
02061
02062
02063
02064
02065
02066
02067
02068 void H1FitChisquareFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
02069 {
02070
02071
02072
02073 TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
02074 hFitter->FitChisquare(npar, gin, f, u, flag);
02075 }
02076
02077
02078 void H1FitLikelihoodFumili(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
02079 {
02080
02081
02082
02083
02084
02085
02086
02087
02088
02089
02090
02091
02092
02093 TFumili *hFitter = (TFumili*)TVirtualFitter::GetFitter();
02094 hFitter->FitLikelihood(npar, gin, f, u, flag);
02095 }
02096
02097
02098
02099
02100 void GraphFitChisquareFumili(Int_t &npar, Double_t * gin, Double_t &f,
02101 Double_t *u, Int_t flag)
02102 {
02103
02104
02105
02106
02107
02108
02109
02110
02111
02112
02113
02114
02115
02116
02117
02118
02119
02120
02121
02122
02123
02124
02125
02126
02127
02128
02129
02130
02131 Double_t cu,eu,exl,exh,ey,eux,fu,fsum;
02132 Double_t x[1];
02133 Int_t i, bin, npfits=0;
02134
02135 TFumili *grFitter = (TFumili*)TVirtualFitter::GetFitter();
02136 TGraph *gr = (TGraph*)grFitter->GetObjectFit();
02137 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
02138 Foption_t fitOption = grFitter->GetFitOption();
02139
02140 Int_t n = gr->GetN();
02141 Double_t *gx = gr->GetX();
02142 Double_t *gy = gr->GetY();
02143 npar = f1->GetNpar();
02144
02145 grFitter->SetParNumber(npar);
02146
02147 if(flag == 9) return;
02148 Double_t *zik = grFitter->GetZ();
02149 Double_t *pl0 = grFitter->GetPL0();
02150 Double_t *df = new Double_t[npar];
02151
02152
02153 f1->InitArgs(x,u);
02154 f = 0;
02155 for (bin=0;bin<n;bin++) {
02156 x[0] = gx[bin];
02157 if (!f1->IsInside(x)) continue;
02158 cu = gy[bin];
02159 TF1::RejectPoint(kFALSE);
02160 fu = f1->EvalPar(x,u);
02161 if (TF1::RejectedPoint()) continue;
02162 npfits++;
02163 Double_t eusq=1.;
02164 if (fitOption.W1) {
02165 eu = 1.;
02166 } else {
02167 exh = gr->GetErrorXhigh(bin);
02168 exl = gr->GetErrorXlow(bin);
02169 ey = gr->GetErrorY(bin);
02170 if (exl < 0) exl = 0;
02171 if (exh < 0) exh = 0;
02172 if (ey < 0) ey = 0;
02173 if (exh > 0 && exl > 0) {
02174
02175 eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
02176 } else
02177 eux = 0.;
02178 eu = ey*ey+eux*eux;
02179 if (eu <= 0) eu = 1;
02180 eusq = TMath::Sqrt(eu);
02181 }
02182 grFitter->Derivatives(df,x);
02183 n = 0;
02184 fsum = (fu-cu)/eusq;
02185 for (i=0;i<npar;i++) {
02186 if (pl0[i]>0){
02187 df[n] = df[i]/eusq;
02188
02189 gin[i] += df[n]*fsum;
02190 n++;
02191 }
02192 }
02193 Int_t l = 0;
02194 for (i=0;i<n;i++)
02195 for (Int_t j=0;j<=i;j++)
02196 zik[l++] += df[i]*df[j];
02197 f += .5*fsum*fsum;
02198
02199 }
02200 delete [] df;
02201 f1->SetNumberFitPoints(npfits);
02202 }
02203