61 ClassImp(RooFitResult);
69 RooFitResult::RooFitResult(
const char* name,
const char* title) :
70 TNamed(name,title), _constPars(0), _initPars(0), _finalPars(0), _globalCorr(0), _randomPars(0), _Lt(0),
71 _CM(0), _VM(0), _GC(0)
73 if (name) appendToDir(
this,kTRUE) ;
80 RooFitResult::RooFitResult(
const RooFitResult& other) :
84 _status(other._status),
85 _covQual(other._covQual),
86 _numBadNLL(other._numBadNLL),
87 _minNLL(other._minNLL),
95 _statusHistory(other._statusHistory)
97 _constPars = (RooArgList*) other._constPars->snapshot() ;
98 _initPars = (RooArgList*) other._initPars->snapshot() ;
99 _finalPars = (RooArgList*) other._finalPars->snapshot() ;
100 if (other._randomPars) _randomPars = (RooArgList*) other._randomPars->snapshot() ;
101 if (other._Lt) _Lt =
new TMatrix(*other._Lt);
102 if (other._VM) _VM =
new TMatrixDSym(*other._VM) ;
103 if (other._CM) _CM =
new TMatrixDSym(*other._CM) ;
104 if (other._GC) _GC =
new TVectorD(*other._GC) ;
107 appendToDir(
this, kTRUE);
115 RooFitResult::~RooFitResult()
117 if (_constPars)
delete _constPars ;
118 if (_initPars)
delete _initPars ;
119 if (_finalPars)
delete _finalPars ;
120 if (_globalCorr)
delete _globalCorr;
121 if (_randomPars)
delete _randomPars;
123 if (_CM)
delete _CM ;
124 if (_VM)
delete _VM ;
125 if (_GC)
delete _GC ;
127 _corrMatrix.RemoveAll();
128 _corrMatrix.Delete();
130 removeFromDir(
this) ;
137 void RooFitResult::setConstParList(
const RooArgList& list)
139 if (_constPars)
delete _constPars ;
140 _constPars = (RooArgList*) list.snapshot() ;
141 TIterator* iter = _constPars->createIterator() ;
143 while((arg=(RooAbsArg*)iter->Next())) {
144 RooRealVar* rrv =
dynamic_cast<RooRealVar*
>(arg) ;
146 rrv->deleteSharedProperties() ;
157 void RooFitResult::setInitParList(
const RooArgList& list)
159 if (_initPars)
delete _initPars ;
160 _initPars = (RooArgList*) list.snapshot() ;
161 TIterator* iter = _initPars->createIterator() ;
163 while((arg=(RooAbsArg*)iter->Next())) {
164 RooRealVar* rrv =
dynamic_cast<RooRealVar*
>(arg) ;
166 rrv->deleteSharedProperties() ;
177 void RooFitResult::setFinalParList(
const RooArgList& list)
179 if (_finalPars)
delete _finalPars ;
180 _finalPars = (RooArgList*) list.snapshot() ;
182 TIterator* iter = _finalPars->createIterator() ;
184 while((arg=(RooAbsArg*)iter->Next())) {
185 RooRealVar* rrv =
dynamic_cast<RooRealVar*
>(arg) ;
187 rrv->deleteSharedProperties() ;
197 Int_t RooFitResult::statusCodeHistory(UInt_t icycle)
const
199 if (icycle>=_statusHistory.size()) {
200 coutE(InputArguments) <<
"RooFitResult::statusCodeHistory(" << GetName()
201 <<
" ERROR request for status history slot "
202 << icycle <<
" exceeds history count of " << _statusHistory.size() << endl ;
204 return _statusHistory[icycle].second ;
211 const char* RooFitResult::statusLabelHistory(UInt_t icycle)
const
213 if (icycle>=_statusHistory.size()) {
214 coutE(InputArguments) <<
"RooFitResult::statusLabelHistory(" << GetName()
215 <<
" ERROR request for status history slot "
216 << icycle <<
" exceeds history count of " << _statusHistory.size() << endl ;
218 return _statusHistory[icycle].first.c_str() ;
248 RooPlot *RooFitResult::plotOn(RooPlot *frame,
const char *parName1,
const char *parName2,
249 const char *options)
const
252 const RooRealVar *par1=
dynamic_cast<const RooRealVar*
>(floatParsFinal().find(parName1));
254 coutE(InputArguments) <<
"RooFitResult::correlationPlot: parameter not floated in fit: " << parName1 << endl;
257 const RooRealVar *par2=
dynamic_cast<const RooRealVar*
>(floatParsFinal().find(parName2));
259 coutE(InputArguments) <<
"RooFitResult::correlationPlot: parameter not floated in fit: " << parName2 << endl;
264 TString opt(options);
268 Double_t x1= par1->getVal();
269 Double_t x2= par2->getVal();
270 Double_t s1= par1->getError();
271 Double_t s2= par2->getError();
272 Double_t rho= correlation(parName1, parName2);
275 if(opt.Contains(
"E")) {
276 RooEllipse *contour=
new RooEllipse(
"contour",x1,x2,s1,s2,rho);
277 contour->SetLineWidth(2) ;
278 frame->addPlotable(contour);
282 if(opt.Contains(
"1")) {
283 TLine *hline=
new TLine(x1-s1,x2,x1+s1,x2);
284 hline->SetLineColor(kRed);
285 frame->addObject(hline);
288 if(opt.Contains(
"2")) {
289 TLine *vline=
new TLine(x1,x2-s2,x1,x2+s2);
290 vline->SetLineColor(kRed);
291 frame->addObject(vline);
294 if(opt.Contains(
"B")) {
295 TBox *box=
new TBox(x1-s1,x2-s2,x1+s1,x2+s2);
296 box->SetLineStyle(kDashed);
297 box->SetLineColor(kRed);
298 box->SetFillStyle(0);
299 frame->addObject(box);
302 if(opt.Contains(
"H")) {
303 TLine *line=
new TLine(x1-rho*s1,x2-s2,x1+rho*s1,x2+s2);
304 line->SetLineStyle(kDashed);
305 line->SetLineColor(kBlue);
306 line->SetLineWidth(2) ;
307 frame->addObject(line);
308 if(opt.Contains(
"A")) {
309 TGaxis *axis=
new TGaxis(x1-s1,x2-s2,x1+s1,x2-s2,-1.,+1.,502,
"-=");
310 axis->SetLineColor(kBlue);
311 frame->addObject(axis);
315 if(opt.Contains(
"V")) {
316 TLine *line=
new TLine(x1-s1,x2-rho*s2,x1+s1,x2+rho*s2);
317 line->SetLineStyle(kDashed);
318 line->SetLineColor(kBlue);
319 line->SetLineWidth(2) ;
320 frame->addObject(line);
321 if(opt.Contains(
"A")) {
322 TGaxis *axis=
new TGaxis(x1-s1,x2-s2,x1-s1,x2+s2,-1.,+1.,502,
"-=");
323 axis->SetLineColor(kBlue);
324 frame->addObject(axis);
329 if(opt.Contains(
"M")) {
330 TMarker *marker=
new TMarker(x1,x2,20);
331 marker->SetMarkerColor(kBlack);
332 frame->addObject(marker);
346 const RooArgList& RooFitResult::randomizePars()
const
348 Int_t nPar= _finalPars->getSize();
349 if(0 == _randomPars) {
350 assert(0 != _finalPars);
352 _randomPars= (RooArgList*)_finalPars->snapshot();
355 TMatrix L(nPar,nPar);
356 for(Int_t iPar= 0; iPar < nPar; iPar++) {
358 L(iPar,iPar)= covariance(iPar,iPar);
359 for(Int_t k= 0; k < iPar; k++) {
360 Double_t tmp= L(k,iPar);
361 L(iPar,iPar)-= tmp*tmp;
363 L(iPar,iPar)= sqrt(L(iPar,iPar));
365 for(Int_t jPar= iPar+1; jPar < nPar; jPar++) {
366 L(iPar,jPar)= covariance(iPar,jPar);
367 for(Int_t k= 0; k < iPar; k++) {
368 L(iPar,jPar)-= L(k,iPar)*L(k,jPar);
370 L(iPar,jPar)/= L(iPar,iPar);
374 _Lt=
new TMatrix(TMatrix::kTransposed,L);
378 *_randomPars= *_finalPars;
383 for(Int_t k= 0; k < nPar; k++) g(k)= RooRandom::gaussian();
387 TIterator *iter= _randomPars->createIterator();
390 while((0 != (par= (RooRealVar*)iter->Next()))) {
391 par->setVal(par->getVal() + g(index++));
402 Double_t RooFitResult::correlation(
const char* parname1,
const char* parname2)
const
404 Int_t idx1 = _finalPars->index(parname1) ;
405 Int_t idx2 = _finalPars->index(parname2) ;
407 coutE(InputArguments) <<
"RooFitResult::correlation(" << GetName() <<
") parameter " << parname1 <<
" is not a floating fit parameter" << endl ;
411 coutE(InputArguments) <<
"RooFitResult::correlation(" << GetName() <<
") parameter " << parname2 <<
" is not a floating fit parameter" << endl ;
414 return correlation(idx1,idx2) ;
423 const RooArgList* RooFitResult::correlation(
const char* parname)
const
425 if (_globalCorr==0) {
426 fillLegacyCorrMatrix() ;
429 RooAbsArg* arg = _initPars->find(parname) ;
431 coutE(InputArguments) <<
"RooFitResult::correlation: variable " << parname <<
" not a floating parameter in fit" << endl ;
434 return (RooArgList*)_corrMatrix.At(_initPars->index(arg)) ;
442 Double_t RooFitResult::globalCorr(
const char* parname)
444 if (_globalCorr==0) {
445 fillLegacyCorrMatrix() ;
448 RooAbsArg* arg = _initPars->find(parname) ;
450 coutE(InputArguments) <<
"RooFitResult::globalCorr: variable " << parname <<
" not a floating parameter in fit" << endl ;
455 return ((RooAbsReal*)_globalCorr->at(_initPars->index(arg)))->getVal() ;
466 const RooArgList* RooFitResult::globalCorr()
468 if (_globalCorr==0) {
469 fillLegacyCorrMatrix() ;
480 Double_t RooFitResult::correlation(Int_t row, Int_t col)
const
482 return (*_CM)(row,col) ;
489 Double_t RooFitResult::covariance(Int_t row, Int_t col)
const
491 return (*_VM)(row,col) ;
501 void RooFitResult::printMultiline(ostream& os, Int_t , Bool_t verbose, TString indent)
const
505 << indent <<
" RooFitResult: minimized FCN value: " << _minNLL <<
", estimated distance to minimum: " << _edm << endl
506 << indent <<
" covariance matrix quality: " ;
508 case -1 : os <<
"Unknown, matrix was externally provided" ; break ;
509 case 0 : os <<
"Not calculated at all" ; break ;
510 case 1 : os <<
"Approximation only, not accurate" ; break ;
511 case 2 : os <<
"Full matrix, but forced positive-definite" ; break ;
512 case 3 : os <<
"Full, accurate covariance matrix" ; break ;
515 os << indent <<
" Status : " ;
516 for (vector<pair<string,int> >::const_iterator iter = _statusHistory.begin() ; iter != _statusHistory.end() ; ++iter) {
517 os << iter->first <<
"=" << iter->second <<
" " ;
519 os << endl << endl ;;
523 if (_constPars->getSize()>0) {
524 os << indent <<
" Constant Parameter Value " << endl
525 << indent <<
" -------------------- ------------" << endl ;
527 for (i=0 ; i<_constPars->getSize() ; i++) {
528 os << indent <<
" " << setw(20) << ((RooAbsArg*)_constPars->at(i))->GetName()
529 <<
" " << setw(12) << Form(
"%12.4e",((RooRealVar*)_constPars->at(i))->getVal())
537 Bool_t doAsymErr(kFALSE) ;
538 for (i=0 ; i<_finalPars->getSize() ; i++) {
539 if (((RooRealVar*)_finalPars->at(i))->hasAsymError()) {
546 os << indent <<
" Floating Parameter InitialValue FinalValue (+HiError,-LoError) GblCorr." << endl
547 << indent <<
" -------------------- ------------ ---------------------------------- --------" << endl ;
549 os << indent <<
" Floating Parameter InitialValue FinalValue +/- Error GblCorr." << endl
550 << indent <<
" -------------------- ------------ -------------------------- --------" << endl ;
553 for (i=0 ; i<_finalPars->getSize() ; i++) {
554 os << indent <<
" " << setw(20) << ((RooAbsArg*)_finalPars->at(i))->GetName() ;
555 os << indent <<
" " << setw(12) << Form(
"%12.4e",((RooRealVar*)_initPars->at(i))->getVal())
556 << indent <<
" " << setw(12) << Form(
"%12.4e",((RooRealVar*)_finalPars->at(i))->getVal()) ;
558 if (((RooRealVar*)_finalPars->at(i))->hasAsymError()) {
559 os << setw(21) << Form(
" (+%8.2e,-%8.2e)",((RooRealVar*)_finalPars->at(i))->getAsymErrorHi(),
560 -1*((RooRealVar*)_finalPars->at(i))->getAsymErrorLo()) ;
562 Double_t err = ((RooRealVar*)_finalPars->at(i))->getError() ;
563 os << (doAsymErr?
" ":
"") <<
" +/- " << setw(9) << Form(
"%9.2e",err) ;
567 os <<
" " << setw(8) << Form(
"%8.6f" ,((RooRealVar*)_globalCorr->at(i))->getVal()) ;
576 os << indent <<
" Floating Parameter FinalValue +/- Error " << endl
577 << indent <<
" -------------------- --------------------------" << endl ;
579 for (i=0 ; i<_finalPars->getSize() ; i++) {
580 Double_t err = ((RooRealVar*)_finalPars->at(i))->getError() ;
581 os << indent <<
" " << setw(20) << ((RooAbsArg*)_finalPars->at(i))->GetName()
582 <<
" " << setw(12) << Form(
"%12.4e",((RooRealVar*)_finalPars->at(i))->getVal())
583 <<
" +/- " << setw(9) << Form(
"%9.2e",err)
596 void RooFitResult::fillCorrMatrix(
const std::vector<double>& globalCC,
const TMatrixDSym& corrs,
const TMatrixDSym& covs)
599 if (globalCC.empty() || corrs.GetNoElements() < 1 || covs.GetNoElements() < 1) {
600 coutI(Minimization) <<
"RooFitResult::fillCorrMatrix: number of floating parameters is zero, correlation matrix not filled" << endl ;
605 coutE(Minimization) <<
"RooFitResult::fillCorrMatrix: ERROR: list of initial parameters must be filled first" << endl ;
610 if (_CM)
delete _CM ;
611 if (_VM)
delete _VM ;
612 if (_GC)
delete _GC ;
615 _CM =
new TMatrixDSym(corrs) ;
616 _VM =
new TMatrixDSym(covs) ;
617 _GC =
new TVectorD(_CM->GetNcols()) ;
618 for(
int i=0 ; i<_CM->GetNcols() ; i++) {
619 (*_GC)[i] = globalCC[i] ;
631 void RooFitResult::fillLegacyCorrMatrix()
const
636 if (_globalCorr)
delete _globalCorr ;
637 _corrMatrix.Delete();
640 _globalCorr =
new RooArgList(
"globalCorrelations") ;
642 TIterator* vIter = _initPars->createIterator() ;
645 while((arg=(RooAbsArg*)vIter->Next())) {
647 TString gcName(
"GC[") ;
648 gcName.Append(arg->GetName()) ;
650 TString gcTitle(arg->GetTitle()) ;
651 gcTitle.Append(
" Global Correlation") ;
652 _globalCorr->addOwned(*(
new RooRealVar(gcName.Data(),gcTitle.Data(),0.))) ;
656 name.Append(arg->GetName()) ;
658 RooArgList* corrMatrixRow =
new RooArgList(name.Data()) ;
659 _corrMatrix.Add(corrMatrixRow) ;
660 TIterator* vIter2 = _initPars->createIterator() ;
662 while((arg2=(RooAbsArg*)vIter2->Next())) {
664 TString cName(
"C[") ;
665 cName.Append(arg->GetName()) ;
667 cName.Append(arg2->GetName()) ;
669 TString cTitle(
"Correlation between ") ;
670 cTitle.Append(arg->GetName()) ;
671 cTitle.Append(
" and ") ;
672 cTitle.Append(arg2->GetName()) ;
673 corrMatrixRow->addOwned(*(
new RooRealVar(cName.Data(),cTitle.Data(),0.))) ;
680 TIterator *gcIter = _globalCorr->createIterator() ;
681 TIterator *parIter = _finalPars->createIterator() ;
682 RooRealVar* gcVal = 0;
683 for (
unsigned int i = 0; i < (
unsigned int)_CM->GetNcols() ; ++i) {
686 gcVal = (RooRealVar*) gcIter->Next() ;
687 gcVal->setVal((*_GC)(i)) ;
690 TIterator* cIter = ((RooArgList*)_corrMatrix.At(i))->createIterator() ;
691 for (
unsigned int it = 0; it < (
unsigned int)_CM->GetNcols() ; ++it) {
692 RooRealVar* cVal = (RooRealVar*) cIter->Next() ;
693 double value = (*_CM)(i,it) ;
695 (*_CM)(i,it) = value;
714 void RooFitResult::fillCorrMatrix()
717 if (gMinuit->fNpar < 1) {
718 coutI(Minimization) <<
"RooFitResult::fillCorrMatrix: number of floating parameters is zero, correlation matrix not filled" << endl ;
723 coutE(Minimization) <<
"RooFitResult::fillCorrMatrix: ERROR: list of initial parameters must be filled first" << endl ;
728 if (_CM)
delete _CM ;
729 if (_VM)
delete _VM ;
730 if (_GC)
delete _GC ;
733 _CM =
new TMatrixDSym(_initPars->getSize()) ;
734 _VM =
new TMatrixDSym(_initPars->getSize()) ;
735 _GC =
new TVectorD(_initPars->getSize()) ;
741 Int_t ndex, i, j, m, n, it ;
743 for (i = 1; i <= gMinuit->fNpar; ++i) {
745 for (j = 1; j <= gMinuit->fNpar; ++j) {
748 ndex = m*(m-1) / 2 + n;
750 gMinuit->fMATUvline[j-1] = gMinuit->fVhmat[ndex-1] / TMath::Sqrt(TMath::Abs(gMinuit->fVhmat[ndi-1]*gMinuit->fVhmat[ndj-1]));
753 (*_GC)(i-1) = gMinuit->fGlobcc[i-1] ;
756 for (it = 1; it <= gMinuit->fNpar ; ++it) {
757 (*_CM)(i-1,it-1) = gMinuit->fMATUvline[it-1] ;
761 for (
int ii=0 ; ii<_finalPars->getSize() ; ii++) {
762 for (
int jj=0 ; jj<_finalPars->getSize() ; jj++) {
763 (*_VM)(ii,jj) = (*_CM)(ii,jj) * ((RooRealVar*)_finalPars->at(ii))->getError() * ((RooRealVar*)_finalPars->at(jj))->getError() ;
770 void RooFitResult::fillPrefitCorrMatrix()
782 _CM =
new TMatrixDSym(_initPars->getSize());
783 _VM =
new TMatrixDSym(_initPars->getSize());
784 _GC =
new TVectorD(_initPars->getSize());
786 for (
int ii = 0; ii < _finalPars->getSize(); ii++) {
788 (*_VM)(ii, ii) = ((RooRealVar *)_finalPars->at(ii))->getError() * ((RooRealVar *)_finalPars->at(ii))->getError();
800 Bool_t RooFitResult::isIdentical(const RooFitResult& other, Double_t tol, Double_t tolCorr, Bool_t )
const
803 auto deviation = [tol](
const double left,
const double right){
805 return fabs((left - right)/right) >= tol;
807 return fabs(left) >= tol;
810 auto errMsg = [](std::string msgHead,
const RooAbsReal* tv,
const RooAbsReal* ov) {
811 cout <<
"RooFitResult::isIdentical: " << msgHead <<
" " << tv->GetName() <<
" differs in value:\t"
812 << tv->getVal() <<
" vs.\t" << ov->getVal()
813 <<
"\t(" << (tv->getVal()-ov->getVal())/ov->getVal() <<
")" << endl;
816 if (deviation(_minNLL, other._minNLL)) {
817 cout <<
"RooFitResult::isIdentical: minimized value of -log(L) is different " << _minNLL <<
" vs. " << other._minNLL << endl ;
821 for (Int_t i=0 ; i<_constPars->getSize() ; i++) {
822 auto tv =
static_cast<const RooAbsReal*
>(_constPars->at(i));
823 auto ov =
static_cast<const RooAbsReal*
>(other._constPars->find(tv->GetName())) ;
825 cout <<
"RooFitResult::isIdentical: cannot find constant parameter " << _constPars->at(i)->GetName() <<
" in reference" << endl ;
828 if (ov && deviation(tv->getVal(), ov->getVal())) {
829 errMsg(
"constant parameter", tv, ov);
834 for (Int_t i=0 ; i<_initPars->getSize() ; i++) {
835 auto ov =
static_cast<const RooAbsReal*
>(other._initPars->find(_initPars->at(i)->GetName())) ;
836 auto tv =
static_cast<const RooAbsReal*
>(_initPars->at(i));
838 cout <<
"RooFitResult::isIdentical: cannot find initial parameter " << _initPars->at(i)->GetName() <<
" in reference" << endl ;
841 if (ov && deviation(tv->getVal(), ov->getVal())) {
842 errMsg(
"initial parameter", tv, ov);
847 for (Int_t i=0 ; i<_finalPars->getSize() ; i++) {
848 auto tv =
static_cast<const RooAbsReal*
>(_finalPars->at(i));
849 auto ov =
static_cast<const RooAbsReal*
>(other._finalPars->find(tv->GetName())) ;
851 cout <<
"RooFitResult::isIdentical: cannot find final parameter " << tv->GetName() <<
" in reference" << endl ;
854 if (ov && deviation(tv->getVal(), ov->getVal())) {
855 errMsg(
"final parameter", tv, ov);
860 auto deviationCorr = [tolCorr](
const double left,
const double right){
861 return fabs(left - right) >= tolCorr;
865 if (_finalPars->getSize()>1) {
867 fillLegacyCorrMatrix() ;
868 other.fillLegacyCorrMatrix() ;
870 for (Int_t i=0 ; i<_globalCorr->getSize() ; i++) {
871 auto tv =
static_cast<const RooAbsReal*
>(_globalCorr->at(i));
872 auto ov =
static_cast<const RooAbsReal*
>(other._globalCorr->find(_globalCorr->at(i)->GetName())) ;
874 cout <<
"RooFitResult::isIdentical: cannot find global correlation coefficient " << tv->GetName() <<
" in reference" << endl ;
877 if (ov && deviationCorr(tv->getVal(), ov->getVal())) {
878 errMsg(
"global correlation coefficient", tv, ov);
883 for (Int_t j=0 ; j<_corrMatrix.GetSize() ; j++) {
884 RooArgList* row = (RooArgList*) _corrMatrix.At(j) ;
885 RooArgList* orow = (RooArgList*) other._corrMatrix.At(j) ;
886 for (Int_t i=0 ; i<row->getSize() ; i++) {
887 auto tv =
static_cast<const RooAbsReal*
>(row->at(i));
888 auto ov =
static_cast<const RooAbsReal*
>(orow->find(tv->GetName())) ;
890 cout <<
"RooFitResult::isIdentical: cannot find correlation coefficient " << tv->GetName() <<
" in reference" << endl ;
893 if (ov && deviationCorr(tv->getVal(), ov->getVal())) {
894 errMsg(
"correlation coefficient", tv, ov);
910 RooFitResult* RooFitResult::lastMinuitFit(
const RooArgList& varList)
913 if (varList.getSize()>0 && varList.getSize()!=gMinuit->fNu) {
914 oocoutE((TObject*)0,InputArguments) <<
"RooFitResult::lastMinuitFit: ERROR: supplied variable list must be either empty " << endl
915 <<
" or match the number of variables of the last fit (" << gMinuit->fNu <<
")" << endl ;
920 TIterator* iter = varList.createIterator() ;
922 while((arg=(RooAbsArg*)iter->Next())) {
923 if (!dynamic_cast<RooRealVar*>(arg)) {
924 oocoutE((TObject*)0,InputArguments) <<
"RooFitResult::lastMinuitFit: ERROR: variable '" << arg->GetName() <<
"' is not of type RooRealVar" << endl ;
930 RooFitResult* r =
new RooFitResult(
"lastMinuitFit",
"Last MINUIT fit") ;
934 RooArgList constPars(
"constPars") ;
935 RooArgList floatPars(
"floatPars") ;
938 for (i = 1; i <= gMinuit->fNu; ++i) {
939 if (gMinuit->fNvarl[i-1] < 0)
continue;
940 Int_t l = gMinuit->fNiofex[i-1];
941 TString varName(gMinuit->fCpnam[i-1]) ;
942 Bool_t isConst(l==0) ;
944 Double_t xlo = gMinuit->fAlim[i-1];
945 Double_t xhi = gMinuit->fBlim[i-1];
946 Double_t xerr = gMinuit->fWerr[l-1];
947 Double_t xval = gMinuit->fU[i-1] ;
950 if (varList.getSize()==0) {
952 if ((xlo<xhi) && !isConst) {
953 var =
new RooRealVar(varName,varName,xval,xlo,xhi) ;
955 var =
new RooRealVar(varName,varName,xval) ;
957 var->setConstant(isConst) ;
960 var = (RooRealVar*) varList.at(i-1)->Clone() ;
961 var->setConstant(isConst) ;
964 var->setRange(xlo,xhi) ;
966 if (varName.CompareTo(var->GetName())) {
967 oocoutI((TObject*)0,Eval) <<
"RooFitResult::lastMinuitFit: fit parameter '" << varName
968 <<
"' stored in variable '" << var->GetName() <<
"'" << endl ;
974 constPars.addOwned(*var) ;
976 var->setError(xerr) ;
977 floatPars.addOwned(*var) ;
981 Int_t icode,npari,nparx ;
982 Double_t fmin,edm,errdef ;
983 gMinuit->mnstat(fmin,edm,errdef,npari,nparx,icode) ;
985 r->setConstParList(constPars) ;
986 r->setInitParList(floatPars) ;
987 r->setFinalParList(floatPars) ;
990 r->setCovQual(icode) ;
991 r->setStatus(gMinuit->fStatus) ;
992 r->fillCorrMatrix() ;
1003 RooFitResult *RooFitResult::prefitResult(
const RooArgList ¶mList)
1006 TIterator *iter = paramList.createIterator();
1008 while ((arg = (RooAbsArg *)iter->Next())) {
1009 if (!dynamic_cast<RooRealVar *>(arg)) {
1010 oocoutE((TObject *)0, InputArguments) <<
"RooFitResult::lastMinuitFit: ERROR: variable '" << arg->GetName()
1011 <<
"' is not of type RooRealVar" << endl;
1016 RooFitResult *r =
new RooFitResult(
"lastMinuitFit",
"Last MINUIT fit");
1020 RooArgList constPars(
"constPars");
1021 RooArgList floatPars(
"floatPars");
1024 while ((arg = (RooAbsArg *)iter->Next())) {
1025 if (arg->isConstant()) {
1026 constPars.addClone(*arg);
1028 floatPars.addClone(*arg);
1033 r->setConstParList(constPars);
1034 r->setInitParList(floatPars);
1035 r->setFinalParList(floatPars);
1040 r->fillPrefitCorrMatrix();
1048 void RooFitResult::setCovarianceMatrix(TMatrixDSym& V)
1059 _VM = (TMatrixDSym*) V.Clone() ;
1062 _CM = (TMatrixDSym*) _VM->Clone() ;
1063 for (Int_t i=0 ; i<_CM->GetNrows() ; i++) {
1064 for (Int_t j=0 ; j<_CM->GetNcols() ; j++) {
1066 (*_CM)(i,j) = (*_CM)(i,j) / sqrt((*_CM)(i,i)*(*_CM)(j,j)) ;
1070 for (Int_t i=0 ; i<_CM->GetNrows() ; i++) {
1082 TH2* RooFitResult::correlationHist(
const char* name)
const
1084 Int_t n = _CM->GetNcols() ;
1086 TH2D* hh =
new TH2D(name,name,n,0,n,n,0,n) ;
1088 for (Int_t i = 0 ; i<n ; i++) {
1089 for (Int_t j = 0 ; j<n; j++) {
1090 hh->Fill(i+0.5,n-j-0.5,(*_CM)(i,j)) ;
1092 hh->GetXaxis()->SetBinLabel(i+1,_finalPars->at(i)->GetName()) ;
1093 hh->GetYaxis()->SetBinLabel(n-i,_finalPars->at(i)->GetName()) ;
1095 hh->SetMinimum(-1) ;
1096 hh->SetMaximum(+1) ;
1108 const TMatrixDSym& RooFitResult::covarianceMatrix()
const
1120 TMatrixDSym RooFitResult::reducedCovarianceMatrix(
const RooArgList& params)
const
1122 const TMatrixDSym& V = covarianceMatrix() ;
1126 RooArgList params2 ;
1127 TIterator* iter = params.createIterator() ;
1129 while((arg=(RooAbsArg*)iter->Next())) {
1130 if (_finalPars->find(arg->GetName())) {
1133 coutW(InputArguments) <<
"RooFitResult::reducedCovarianceMatrix(" << GetName() <<
") WARNING input variable "
1134 << arg->GetName() <<
" was not a floating parameters in fit result and is ignored" << endl ;
1141 vector<int> indexMap(params2.getSize());
1142 for (
int i=0 ; i<params2.getSize() ; i++) {
1143 indexMap[i] = _finalPars->index(params2[i].GetName());
1144 assert(indexMap[i] < V.GetNrows());
1147 TMatrixDSym Vred(indexMap.size());
1148 for (
int i = 0; i < Vred.GetNrows(); ++i) {
1149 for (
int j = 0; j < Vred.GetNcols(); ++j) {
1150 Vred(i,j) = V( indexMap[i], indexMap[j]);
1169 TMatrixDSym RooFitResult::conditionalCovarianceMatrix(
const RooArgList& params)
const
1171 const TMatrixDSym& V = covarianceMatrix() ;
1174 if (V.GetNcols()==params.getSize()) {
1178 Double_t det = V.Determinant() ;
1181 coutE(Eval) <<
"RooFitResult::conditionalCovarianceMatrix(" << GetName() <<
") ERROR: covariance matrix is not positive definite (|V|="
1182 << det <<
") cannot reduce it" << endl ;
1183 throw string(
"RooFitResult::conditionalCovarianceMatrix() ERROR, input covariance matrix is not positive definite") ;
1187 RooArgList params2 ;
1188 TIterator* iter = params.createIterator() ;
1190 while((arg=(RooAbsArg*)iter->Next())) {
1191 if (_finalPars->find(arg->GetName())) {
1194 coutW(InputArguments) <<
"RooFitResult::conditionalCovarianceMatrix(" << GetName() <<
") WARNING input variable "
1195 << arg->GetName() <<
" was not a floating parameters in fit result and is ignored" << endl ;
1201 RooArgList params3 ;
1202 iter = _finalPars->createIterator() ;
1203 while((arg=(RooAbsArg*)iter->Next())) {
1204 if (params2.find(arg->GetName())) {
1211 vector<int> map1, map2 ;
1212 for (
int i=0 ; i<_finalPars->getSize() ; i++) {
1213 if (params3.find(_finalPars->at(i)->GetName())) {
1222 TMatrixDSym S11, S22 ;
1224 RooMultiVarGaussian::blockDecompose(V,map1,map2,S11,S12,S21,S22) ;
1230 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
1231 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
1234 TMatrixDSym Vred(S22bar.GetNcols()) ;
1235 for (
int i=0 ; i<Vred.GetNcols() ; i++) {
1236 for (
int j=i ; j<Vred.GetNcols() ; j++) {
1237 Vred(i,j) = (S22bar(i,j) + S22bar(j,i))/2 ;
1238 Vred(j,i) = Vred(i,j) ;
1250 const TMatrixDSym& RooFitResult::correlationMatrix()
const
1261 RooAbsPdf* RooFitResult::createHessePdf(
const RooArgSet& params)
const
1263 const TMatrixDSym& V = covarianceMatrix() ;
1264 Double_t det = V.Determinant() ;
1267 coutE(Eval) <<
"RooFitResult::createHessePdf(" << GetName() <<
") ERROR: covariance matrix is not positive definite (|V|="
1268 << det <<
") cannot construct p.d.f" << endl ;
1273 RooArgList params2 ;
1274 TIterator* iter = params.createIterator() ;
1276 while((arg=(RooAbsArg*)iter->Next())) {
1277 if (_finalPars->find(arg->GetName())) {
1280 coutW(InputArguments) <<
"RooFitResult::createHessePdf(" << GetName() <<
") WARNING input variable "
1281 << arg->GetName() <<
" was not a floating parameters in fit result and is ignored" << endl ;
1287 RooArgList params3 ;
1288 iter = _finalPars->createIterator() ;
1289 while((arg=(RooAbsArg*)iter->Next())) {
1290 if (params2.find(arg->GetName())) {
1298 if (params3.getSize()==_finalPars->getSize()) {
1301 for (Int_t i=0 ; i<_finalPars->getSize() ; i++) {
1302 RooRealVar* parclone = (RooRealVar*) _finalPars->at(i)->Clone(Form(
"%s_centralvalue",_finalPars->at(i)->GetName())) ;
1303 parclone->setConstant(kTRUE) ;
1307 string name = Form(
"pdf_%s",GetName()) ;
1308 string title = Form(
"P.d.f of %s",GetTitle()) ;
1311 RooAbsPdf* mvg =
new RooMultiVarGaussian(name.c_str(),title.c_str(),params3,mu,V) ;
1312 mvg->addOwnedComponents(mu) ;
1320 vector<int> map1, map2 ;
1321 for (
int i=0 ; i<_finalPars->getSize() ; i++) {
1322 if (params3.find(_finalPars->at(i)->GetName())) {
1331 TMatrixDSym S11, S22 ;
1333 RooMultiVarGaussian::blockDecompose(V,map1,map2,S11,S12,S21,S22) ;
1337 for (UInt_t i=0 ; i<map1.size() ; i++) {
1338 RooRealVar* parclone = (RooRealVar*) _finalPars->at(map1[i])->Clone(Form(
"%s_centralvalue",_finalPars->at(map1[i])->GetName())) ;
1339 parclone->setConstant(kTRUE) ;
1340 mu1.add(*parclone) ;
1347 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
1348 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
1351 TMatrixDSym Vred(S22bar.GetNcols()) ;
1352 for (
int i=0 ; i<Vred.GetNcols() ; i++) {
1353 for (
int j=i ; j<Vred.GetNcols() ; j++) {
1354 Vred(i,j) = (S22bar(i,j) + S22bar(j,i))/2 ;
1355 Vred(j,i) = Vred(i,j) ;
1358 string name = Form(
"pdf_%s",GetName()) ;
1359 string title = Form(
"P.d.f of %s",GetTitle()) ;
1362 RooAbsPdf* ret =
new RooMultiVarGaussian(name.c_str(),title.c_str(),params3,mu1,Vred) ;
1363 ret->addOwnedComponents(mu1) ;
1372 void RooFitResult::SetName(
const char *name)
1374 if (_dir) _dir->GetList()->Remove(
this);
1375 TNamed::SetName(name) ;
1376 if (_dir) _dir->GetList()->Add(
this);
1383 void RooFitResult::SetNameTitle(
const char *name,
const char* title)
1385 if (_dir) _dir->GetList()->Remove(
this);
1386 TNamed::SetNameTitle(name,title) ;
1387 if (_dir) _dir->GetList()->Add(
this);
1394 void RooFitResult::printName(ostream& os)
const
1403 void RooFitResult::printTitle(ostream& os)
const
1412 void RooFitResult::printClassName(ostream& os)
const
1414 os << IsA()->GetName() ;
1421 void RooFitResult::printArgs(ostream& os)
const
1423 os <<
"[constPars=" << *_constPars <<
",floatPars=" << *_finalPars <<
"]" ;
1431 void RooFitResult::printValue(ostream& os)
const
1433 os <<
"(status=" << _status <<
",FCNmin=" << _minNLL <<
",EDM=" << _edm <<
",covQual=" << _covQual <<
")" ;
1440 Int_t RooFitResult::defaultPrintContents(Option_t* )
const
1442 return kName|kClassName|kArgs|kValue ;
1449 RooPrintable::StyleOption RooFitResult::defaultPrintStyle(Option_t* opt)
const
1451 if (!opt || strlen(opt)==0) {
1454 return RooPrintable::defaultPrintStyle(opt) ;
1461 void RooFitResult::Streamer(TBuffer &R__b)
1463 if (R__b.IsReading()) {
1465 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
1467 R__b.ReadClassBuffer(RooFitResult::Class(),
this,R__v,R__s,R__c);
1468 RooAbsArg::ioStreamerPass2Finalize();
1469 _corrMatrix.SetOwner();
1472 TNamed::Streamer(R__b);
1473 RooPrintable::Streamer(R__b);
1474 RooDirItem::Streamer(R__b);
1483 R__b >> _globalCorr;
1484 _corrMatrix.Streamer(R__b);
1485 R__b.CheckByteCount(R__s, R__c, RooFitResult::IsA());
1489 _CM =
new TMatrixDSym(_finalPars->getSize()) ;
1490 _VM =
new TMatrixDSym(_CM->GetNcols()) ;
1491 _GC =
new TVectorD(_CM->GetNcols()) ;
1493 TIterator *gcIter = _globalCorr->createIterator() ;
1494 TIterator *parIter = _finalPars->createIterator() ;
1495 RooRealVar* gcVal = 0;
1496 for (
unsigned int i = 0; i < (
unsigned int)_CM->GetNcols() ; ++i) {
1499 gcVal = (RooRealVar*) gcIter->Next() ;
1500 (*_GC)(i) = gcVal->getVal() ;
1503 TIterator* cIter = ((RooArgList*)_corrMatrix.At(i))->createIterator() ;
1504 for (
unsigned int it = 0; it < (
unsigned int)_CM->GetNcols() ; ++it) {
1505 RooRealVar* cVal = (RooRealVar*) cIter->Next() ;
1506 double value = cVal->getVal() ;
1507 (*_CM)(it,i) = value ;
1508 (*_CM)(i,it) = value;
1509 (*_VM)(it,i) = value*((RooRealVar*)_finalPars->at(i))->getError()*((RooRealVar*)_finalPars->at(it))->getError() ;
1510 (*_VM)(i,it) = (*_VM)(it,i) ;
1520 R__b.WriteClassBuffer(RooFitResult::Class(),
this);