00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include "TMinuit.h"
00012 #include "TFitter.h"
00013 #include "TH1.h"
00014 #include "TF1.h"
00015 #include "TF2.h"
00016 #include "TF3.h"
00017 #include "TList.h"
00018 #include "TGraph.h"
00019 #include "TGraph2D.h"
00020 #include "TMultiGraph.h"
00021 #include "TMath.h"
00022
00023 extern void H1FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00024 extern void H1FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00025 extern void GraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00026 extern void Graph2DFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00027 extern void MultiGraphFitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00028 extern void F2Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00029 extern void F3Fit(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00030
00031 ClassImp(TFitter)
00032
00033
00034 TFitter::TFitter(Int_t maxpar)
00035 {
00036
00037
00038
00039 fMinuit = new TMinuit(maxpar);
00040 fNlog = 0;
00041 fSumLog = 0;
00042 fCovar = 0;
00043 SetName("MinuitFitter");
00044 }
00045
00046
00047 TFitter::~TFitter()
00048 {
00049
00050
00051
00052 if (fCovar) delete [] fCovar;
00053 if (fSumLog) delete [] fSumLog;
00054 delete fMinuit;
00055 }
00056
00057
00058 Double_t TFitter::Chisquare(Int_t npar, Double_t *params) const
00059 {
00060
00061
00062 Double_t amin = 0;
00063 H1FitChisquare(npar,params,amin,params,1);
00064 return amin;
00065 }
00066
00067
00068 void TFitter::Clear(Option_t *)
00069 {
00070
00071
00072 if (fCovar) {delete [] fCovar; fCovar = 0;}
00073 fMinuit->mncler();
00074
00075
00076 Double_t val = 3;
00077 Int_t inseed = 12345;
00078 fMinuit->mnrn15(val,inseed);
00079 }
00080
00081
00082 Int_t TFitter::ExecuteCommand(const char *command, Double_t *args, Int_t nargs)
00083 {
00084
00085
00086
00087
00088 if (fCovar) {delete [] fCovar; fCovar = 0;}
00089 Int_t ierr = 0;
00090 fMinuit->mnexcm(command,args,nargs,ierr);
00091 return ierr;
00092 }
00093
00094
00095 void TFitter::FixParameter(Int_t ipar)
00096 {
00097
00098
00099 if (fCovar) {delete [] fCovar; fCovar = 0;}
00100 fMinuit->FixParameter(ipar);
00101 }
00102
00103
00104 void TFitter::GetConfidenceIntervals(Int_t n, Int_t ndim, const Double_t *x, Double_t *ci, Double_t cl)
00105 {
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116 TF1 *f = (TF1*)fUserFunc;
00117 Int_t npar = f->GetNumberFreeParameters();
00118 Int_t npar_real = f->GetNpar();
00119 Double_t *grad = new Double_t[npar_real];
00120 Double_t *sum_vector = new Double_t[npar];
00121 Bool_t *fixed=0;
00122 Double_t al, bl;
00123 if (npar_real != npar){
00124 fixed = new Bool_t[npar_real];
00125 memset(fixed,0,npar_real*sizeof(Bool_t));
00126
00127 for (Int_t ipar=0; ipar<npar_real; ipar++){
00128 fixed[ipar]=0;
00129 f->GetParLimits(ipar,al,bl);
00130 if (al*bl != 0 && al >= bl) {
00131
00132 fixed[ipar]=1;
00133 }
00134 }
00135 }
00136 Double_t c=0;
00137
00138 Double_t *matr = GetCovarianceMatrix();
00139 if (!matr)
00140 return;
00141 Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
00142 Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
00143 Int_t igrad, ifree=0;
00144 for (Int_t ipoint=0; ipoint<n; ipoint++){
00145 c=0;
00146 f->GradientPar(x+ndim*ipoint, grad);
00147
00148 for (Int_t irow=0; irow<npar; irow++){
00149 sum_vector[irow]=0;
00150 igrad = 0;
00151 for (Int_t icol=0; icol<npar; icol++){
00152 igrad = 0;
00153 ifree=0;
00154 if (fixed) {
00155
00156 while (ifree<icol+1){
00157 if (fixed[igrad]==0) ifree++;
00158 igrad++;
00159 }
00160 igrad--;
00161
00162 } else {
00163 igrad = icol;
00164 }
00165 sum_vector[irow]+=matr[irow*npar_real+icol]*grad[igrad];
00166 }
00167 }
00168 igrad = 0;
00169 for (Int_t i=0; i<npar; i++){
00170 igrad = 0; ifree=0;
00171 if (fixed) {
00172
00173 while (ifree<i+1){
00174 if (fixed[igrad]==0) ifree++;
00175 igrad++;
00176 }
00177 igrad--;
00178 } else {
00179 igrad = i;
00180 }
00181 c+=grad[igrad]*sum_vector[i];
00182 }
00183
00184 c=TMath::Sqrt(c);
00185 ci[ipoint]=c*t*chidf;
00186 }
00187
00188 delete [] grad;
00189 delete [] sum_vector;
00190 if (fixed)
00191 delete [] fixed;
00192 }
00193
00194
00195 void TFitter::GetConfidenceIntervals(TObject *obj, Double_t cl)
00196 {
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218 if (obj->InheritsFrom(TGraph::Class())) {
00219 TGraph *gr = (TGraph*)obj;
00220 if (!gr->GetEY()){
00221 Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a graph");
00222 return;
00223 }
00224 if (fObjectFit->InheritsFrom(TGraph2D::Class())){
00225 Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a graph");
00226 return;
00227 }
00228 if (fObjectFit->InheritsFrom(TH1::Class())){
00229 if (((TH1*)(fObjectFit))->GetDimension()>1){
00230 Error("GetConfidenceIntervals", "A TGraph2DErrors or a TH23 should be passed instead of a graph");
00231 return;
00232 }
00233 }
00234 GetConfidenceIntervals(gr->GetN(),1,gr->GetX(), gr->GetEY(), cl);
00235 for (Int_t i=0; i<gr->GetN(); i++)
00236 gr->SetPoint(i, gr->GetX()[i], ((TF1*)(fUserFunc))->Eval(gr->GetX()[i]));
00237 }
00238
00239
00240 else if (obj->InheritsFrom(TGraph2D::Class())) {
00241 TGraph2D *gr2 = (TGraph2D*)obj;
00242 if (!gr2->GetEZ()){
00243 Error("GetConfidenceIntervals", "A TGraph2DErrors should be passed instead of a TGraph2D");
00244 return;
00245 }
00246 if (fObjectFit->InheritsFrom(TGraph::Class())){
00247 Error("GetConfidenceIntervals", "A TGraphErrors should be passed instead of a TGraph2D");
00248 return;
00249 }
00250 if (fObjectFit->InheritsFrom(TH1::Class())){
00251 if (((TH1*)(fObjectFit))->GetDimension()==1){
00252 Error("GetConfidenceIntervals", "A TGraphErrors or a TH1 should be passed instead of a graph");
00253 return;
00254 }
00255 }
00256 TF2 *f=(TF2*)fUserFunc;
00257 Double_t xy[2];
00258 Int_t np = gr2->GetN();
00259 Int_t npar = f->GetNpar();
00260 Double_t *grad = new Double_t[npar];
00261 Double_t *sum_vector = new Double_t[npar];
00262 Double_t *x = gr2->GetX();
00263 Double_t *y = gr2->GetY();
00264 Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
00265 Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
00266 Double_t *matr=GetCovarianceMatrix();
00267 Double_t c = 0;
00268 for (Int_t ipoint=0; ipoint<np; ipoint++){
00269 xy[0]=x[ipoint];
00270 xy[1]=y[ipoint];
00271 f->GradientPar(xy, grad);
00272 for (Int_t irow=0; irow<f->GetNpar(); irow++){
00273 sum_vector[irow]=0;
00274 for (Int_t icol=0; icol<npar; icol++)
00275 sum_vector[irow]+=matr[irow*npar+icol]*grad[icol];
00276 }
00277 c = 0;
00278 for (Int_t i=0; i<npar; i++)
00279 c+=grad[i]*sum_vector[i];
00280 c=TMath::Sqrt(c);
00281 gr2->SetPoint(ipoint, xy[0], xy[1], f->EvalPar(xy));
00282 gr2->GetEZ()[ipoint]=c*t*chidf;
00283
00284 }
00285 delete [] grad;
00286 delete [] sum_vector;
00287 }
00288
00289
00290 else if (obj->InheritsFrom(TH1::Class())) {
00291 if (fObjectFit->InheritsFrom(TGraph::Class())){
00292 if (((TH1*)obj)->GetDimension()>1){
00293 Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
00294 return;
00295 }
00296 }
00297 if (fObjectFit->InheritsFrom(TGraph2D::Class())){
00298 if (((TH1*)obj)->GetDimension()!=2){
00299 Error("GetConfidenceIntervals", "Fitted graph and passed histogram have different number of dimensions");
00300 return;
00301 }
00302 }
00303 if (fObjectFit->InheritsFrom(TH1::Class())){
00304 if (((TH1*)(fObjectFit))->GetDimension()!=((TH1*)(obj))->GetDimension()){
00305 Error("GetConfidenceIntervals", "Fitted and passed histograms have different number of dimensions");
00306 return;
00307 }
00308 }
00309
00310
00311 TH1 *hfit = (TH1*)obj;
00312 TF1 *f = (TF1*)GetUserFunc();
00313 Int_t npar = f->GetNpar();
00314 Double_t *grad = new Double_t[npar];
00315 Double_t *sum_vector = new Double_t[npar];
00316 Double_t x[3];
00317
00318 Int_t hxfirst = hfit->GetXaxis()->GetFirst();
00319 Int_t hxlast = hfit->GetXaxis()->GetLast();
00320 Int_t hyfirst = hfit->GetYaxis()->GetFirst();
00321 Int_t hylast = hfit->GetYaxis()->GetLast();
00322 Int_t hzfirst = hfit->GetZaxis()->GetFirst();
00323 Int_t hzlast = hfit->GetZaxis()->GetLast();
00324
00325 TAxis *xaxis = hfit->GetXaxis();
00326 TAxis *yaxis = hfit->GetYaxis();
00327 TAxis *zaxis = hfit->GetZaxis();
00328 Double_t t = TMath::StudentQuantile(0.5 + cl/2, f->GetNDF());
00329 Double_t chidf = TMath::Sqrt(f->GetChisquare()/f->GetNDF());
00330 Double_t *matr=GetCovarianceMatrix();
00331 Double_t c=0;
00332 for (Int_t binz=hzfirst; binz<=hzlast; binz++){
00333 x[2]=zaxis->GetBinCenter(binz);
00334 for (Int_t biny=hyfirst; biny<=hylast; biny++) {
00335 x[1]=yaxis->GetBinCenter(biny);
00336 for (Int_t binx=hxfirst; binx<=hxlast; binx++) {
00337 x[0]=xaxis->GetBinCenter(binx);
00338 f->GradientPar(x, grad);
00339 for (Int_t irow=0; irow<npar; irow++){
00340 sum_vector[irow]=0;
00341 for (Int_t icol=0; icol<npar; icol++)
00342 sum_vector[irow]+=matr[irow*npar+icol]*grad[icol];
00343 }
00344 c = 0;
00345 for (Int_t i=0; i<npar; i++)
00346 c+=grad[i]*sum_vector[i];
00347 c=TMath::Sqrt(c);
00348 hfit->SetBinContent(binx, biny, binz, f->EvalPar(x));
00349 hfit->SetBinError(binx, biny, binz, c*t*chidf);
00350 }
00351 }
00352 }
00353 delete [] grad;
00354 delete [] sum_vector;
00355 }
00356 else {
00357 Error("GetConfidenceIntervals", "This object type is not supported");
00358 return;
00359 }
00360
00361 }
00362
00363
00364 Double_t *TFitter::GetCovarianceMatrix() const
00365 {
00366
00367
00368 if (fCovar) return fCovar;
00369 Int_t npars = fMinuit->GetNumPars();
00370 ((TFitter*)this)->fCovar = new Double_t[npars*npars];
00371 fMinuit->mnemat(fCovar,npars);
00372 return fCovar;
00373 }
00374
00375
00376 Double_t TFitter::GetCovarianceMatrixElement(Int_t i, Int_t j) const
00377 {
00378
00379
00380 GetCovarianceMatrix();
00381 Int_t npars = fMinuit->GetNumPars();
00382 if (i < 0 || i >= npars || j < 0 || j >= npars) {
00383 Error("GetCovarianceMatrixElement","Illegal arguments i=%d, j=%d",i,j);
00384 return 0;
00385 }
00386 return fCovar[j+npars*i];
00387 }
00388
00389
00390 Int_t TFitter::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const
00391 {
00392
00393
00394
00395
00396
00397
00398
00399
00400 Int_t ierr = 0;
00401 fMinuit->mnerrs(ipar, eplus,eminus,eparab,globcc);
00402 return ierr;
00403 }
00404
00405
00406
00407 Int_t TFitter::GetNumberTotalParameters() const
00408 {
00409
00410
00411 return fMinuit->fNpar + fMinuit->fNpfix;
00412 }
00413
00414
00415 Int_t TFitter::GetNumberFreeParameters() const
00416 {
00417
00418
00419 return fMinuit->fNpar;
00420 }
00421
00422
00423
00424 Double_t TFitter::GetParError(Int_t ipar) const
00425 {
00426
00427
00428 Int_t ierr = 0;
00429 TString pname;
00430 Double_t value,verr,vlow,vhigh;
00431
00432 fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
00433 return verr;
00434 }
00435
00436
00437
00438 Double_t TFitter::GetParameter(Int_t ipar) const
00439 {
00440
00441
00442 Int_t ierr = 0;
00443 TString pname;
00444 Double_t value,verr,vlow,vhigh;
00445
00446 fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
00447 return value;
00448 }
00449
00450
00451 Int_t TFitter::GetParameter(Int_t ipar, char *parname,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const
00452 {
00453
00454
00455
00456
00457
00458
00459
00460
00461
00462 Int_t ierr = 0;
00463 TString pname;
00464 fMinuit->mnpout(ipar, pname,value,verr,vlow,vhigh,ierr);
00465 strcpy(parname,pname.Data());
00466 return ierr;
00467 }
00468
00469
00470 const char *TFitter::GetParName(Int_t ipar) const
00471 {
00472
00473
00474 if (!fMinuit || ipar < 0 || ipar > fMinuit->fNu) return "";
00475 return fMinuit->fCpnam[ipar];
00476 }
00477
00478
00479 Int_t TFitter::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const
00480 {
00481
00482
00483
00484
00485
00486
00487
00488 Int_t ierr = 0;
00489 fMinuit->mnstat(amin,edm,errdef,nvpar,nparx,ierr);
00490 return ierr;
00491 }
00492
00493
00494 Double_t TFitter::GetSumLog(Int_t n)
00495 {
00496
00497
00498
00499 if (n < 0) return 0;
00500 if (n > fNlog) {
00501 if (fSumLog) delete [] fSumLog;
00502 fNlog = 2*n+1000;
00503 fSumLog = new Double_t[fNlog+1];
00504 Double_t fobs = 0;
00505 for (Int_t j=0;j<=fNlog;j++) {
00506 if (j > 1) fobs += TMath::Log(j);
00507 fSumLog[j] = fobs;
00508 }
00509 }
00510 if (fSumLog) return fSumLog[n];
00511 return 0;
00512 }
00513
00514
00515
00516 Bool_t TFitter::IsFixed(Int_t ipar) const
00517 {
00518
00519
00520 if (fMinuit->fNiofex[ipar] == 0 ) return kTRUE;
00521 return kFALSE;
00522 }
00523
00524
00525
00526 void TFitter::PrintResults(Int_t level, Double_t amin) const
00527 {
00528
00529
00530 fMinuit->mnprin(level,amin);
00531 }
00532
00533
00534 void TFitter::ReleaseParameter(Int_t ipar)
00535 {
00536
00537
00538 if (fCovar) {delete [] fCovar; fCovar = 0;}
00539 fMinuit->Release(ipar);
00540 }
00541
00542
00543 void TFitter::SetFCN(void *fcn)
00544 {
00545
00546
00547 if (fCovar) {delete [] fCovar; fCovar = 0;}
00548 TVirtualFitter::SetFCN(fcn);
00549 fMinuit->SetFCN(fcn);
00550
00551 }
00552
00553
00554 void TFitter::SetFCN(void (*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
00555 {
00556
00557
00558 if (fCovar) {delete [] fCovar; fCovar = 0;}
00559 TVirtualFitter::SetFCN(fcn);
00560 fMinuit->SetFCN(fcn);
00561 }
00562
00563
00564 void TFitter::SetFitMethod(const char *name)
00565 {
00566
00567
00568 if (fCovar) {delete [] fCovar; fCovar = 0;}
00569 if (!strcmp(name,"H1FitChisquare")) SetFCN(H1FitChisquare);
00570 if (!strcmp(name,"H1FitLikelihood")) SetFCN(H1FitLikelihood);
00571 if (!strcmp(name,"GraphFitChisquare")) SetFCN(GraphFitChisquare);
00572 if (!strcmp(name,"Graph2DFitChisquare")) SetFCN(Graph2DFitChisquare);
00573 if (!strcmp(name,"MultiGraphFitChisquare")) SetFCN(MultiGraphFitChisquare);
00574 if (!strcmp(name,"F2Minimizer")) SetFCN(F2Fit);
00575 if (!strcmp(name,"F3Minimizer")) SetFCN(F3Fit);
00576 }
00577
00578
00579 Int_t TFitter::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh)
00580 {
00581
00582
00583
00584
00585
00586
00587
00588
00589 if (fCovar) {delete [] fCovar; fCovar = 0;}
00590 Int_t ierr = 0;
00591 fMinuit->mnparm(ipar,parname,value,verr,vlow,vhigh,ierr);
00592 return ierr;
00593 }
00594
00595
00596 void TFitter::FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00597 {
00598
00599
00600
00601
00602
00603
00604
00605 Foption_t fitOption = GetFitOption();
00606 if (fitOption.Integral) {
00607 FitChisquareI(npar,gin,f,u,flag);
00608 return;
00609 }
00610 Double_t cu,eu,fu,fsum;
00611 Double_t dersum[100], grad[100];
00612 memset(grad,0,800);
00613 Double_t x[3];
00614
00615 TH1 *hfit = (TH1*)GetObjectFit();
00616 TF1 *f1 = (TF1*)GetUserFunc();
00617 Int_t nd = hfit->GetDimension();
00618 Int_t j;
00619
00620 f1->InitArgs(x,u);
00621 npar = f1->GetNpar();
00622 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
00623 f = 0;
00624
00625 Int_t npfit = 0;
00626 Double_t *cache = fCache;
00627 for (Int_t i=0;i<fNpoints;i++) {
00628 if (nd > 2) x[2] = cache[4];
00629 if (nd > 1) x[1] = cache[3];
00630 x[0] = cache[2];
00631 cu = cache[0];
00632 TF1::RejectPoint(kFALSE);
00633 fu = f1->EvalPar(x,u);
00634 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
00635 eu = cache[1];
00636 if (flag == 2) {
00637 for (j=0;j<npar;j++) dersum[j] += 1;
00638 for (j=0;j<npar;j++) grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
00639 }
00640 fsum = (cu-fu)/eu;
00641 f += fsum*fsum;
00642 npfit++;
00643 cache += fPointSize;
00644 }
00645 f1->SetNumberFitPoints(npfit);
00646 }
00647
00648
00649 void TFitter::FitChisquareI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00650 {
00651
00652
00653
00654
00655
00656
00657
00658 Double_t cu,eu,fu,fsum;
00659 Double_t dersum[100], grad[100];
00660 memset(grad,0,800);
00661 Double_t x[3];
00662
00663 TH1 *hfit = (TH1*)GetObjectFit();
00664 TF1 *f1 = (TF1*)GetUserFunc();
00665 Int_t nd = hfit->GetDimension();
00666 Int_t j;
00667
00668 f1->InitArgs(x,u);
00669 npar = f1->GetNpar();
00670 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
00671 f = 0;
00672
00673 Int_t npfit = 0;
00674 Double_t *cache = fCache;
00675 for (Int_t i=0;i<fNpoints;i++) {
00676 cu = cache[0];
00677 TF1::RejectPoint(kFALSE);
00678 f1->SetParameters(u);
00679 if (nd < 2) {
00680 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
00681 } else if (nd < 3) {
00682 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]);
00683 } else {
00684 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]);
00685 }
00686 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
00687 eu = cache[1];
00688 if (flag == 2) {
00689 for (j=0;j<npar;j++) dersum[j] += 1;
00690 for (j=0;j<npar;j++) grad[j] += dersum[j]*(fu-cu)/eu; dersum[j] = 0;
00691 }
00692 fsum = (cu-fu)/eu;
00693 f += fsum*fsum;
00694 npfit++;
00695 cache += fPointSize;
00696 }
00697 f1->SetNumberFitPoints(npfit);
00698 }
00699
00700
00701
00702 void TFitter::FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00703 {
00704
00705
00706
00707
00708
00709
00710
00711
00712
00713
00714
00715 Foption_t fitOption = GetFitOption();
00716 if (fitOption.Integral) {
00717 FitLikelihoodI(npar,gin,f,u,flag);
00718 return;
00719 }
00720 Double_t cu,fu,fobs,fsub;
00721 Double_t dersum[100];
00722 Double_t x[3];
00723 Int_t icu;
00724
00725 TH1 *hfit = (TH1*)GetObjectFit();
00726 TF1 *f1 = (TF1*)GetUserFunc();
00727 Int_t nd = hfit->GetDimension();
00728 Int_t j;
00729
00730 f1->InitArgs(x,u);
00731 npar = f1->GetNpar();
00732 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
00733 f = 0;
00734
00735 Int_t npfit = 0;
00736 Double_t *cache = fCache;
00737 for (Int_t i=0;i<fNpoints;i++) {
00738 if (nd > 2) x[2] = cache[4];
00739 if (nd > 1) x[1] = cache[3];
00740 x[0] = cache[2];
00741 cu = cache[0];
00742 TF1::RejectPoint(kFALSE);
00743 fu = f1->EvalPar(x,u);
00744 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
00745 if (flag == 2) {
00746 for (j=0;j<npar;j++) {
00747 dersum[j] += 1;
00748
00749 }
00750 }
00751 if (fu < 1.e-9) fu = 1.e-9;
00752 if (fitOption.Like == 1) {
00753 icu = Int_t(cu);
00754 fsub = -fu +cu*TMath::Log(fu);
00755 if (icu <10000) fobs = GetSumLog(icu);
00756 else fobs = TMath::LnGamma(cu+1);
00757 } else {
00758 fsub = -fu +cu*TMath::Log(fu);
00759 fobs = TMath::LnGamma(cu+1);
00760 }
00761 fsub -= fobs;
00762 f -= fsub;
00763 npfit++;
00764 cache += fPointSize;
00765 }
00766 f *= 2;
00767 f1->SetNumberFitPoints(npfit);
00768 }
00769
00770
00771
00772 void TFitter::FitLikelihoodI(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00773 {
00774
00775
00776
00777
00778
00779
00780
00781
00782
00783
00784
00785 Double_t cu,fu,fobs,fsub;
00786 Double_t dersum[100];
00787 Double_t x[3];
00788 Int_t icu;
00789
00790 TH1 *hfit = (TH1*)GetObjectFit();
00791 TF1 *f1 = (TF1*)GetUserFunc();
00792 Foption_t fitOption = GetFitOption();
00793 Int_t nd = hfit->GetDimension();
00794 Int_t j;
00795
00796 f1->InitArgs(x,u);
00797 npar = f1->GetNpar();
00798 if (flag == 2) for (j=0;j<npar;j++) dersum[j] = gin[j] = 0;
00799 f = 0;
00800
00801 Int_t npfit = 0;
00802 Double_t *cache = fCache;
00803 for (Int_t i=0;i<fNpoints;i++) {
00804 if (nd > 2) x[2] = cache[6];
00805 if (nd > 1) x[1] = cache[4];
00806 x[0] = cache[2];
00807 cu = cache[0];
00808 TF1::RejectPoint(kFALSE);
00809 f1->SetParameters(u);
00810 if (nd < 2) {
00811 fu = f1->Integral(cache[2] - 0.5*cache[3],cache[2] + 0.5*cache[3],u)/cache[3];
00812 } else if (nd < 3) {
00813 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]);
00814 } else {
00815 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]);
00816 }
00817 if (TF1::RejectedPoint()) {cache += fPointSize; continue;}
00818 if (flag == 2) {
00819 for (j=0;j<npar;j++) {
00820 dersum[j] += 1;
00821
00822 }
00823 }
00824 if (fu < 1.e-9) fu = 1.e-9;
00825 if (fitOption.Like == 1) {
00826 icu = Int_t(cu);
00827 fsub = -fu +cu*TMath::Log(fu);
00828 if (icu <10000) fobs = GetSumLog(icu);
00829 else fobs = TMath::LnGamma(cu+1);
00830 } else {
00831 fsub = -fu +cu*TMath::Log(fu);
00832 fobs = TMath::LnGamma(cu+1);
00833 }
00834 fsub -= fobs;
00835 f -= fsub;
00836 npfit++;
00837 cache += fPointSize;
00838 }
00839 f *= 2;
00840 f1->SetNumberFitPoints(npfit);
00841 }
00842
00843
00844
00845
00846 void H1FitChisquare(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00847 {
00848
00849
00850
00851 TFitter *hFitter = (TFitter*)TVirtualFitter::GetFitter();
00852 hFitter->FitChisquare(npar, gin, f, u, flag);
00853 }
00854
00855
00856 void H1FitLikelihood(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
00857 {
00858
00859
00860
00861
00862
00863
00864
00865 TFitter *hFitter = (TFitter*)TVirtualFitter::GetFitter();
00866 hFitter->FitLikelihood(npar, gin, f, u, flag);
00867 }
00868
00869
00870 void GraphFitChisquare(Int_t &npar, Double_t * , Double_t &f,
00871 Double_t *u, Int_t )
00872 {
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905
00906
00907
00908
00909 Double_t cu,eu,exh,exl,ey,eux,fu,fsum;
00910 Double_t x[1];
00911
00912 Int_t bin, npfits=0;
00913
00914 TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
00915 TGraph *gr = (TGraph*)grFitter->GetObjectFit();
00916 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
00917 Foption_t fitOption = grFitter->GetFitOption();
00918 Int_t n = gr->GetN();
00919 Double_t *gx = gr->GetX();
00920 Double_t *gy = gr->GetY();
00921
00922
00923 npar = f1->GetNpar();
00924
00925 f = 0;
00926 for (bin=0;bin<n;bin++) {
00927 f1->InitArgs(x,u);
00928 x[0] = gx[bin];
00929 if (!f1->IsInside(x)) continue;
00930 cu = gy[bin];
00931 TF1::RejectPoint(kFALSE);
00932 fu = f1->EvalPar(x,u);
00933 if (TF1::RejectedPoint()) continue;
00934 fsum = (cu-fu);
00935 npfits++;
00936 if (fitOption.W1) {
00937 f += fsum*fsum;
00938 continue;
00939 }
00940
00941 exh = gr->GetErrorXhigh(bin);
00942 exl = gr->GetErrorXlow(bin);
00943 if (fsum < 0)
00944 ey = gr->GetErrorYhigh(bin);
00945 else
00946 ey = gr->GetErrorYlow(bin);
00947 if (exl < 0) exl = 0;
00948 if (exh < 0) exh = 0;
00949 if (ey < 0) ey = 0;
00950 if (exh > 0 || exl > 0) {
00951
00952
00953
00954
00955
00956
00957
00958
00959
00960
00961
00962 eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
00963 } else
00964 eux = 0.;
00965 eu = ey*ey+eux*eux;
00966 if (eu <= 0) eu = 1;
00967 f += fsum*fsum/eu;
00968 }
00969 f1->SetNumberFitPoints(npfits);
00970
00971 }
00972
00973
00974
00975 void Graph2DFitChisquare(Int_t &npar, Double_t * , Double_t &f,
00976 Double_t *u, Int_t )
00977 {
00978
00979
00980
00981 Double_t cu,eu,ex,ey,ez,eux,euy,fu,fsum,fm,fp;
00982 Double_t x[2];
00983 Double_t xm,xp,ym,yp;
00984 Int_t bin, npfits=0;
00985
00986 TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
00987 TGraph2D *gr = (TGraph2D*)grFitter->GetObjectFit();
00988 TF2 *f2 = (TF2*)grFitter->GetUserFunc();
00989 Foption_t fitOption = grFitter->GetFitOption();
00990
00991 Int_t n = gr->GetN();
00992 Double_t *gx = gr->GetX();
00993 Double_t *gy = gr->GetY();
00994 Double_t *gz = gr->GetZ();
00995 Double_t fxmin = f2->GetXmin();
00996 Double_t fxmax = f2->GetXmax();
00997 Double_t fymin = f2->GetYmin();
00998 Double_t fymax = f2->GetYmax();
00999 npar = f2->GetNpar();
01000 f = 0;
01001 for (bin=0;bin<n;bin++) {
01002 f2->InitArgs(x,u);
01003 x[0] = gx[bin];
01004 x[1] = gy[bin];
01005 if (!f2->IsInside(x)) continue;
01006 cu = gz[bin];
01007 TF2::RejectPoint(kFALSE);
01008 fu = f2->EvalPar(x,u);
01009 if (TF2::RejectedPoint()) continue;
01010 fsum = (cu-fu);
01011 npfits++;
01012 if (fitOption.W1) {
01013 f += fsum*fsum;
01014 continue;
01015 }
01016 ex = gr->GetErrorX(bin);
01017 ey = gr->GetErrorY(bin);
01018 ez = gr->GetErrorZ(bin);
01019 if (ex < 0) ex = 0;
01020 if (ey < 0) ey = 0;
01021 if (ez < 0) ez = 0;
01022 eux = euy = 0;
01023 if (ex > 0) {
01024 xm = x[0] - ex; if (xm < fxmin) xm = fxmin;
01025 xp = x[0] + ex; if (xp > fxmax) xp = fxmax;
01026 x[0] = xm; fm = f2->EvalPar(x,u);
01027 x[0] = xp; fp = f2->EvalPar(x,u);
01028 eux = fp-fm;
01029 }
01030 if (ey > 0) {
01031 x[0] = gx[bin];
01032 ym = x[1] - ey; if (ym < fymin) ym = fxmin;
01033 yp = x[1] + ey; if (yp > fymax) yp = fymax;
01034 x[1] = ym; fm = f2->EvalPar(x,u);
01035 x[1] = yp; fp = f2->EvalPar(x,u);
01036 euy = fp-fm;
01037 }
01038 eu = ez*ez+eux*eux+euy*euy;
01039 if (eu <= 0) eu = 1;
01040 f += fsum*fsum/eu;
01041 }
01042 f2->SetNumberFitPoints(npfits);
01043 }
01044
01045
01046 void MultiGraphFitChisquare(Int_t &npar, Double_t * , Double_t &f,
01047 Double_t *u, Int_t )
01048 {
01049
01050 Double_t cu,eu,exh,exl,ey,eux,fu,fsum;
01051 Double_t x[1];
01052
01053 Int_t bin, npfits=0;
01054
01055 TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
01056 TMultiGraph *mg = (TMultiGraph*)grFitter->GetObjectFit();
01057 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
01058 Foption_t fitOption = grFitter->GetFitOption();
01059 TGraph *gr;
01060 TIter next(mg->GetListOfGraphs());
01061
01062 Int_t n;
01063 Double_t *gx;
01064 Double_t *gy;
01065
01066
01067 npar = f1->GetNpar();
01068
01069 f = 0;
01070
01071 while ((gr = (TGraph*) next())) {
01072 n = gr->GetN();
01073 gx = gr->GetX();
01074 gy = gr->GetY();
01075 for (bin=0;bin<n;bin++) {
01076 f1->InitArgs(x,u);
01077 x[0] = gx[bin];
01078 if (!f1->IsInside(x)) continue;
01079 cu = gy[bin];
01080 TF1::RejectPoint(kFALSE);
01081 fu = f1->EvalPar(x,u);
01082 if (TF1::RejectedPoint()) continue;
01083 fsum = (cu-fu);
01084 npfits++;
01085 if (fitOption.W1) {
01086 f += fsum*fsum;
01087 continue;
01088 }
01089 exh = gr->GetErrorXhigh(bin);
01090 exl = gr->GetErrorXlow(bin);
01091 ey = gr->GetErrorY(bin);
01092 if (exl < 0) exl = 0;
01093 if (exh < 0) exh = 0;
01094 if (ey < 0) ey = 0;
01095 if (exh > 0 && exl > 0) {
01096
01097
01098
01099
01100
01101
01102
01103
01104
01105
01106
01107 eux = 0.5*(exl + exh)*f1->Derivative(x[0], u);
01108 } else
01109 eux = 0.;
01110 eu = ey*ey+eux*eux;
01111 if (eu <= 0) eu = 1;
01112 f += fsum*fsum/eu;
01113 }
01114 }
01115 f1->SetNumberFitPoints(npfits);
01116 }
01117
01118
01119 void F2Fit(Int_t &, Double_t * , Double_t &f,Double_t *u, Int_t )
01120 {
01121 TVirtualFitter *fitter = TVirtualFitter::GetFitter();
01122 TF2 *f2 = (TF2*)fitter->GetObjectFit();
01123 f2->InitArgs(u, f2->GetParameters() );
01124 f = f2->EvalPar(u);
01125 }
01126
01127
01128 void F3Fit(Int_t &, Double_t * , Double_t &f,Double_t *u, Int_t )
01129 {
01130 TVirtualFitter *fitter = TVirtualFitter::GetFitter();
01131 TF3 *f3 = (TF3*)fitter->GetObjectFit();
01132 f3->InitArgs(u, f3->GetParameters() );
01133 f = f3->EvalPar(u);
01134 }