Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooFitResult.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/// \class RooFitResult
19/// RooFitResult is a container class to hold the input and output
20/// of a PDF fit to a dataset. It contains:
21///
22/// * Values of all constant parameters
23/// * Initial and final values of floating parameters with error
24/// * Correlation matrix and global correlation coefficients
25/// * NLL and EDM at minimum
26///
27/// No references to the fitted PDF and dataset are stored
28///
29
30#include <iostream>
31#include <iomanip>
32
33#include "TBuffer.h"
34#include "TMinuit.h"
35#include "TMath.h"
36#include "TMarker.h"
37#include "TLine.h"
38#include "TBox.h"
39#include "TGaxis.h"
40#include "TMatrix.h"
41#include "TVector.h"
42#include "TDirectory.h"
43#include "TClass.h"
44#include "RooFitResult.h"
45#include "RooArgSet.h"
46#include "RooArgList.h"
47#include "RooRealVar.h"
48#include "RooPlot.h"
49#include "RooEllipse.h"
50#include "RooRandom.h"
51#include "RooMsgService.h"
52#include "TH2D.h"
53#include "TMatrixDSym.h"
54#include "RooMultiVarGaussian.h"
55
56
57using namespace std;
58
60
61
62
63////////////////////////////////////////////////////////////////////////////////
64/// Constructor with name and title
65
66RooFitResult::RooFitResult(const char *name, const char *title) : TNamed(name, title)
67{
68 if (name)
69 appendToDir(this, true);
70}
71
72
73////////////////////////////////////////////////////////////////////////////////
74/// Copy constructor
75
77 TNamed(other),
78 RooPrintable(other),
79 RooDirItem(other),
80 _status(other._status),
81 _covQual(other._covQual),
82 _numBadNLL(other._numBadNLL),
83 _minNLL(other._minNLL),
84 _edm(other._edm),
85 _statusHistory(other._statusHistory)
86{
93 if (other._randomPars) {
96 }
97 if (other._Lt) _Lt = new TMatrix(*other._Lt);
98 if (other._VM) _VM = new TMatrixDSym(*other._VM) ;
99 if (other._CM) _CM = new TMatrixDSym(*other._CM) ;
100 if (other._GC) _GC = new TVectorD(*other._GC) ;
101
102 if (GetName())
103 appendToDir(this, true);
104}
105
106
107
108////////////////////////////////////////////////////////////////////////////////
109/// Destructor
110
112{
113 if (_constPars) delete _constPars ;
114 if (_initPars) delete _initPars ;
115 if (_finalPars) delete _finalPars ;
116 if (_globalCorr) delete _globalCorr;
117 if (_randomPars) delete _randomPars;
118 if (_Lt) delete _Lt;
119 if (_CM) delete _CM ;
120 if (_VM) delete _VM ;
121 if (_GC) delete _GC ;
122
125
126 removeFromDir(this) ;
127}
128
129
130////////////////////////////////////////////////////////////////////////////////
131/// Fill the list of constant parameters
132
134{
135 if (_constPars) delete _constPars ;
137 list.snapshot(*_constPars);
138 for(auto* rrv : dynamic_range_cast<RooRealVar*>(*_constPars)) {
139 if (rrv) {
140 rrv->deleteSharedProperties() ;
141 }
142 }
143}
144
145
146
147////////////////////////////////////////////////////////////////////////////////
148/// Fill the list of initial values of the floating parameters
149
151{
152 if (_initPars) delete _initPars ;
153 _initPars = new RooArgList;
154 list.snapshot(*_initPars);
155 for(auto* rrv : dynamic_range_cast<RooRealVar*>(*_initPars)) {
156 if (rrv) {
157 rrv->deleteSharedProperties() ;
158 }
159 }
160}
161
162
163
164////////////////////////////////////////////////////////////////////////////////
165/// Fill the list of final values of the floating parameters
166
168{
169 if (_finalPars) delete _finalPars ;
171 list.snapshot(*_finalPars);
172
173 for(auto* rrv : dynamic_range_cast<RooRealVar*>(*_finalPars)) {
174 if (rrv) {
175 rrv->deleteSharedProperties() ;
176 }
177 }
178}
179
180
181
182////////////////////////////////////////////////////////////////////////////////
183
185{
186 if (icycle>=_statusHistory.size()) {
187 coutE(InputArguments) << "RooFitResult::statusCodeHistory(" << GetName()
188 << " ERROR request for status history slot "
189 << icycle << " exceeds history count of " << _statusHistory.size() << endl ;
190 }
191 return _statusHistory[icycle].second ;
192}
193
194
195
196////////////////////////////////////////////////////////////////////////////////
197
199{
200 if (icycle>=_statusHistory.size()) {
201 coutE(InputArguments) << "RooFitResult::statusLabelHistory(" << GetName()
202 << " ERROR request for status history slot "
203 << icycle << " exceeds history count of " << _statusHistory.size() << endl ;
204 }
205 return _statusHistory[icycle].first.c_str() ;
206}
207
208
209
210////////////////////////////////////////////////////////////////////////////////
211/// Add objects to a 2D plot that represent the fit results for the
212/// two named parameters. The input frame with the objects added is
213/// returned, or zero in case of an error. Which objects are added
214/// are determined by the options string which should be a concatenation
215/// of the following (not case sensitive):
216///
217/// * M - a marker at the best fit result
218/// * E - an error ellipse calculated at 1-sigma using the error matrix at the minimum
219/// * 1 - the 1-sigma error bar for parameter 1
220/// * 2 - the 1-sigma error bar for parameter 2
221/// * B - the bounding box for the error ellipse
222/// * H - a line and horizontal axis for reading off the correlation coefficient
223/// * V - a line and vertical axis for reading off the correlation coefficient
224/// * A - draw axes for reading off the correlation coefficients with the H or V options
225///
226/// You can change the attributes of objects in the returned RooPlot using the
227/// various `RooPlot::getAttXxx(name)` member functions, e.g.
228/// ```
229/// plot->getAttLine("contour")->SetLineStyle(kDashed);
230/// ```
231/// Use `plot->Print()` for a list of all objects and their names (unfortunately most
232/// of the ROOT builtin graphics objects like TLine are unnamed). Drag the left mouse
233/// button along the labels of either axis button to interactively zoom in a plot.
234
235RooPlot *RooFitResult::plotOn(RooPlot *frame, const char *parName1, const char *parName2,
236 const char *options) const
237{
238 // lookup the input parameters by name: we require that they were floated in our fit
239 const RooRealVar *par1= dynamic_cast<const RooRealVar*>(floatParsFinal().find(parName1));
240 if(nullptr == par1) {
241 coutE(InputArguments) << "RooFitResult::correlationPlot: parameter not floated in fit: " << parName1 << endl;
242 return nullptr;
243 }
244 const RooRealVar *par2= dynamic_cast<const RooRealVar*>(floatParsFinal().find(parName2));
245 if(nullptr == par2) {
246 coutE(InputArguments) << "RooFitResult::correlationPlot: parameter not floated in fit: " << parName2 << endl;
247 return nullptr;
248 }
249
250 // options are not case sensitive
251 TString opt(options);
252 opt.ToUpper();
253
254 // lookup the 2x2 covariance matrix elements for these variables
255 double x1= par1->getVal();
256 double x2= par2->getVal();
257 double s1= par1->getError();
258 double s2= par2->getError();
259 double rho= correlation(parName1, parName2);
260
261 // add a 1-sigma error ellipse, if requested
262 if(opt.Contains("E")) {
263 RooEllipse *contour= new RooEllipse("contour",x1,x2,s1,s2,rho);
264 contour->SetLineWidth(2) ;
265 frame->addPlotable(contour);
266 }
267
268 // add the error bar for parameter 1, if requested
269 if(opt.Contains("1")) {
270 TLine *hline= new TLine(x1-s1,x2,x1+s1,x2);
271 hline->SetLineColor(kRed);
272 frame->addObject(hline);
273 }
274
275 if(opt.Contains("2")) {
276 TLine *vline= new TLine(x1,x2-s2,x1,x2+s2);
277 vline->SetLineColor(kRed);
278 frame->addObject(vline);
279 }
280
281 if(opt.Contains("B")) {
282 TBox *box= new TBox(x1-s1,x2-s2,x1+s1,x2+s2);
283 box->SetLineStyle(kDashed);
284 box->SetLineColor(kRed);
285 box->SetFillStyle(0);
286 frame->addObject(box);
287 }
288
289 if(opt.Contains("H")) {
290 TLine *line= new TLine(x1-rho*s1,x2-s2,x1+rho*s1,x2+s2);
293 line->SetLineWidth(2) ;
294 frame->addObject(line);
295 if(opt.Contains("A")) {
296 TGaxis *axis= new TGaxis(x1-s1,x2-s2,x1+s1,x2-s2,-1.,+1.,502,"-=");
297 axis->SetLineColor(kBlue);
298 frame->addObject(axis);
299 }
300 }
301
302 if(opt.Contains("V")) {
303 TLine *line= new TLine(x1-s1,x2-rho*s2,x1+s1,x2+rho*s2);
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);
312 }
313 }
314
315 // add a marker at the fitted value, if requested
316 if(opt.Contains("M")) {
317 TMarker *marker= new TMarker(x1,x2,20);
318 marker->SetMarkerColor(kBlack);
319 frame->addObject(marker);
320 }
321
322 return frame;
323}
324
325
326////////////////////////////////////////////////////////////////////////////////
327/// Return a list of floating parameter values that are perturbed from the final
328/// fit values by random amounts sampled from the covariance matrix. The returned
329/// object is overwritten with each call and belongs to the RooFitResult. Uses
330/// the "square root method" to decompose the covariance matrix, which makes inverting
331/// it unnecessary.
332
334{
335 Int_t nPar= _finalPars->getSize();
336 if(nullptr == _randomPars) { // first-time initialization
337 assert(nullptr != _finalPars);
338 // create the list of random values to fill
341 // calculate the elements of the upper-triangular matrix L that gives Lt*L = C
342 // where Lt is the transpose of L (the "square-root method")
343 TMatrix L(nPar,nPar);
344 for(Int_t iPar= 0; iPar < nPar; iPar++) {
345 // calculate the diagonal term first
346 L(iPar,iPar)= covariance(iPar,iPar);
347 for(Int_t k= 0; k < iPar; k++) {
348 double tmp= L(k,iPar);
349 L(iPar,iPar)-= tmp*tmp;
350 }
351 L(iPar,iPar)= sqrt(L(iPar,iPar));
352 // then the off-diagonal terms
353 for(Int_t jPar= iPar+1; jPar < nPar; jPar++) {
354 L(iPar,jPar)= covariance(iPar,jPar);
355 for(Int_t k= 0; k < iPar; k++) {
356 L(iPar,jPar)-= L(k,iPar)*L(k,jPar);
357 }
358 L(iPar,jPar)/= L(iPar,iPar);
359 }
360 }
361 // remember Lt
363 }
364 else {
365 // reset to the final fit values
367 }
368
369 // create a vector of unit Gaussian variables
370 TVector g(nPar);
371 for(Int_t k= 0; k < nPar; k++) g(k)= RooRandom::gaussian();
372 // multiply this vector by Lt to introduce the appropriate correlations
373 g*= (*_Lt);
374 // add the mean value offsets and store the results
375 Int_t index(0);
376 for(auto * par : static_range_cast<RooRealVar*>(*_randomPars)) {
377 par->setVal(par->getVal() + g(index++));
378 }
379
380 return *_randomPars;
381}
382
383
384////////////////////////////////////////////////////////////////////////////////
385/// Return the correlation between parameters 'par1' and 'par2'
386
387double RooFitResult::correlation(const char* parname1, const char* parname2) const
388{
389 Int_t idx1 = _finalPars->index(parname1) ;
390 Int_t idx2 = _finalPars->index(parname2) ;
391 if (idx1<0) {
392 coutE(InputArguments) << "RooFitResult::correlation(" << GetName() << ") parameter " << parname1 << " is not a floating fit parameter" << endl ;
393 return 0 ;
394 }
395 if (idx2<0) {
396 coutE(InputArguments) << "RooFitResult::correlation(" << GetName() << ") parameter " << parname2 << " is not a floating fit parameter" << endl ;
397 return 0 ;
398 }
399 return correlation(idx1,idx2) ;
400}
401
402
403
404////////////////////////////////////////////////////////////////////////////////
405/// Return the set of correlation coefficients of parameter 'par' with
406/// all other floating parameters
407
408const RooArgList* RooFitResult::correlation(const char* parname) const
409{
410 if (_globalCorr==nullptr) {
412 }
413
414 RooAbsArg* arg = _initPars->find(parname) ;
415 if (!arg) {
416 coutE(InputArguments) << "RooFitResult::correlation: variable " << parname << " not a floating parameter in fit" << endl ;
417 return nullptr ;
418 }
419 return (RooArgList*)_corrMatrix.At(_initPars->index(arg)) ;
420}
421
422
423
424////////////////////////////////////////////////////////////////////////////////
425/// Return the global correlation of the named parameter
426
427double RooFitResult::globalCorr(const char* parname)
428{
429 if (_globalCorr==nullptr) {
431 }
432
433 RooAbsArg* arg = _initPars->find(parname) ;
434 if (!arg) {
435 coutE(InputArguments) << "RooFitResult::globalCorr: variable " << parname << " not a floating parameter in fit" << endl ;
436 return 0 ;
437 }
438
439 if (_globalCorr) {
440 return ((RooAbsReal*)_globalCorr->at(_initPars->index(arg)))->getVal() ;
441 } else {
442 return 1.0 ;
443 }
444}
445
446
447
448////////////////////////////////////////////////////////////////////////////////
449/// Return the list of all global correlations
450
452{
453 if (_globalCorr==nullptr) {
455 }
456
457 return _globalCorr ;
458}
459
460
461
462////////////////////////////////////////////////////////////////////////////////
463/// Return a correlation matrix element addressed with numeric indices.
464
466{
467 return (*_CM)(row,col) ;
468}
469
470
471////////////////////////////////////////////////////////////////////////////////
472/// Return the covariance matrix element addressed with numeric indices.
473
475{
476 return (*_VM)(row,col) ;
477}
478
479
480
481////////////////////////////////////////////////////////////////////////////////
482/// Print fit result to stream 'os'. In Verbose mode, the constant parameters and
483/// the initial and final values of the floating parameters are printed.
484/// Standard mode only the final values of the floating parameters are printed
485
486void RooFitResult::printMultiline(ostream& os, Int_t /*contents*/, bool verbose, TString indent) const
487{
488
489 os << endl
490 << indent << " RooFitResult: minimized FCN value: " << _minNLL << ", estimated distance to minimum: " << _edm << endl
491 << indent << " covariance matrix quality: " ;
492 switch(_covQual) {
493 case -1 : os << "Unknown, matrix was externally provided" ; break ;
494 case 0 : os << "Not calculated at all" ; break ;
495 case 1 : os << "Approximation only, not accurate" ; break ;
496 case 2 : os << "Full matrix, but forced positive-definite" ; break ;
497 case 3 : os << "Full, accurate covariance matrix" ; break ;
498 }
499 os << endl ;
500 os << indent << " Status : " ;
501 for (vector<pair<string,int> >::const_iterator iter = _statusHistory.begin() ; iter != _statusHistory.end() ; ++iter) {
502 os << iter->first << "=" << iter->second << " " ;
503 }
504 os << endl << endl ;;
505
506 Int_t i ;
507 if (verbose) {
508 if (_constPars->getSize()>0) {
509 os << indent << " Constant Parameter Value " << endl
510 << indent << " -------------------- ------------" << endl ;
511
512 for (i=0 ; i<_constPars->getSize() ; i++) {
513 os << indent << " " << setw(20) << _constPars->at(i)->GetName() << " " << setw(12);
514 if(RooRealVar* v = dynamic_cast<RooRealVar*>(_constPars->at(i))) {
515 os << TString::Format("%12.4e",v->getVal());
516 } else {
517 _constPars->at(i)->printValue(os); // for anything other than RooRealVar use printValue method to print
518 }
519 os << endl ;
520 }
521
522 os << endl ;
523 }
524
525 // Has any parameter asymmetric errors?
526 bool doAsymErr(false) ;
527 for (i=0 ; i<_finalPars->getSize() ; i++) {
528 if (((RooRealVar*)_finalPars->at(i))->hasAsymError()) {
529 doAsymErr=true ;
530 break ;
531 }
532 }
533
534 if (doAsymErr) {
535 os << indent << " Floating Parameter InitialValue FinalValue (+HiError,-LoError) GblCorr." << endl
536 << indent << " -------------------- ------------ ---------------------------------- --------" << endl ;
537 } else {
538 os << indent << " Floating Parameter InitialValue FinalValue +/- Error GblCorr." << endl
539 << indent << " -------------------- ------------ -------------------------- --------" << endl ;
540 }
541
542 for (i=0 ; i<_finalPars->getSize() ; i++) {
543 os << indent << " " << setw(20) << ((RooAbsArg*)_finalPars->at(i))->GetName() ;
544 os << indent << " " << setw(12) << Form("%12.4e",((RooRealVar*)_initPars->at(i))->getVal())
545 << indent << " " << setw(12) << Form("%12.4e",((RooRealVar*)_finalPars->at(i))->getVal()) ;
546
547 if (((RooRealVar*)_finalPars->at(i))->hasAsymError()) {
548 os << setw(21) << Form(" (+%8.2e,-%8.2e)",((RooRealVar*)_finalPars->at(i))->getAsymErrorHi(),
549 -1*((RooRealVar*)_finalPars->at(i))->getAsymErrorLo()) ;
550 } else {
551 double err = ((RooRealVar*)_finalPars->at(i))->getError() ;
552 os << (doAsymErr?" ":"") << " +/- " << setw(9) << Form("%9.2e",err) ;
553 }
554
555 if (_globalCorr) {
556 os << " " << setw(8) << Form("%8.6f" ,((RooRealVar*)_globalCorr->at(i))->getVal()) ;
557 } else {
558 os << " <none>" ;
559 }
560
561 os << endl ;
562 }
563
564 } else {
565 os << indent << " Floating Parameter FinalValue +/- Error " << endl
566 << indent << " -------------------- --------------------------" << endl ;
567
568 for (i=0 ; i<_finalPars->getSize() ; i++) {
569 double err = ((RooRealVar*)_finalPars->at(i))->getError() ;
570 os << indent << " " << setw(20) << ((RooAbsArg*)_finalPars->at(i))->GetName()
571 << " " << setw(12) << Form("%12.4e",((RooRealVar*)_finalPars->at(i))->getVal())
572 << " +/- " << setw(9) << Form("%9.2e",err)
573 << endl ;
574 }
575 }
576
577
578 os << endl ;
579}
580
581
582////////////////////////////////////////////////////////////////////////////////
583/// Function called by RooMinimizer
584
585void RooFitResult::fillCorrMatrix(const std::vector<double>& globalCC, const TMatrixDSym& corrs, const TMatrixDSym& covs)
586{
587 // Sanity check
588 if (globalCC.empty() || corrs.GetNoElements() < 1 || covs.GetNoElements() < 1) {
589 coutI(Minimization) << "RooFitResult::fillCorrMatrix: number of floating parameters is zero, correlation matrix not filled" << endl ;
590 return ;
591 }
592
593 if (!_initPars) {
594 coutE(Minimization) << "RooFitResult::fillCorrMatrix: ERROR: list of initial parameters must be filled first" << endl ;
595 return ;
596 }
597
598 // Delete eventual previous correlation data holders
599 if (_CM) delete _CM ;
600 if (_VM) delete _VM ;
601 if (_GC) delete _GC ;
602
603 // Build holding arrays for correlation coefficients
604 _CM = new TMatrixDSym(corrs) ;
605 _VM = new TMatrixDSym(covs) ;
606 _GC = new TVectorD(_CM->GetNcols()) ;
607 for(int i=0 ; i<_CM->GetNcols() ; i++) {
608 (*_GC)[i] = globalCC[i] ;
609 }
610 //fillLegacyCorrMatrix() ;
611}
612
613
614
615
616
617////////////////////////////////////////////////////////////////////////////////
618/// Sanity check
619
621{
622 if (!_CM) return ;
623
624 // Delete eventual previous correlation data holders
625 if (_globalCorr) delete _globalCorr ;
627
628 // Build holding arrays for correlation coefficients
629 _globalCorr = new RooArgList("globalCorrelations") ;
630
631 for(RooAbsArg* arg : *_initPars) {
632 // Create global correlation value holder
633 TString gcName("GC[") ;
634 gcName.Append(arg->GetName()) ;
635 gcName.Append("]") ;
636 TString gcTitle(arg->GetTitle()) ;
637 gcTitle.Append(" Global Correlation") ;
638 _globalCorr->addOwned(std::make_unique<RooRealVar>(gcName.Data(),gcTitle.Data(),0.));
639
640 // Create array with correlation holders for this parameter
641 TString name("C[") ;
642 name.Append(arg->GetName()) ;
643 name.Append(",*]") ;
644 RooArgList* corrMatrixRow = new RooArgList(name.Data()) ;
645 _corrMatrix.Add(corrMatrixRow) ;
646 for(RooAbsArg* arg2 : *_initPars) {
647
648 TString cName("C[") ;
649 cName.Append(arg->GetName()) ;
650 cName.Append(",") ;
651 cName.Append(arg2->GetName()) ;
652 cName.Append("]") ;
653 TString cTitle("Correlation between ") ;
654 cTitle.Append(arg->GetName()) ;
655 cTitle.Append(" and ") ;
656 cTitle.Append(arg2->GetName()) ;
657 corrMatrixRow->addOwned(std::make_unique<RooRealVar>(cName.Data(),cTitle.Data(),0.));
658 }
659 }
660
661 for (unsigned int i = 0; i < (unsigned int)_CM->GetNcols() ; ++i) {
662
663 // Find the next global correlation slot to fill, skipping fixed parameters
664 auto& gcVal = static_cast<RooRealVar&>((*_globalCorr)[i]);
665 gcVal.setVal((*_GC)(i)) ; // WVE FIX THIS
666
667 // Fill a row of the correlation matrix
668 auto corrMatrixCol = static_cast<RooArgList const&>(*_corrMatrix.At(i));
669 for (unsigned int it = 0; it < (unsigned int)_CM->GetNcols() ; ++it) {
670 auto& cVal = static_cast<RooRealVar&>(corrMatrixCol[it]);
671 double value = (*_CM)(i,it) ;
672 cVal.setVal(value);
673 (*_CM)(i,it) = value;
674 }
675 }
676}
677
678
679
680
681
682////////////////////////////////////////////////////////////////////////////////
683/// Internal utility method to extract the correlation matrix and the
684/// global correlation coefficients from the MINUIT memory buffer and
685/// fill the internal arrays.
686
688{
689 // Sanity check
690 if (gMinuit->fNpar < 1) {
691 coutI(Minimization) << "RooFitResult::fillCorrMatrix: number of floating parameters is zero, correlation matrix not filled" << endl ;
692 return ;
693 }
694
695 if (!_initPars) {
696 coutE(Minimization) << "RooFitResult::fillCorrMatrix: ERROR: list of initial parameters must be filled first" << endl ;
697 return ;
698 }
699
700 // Delete eventual previous correlation data holders
701 if (_CM) delete _CM ;
702 if (_VM) delete _VM ;
703 if (_GC) delete _GC ;
704
705 // Build holding arrays for correlation coefficients
706 _CM = new TMatrixDSym(_initPars->getSize()) ;
707 _VM = new TMatrixDSym(_initPars->getSize()) ;
708 _GC = new TVectorD(_initPars->getSize()) ;
709
710 // Extract correlation information for MINUIT (code taken from TMinuit::mnmatu() )
711
712 // WVE: This code directly manipulates minuit internal workspace,
713 // if TMinuit code changes this may need updating
714 Int_t ndex, i, j, m, n, it /* nparm,id,ix */ ;
715 Int_t ndi, ndj /*, iso, isw2, isw5*/;
716 for (i = 1; i <= gMinuit->fNpar; ++i) {
717 ndi = i*(i + 1) / 2;
718 for (j = 1; j <= gMinuit->fNpar; ++j) {
719 m = TMath::Max(i,j);
720 n = TMath::Min(i,j);
721 ndex = m*(m-1) / 2 + n;
722 ndj = j*(j + 1) / 2;
723 gMinuit->fMATUvline[j-1] = gMinuit->fVhmat[ndex-1] / TMath::Sqrt(TMath::Abs(gMinuit->fVhmat[ndi-1]*gMinuit->fVhmat[ndj-1]));
724 }
725
726 (*_GC)(i-1) = gMinuit->fGlobcc[i-1] ;
727
728 // Fill a row of the correlation matrix
729 for (it = 1; it <= gMinuit->fNpar ; ++it) {
730 (*_CM)(i-1,it-1) = gMinuit->fMATUvline[it-1] ;
731 }
732 }
733
734 for (int ii=0 ; ii<_finalPars->getSize() ; ii++) {
735 for (int jj=0 ; jj<_finalPars->getSize() ; jj++) {
736 (*_VM)(ii,jj) = (*_CM)(ii,jj) * ((RooRealVar*)_finalPars->at(ii))->getError() * ((RooRealVar*)_finalPars->at(jj))->getError() ;
737 }
738 }
739}
740
741////////////////////////////////////////////////////////////////////////////////
742
744{
745
746 // Delete eventual previous correlation data holders
747 if (_CM)
748 delete _CM;
749 if (_VM)
750 delete _VM;
751 if (_GC)
752 delete _GC;
753
754 // Build holding arrays for correlation coefficients
757 _GC = new TVectorD(_initPars->getSize());
758
759 for (int ii = 0; ii < _finalPars->getSize(); ii++) {
760 (*_CM)(ii, ii) = 1;
761 (*_VM)(ii, ii) = ((RooRealVar *)_finalPars->at(ii))->getError() * ((RooRealVar *)_finalPars->at(ii))->getError();
762 (*_GC)(ii) = 0;
763 }
764}
765
766
767namespace {
768
769void isIdenticalErrMsg(std::string const& msgHead, const RooAbsReal* tv, const RooAbsReal* ov, bool verbose) {
770 if(!verbose) return;
771 std::cout << "RooFitResult::isIdentical: " << msgHead << " " << tv->GetName() << " differs in value:\t"
772 << tv->getVal() << " vs.\t" << ov->getVal()
773 << "\t(" << (tv->getVal()-ov->getVal())/ov->getVal() << ")" << std::endl;
774}
775
776void isErrorIdenticalErrMsg(std::string const& msgHead, const RooRealVar* tv, const RooRealVar* ov, bool verbose) {
777 if(!verbose) return;
778 std::cout << "RooFitResult::isIdentical: " << msgHead << " " << tv->GetName() << " differs in error:\t"
779 << tv->getError() << " vs.\t" << ov->getError()
780 << "\t(" << (tv->getError()-ov->getError())/ov->getError() << ")" << std::endl;
781}
782
783} // namespace
784
785
786////////////////////////////////////////////////////////////////////////////////
787/// Return true if this fit result is identical to other within tolerances, ignoring the correlation matrix.
788/// \param[in] other Fit result to test against.
789/// \param[in] tol **Relative** tolerance for parameters and NLL.
790/// \param[in] tolErr **Relative** tolerance for parameter errors.
791/// \param[in] verbose If this function will log to the standard output when comparisons fail.
792
793bool RooFitResult::isIdenticalNoCov(const RooFitResult& other, double tol, double tolErr, bool verbose) const
794{
795 bool ret = true;
796 auto deviation = [](const double left, const double right, double tolerance){
797 return right != 0. ? std::abs((left - right)/right) >= tolerance : std::abs(left) >= tolerance;
798 };
799
800 auto compare = [&](RooArgList const& pars, RooArgList const& otherpars, std::string const& prefix, bool isVerbose) {
801 bool out = true;
802
803 for (auto * tv : static_range_cast<const RooAbsReal*>(pars)) {
804 auto ov = static_cast<const RooAbsReal*>(otherpars.find(tv->GetName())) ;
805
806 // Check in the parameter is in the other fit result
807 if (!ov) {
808 if(verbose) cout << "RooFitResult::isIdentical: cannot find " << prefix << " " << tv->GetName() << " in reference" << endl ;
809 out = false;
810 }
811
812 // Compare parameter value
813 if (ov && deviation(tv->getVal(), ov->getVal(), tol)) {
814 isIdenticalErrMsg(prefix, tv, ov, isVerbose);
815 out = false;
816 }
817
818 // Compare parameter error if it's a RooRealVar
819 auto * rtv = dynamic_cast<RooRealVar const*>(tv);
820 auto * rov = dynamic_cast<RooRealVar const*>(ov);
821 if(rtv && rov) {
822 if (ov && deviation(rtv->getError(), rov->getError(), tolErr)) {
823 isErrorIdenticalErrMsg(prefix, rtv, rov, isVerbose);
824 out = false;
825 }
826 }
827 }
828
829 return out;
830 };
831
832 if (deviation(_minNLL, other._minNLL, tol)) {
833 if(verbose) std::cout << "RooFitResult::isIdentical: minimized value of -log(L) is different " << _minNLL << " vs. " << other._minNLL << std::endl;
834 ret = false;
835 }
836
837 ret &= compare(*_constPars, *other._constPars, "constant parameter", verbose);
838 ret &= compare(*_initPars, *other._initPars, "initial parameter", verbose);
839 ret &= compare(*_finalPars, *other._finalPars, "final parameter", verbose);
840
841 return ret;
842}
843
844
845////////////////////////////////////////////////////////////////////////////////
846/// Return true if this fit result is identical to other within tolerances.
847/// \param[in] other Fit result to test against.
848/// \param[in] tol **Relative** tolerance for parameters and NLL.
849/// \param[in] tolCorr **absolute** tolerance for correlation coefficients.
850/// \param[in] verbose If this function will log to the standard output when comparisons fail.
851///
852/// As the relative tolerance for the parameter errors, the default value of
853/// `1e-3` will be used.
854
855bool RooFitResult::isIdentical(const RooFitResult& other, double tol, double tolCorr, bool verbose) const
856{
857 bool ret = isIdenticalNoCov(other, tol, 1e-3 /* synced with default parameter*/, verbose);
858
859 auto deviationCorr = [tolCorr](const double left, const double right){
860 return std::abs(left - right) >= tolCorr;
861 };
862
863 // Only examine correlations for cases with >1 floating parameter
864 if (_finalPars->getSize()>1) {
865
867 other.fillLegacyCorrMatrix() ;
868
869 for (Int_t i=0 ; i<_globalCorr->getSize() ; i++) {
870 auto tv = static_cast<const RooAbsReal*>(_globalCorr->at(i));
871 auto ov = static_cast<const RooAbsReal*>(other._globalCorr->find(_globalCorr->at(i)->GetName())) ;
872 if (!ov) {
873 if(verbose) cout << "RooFitResult::isIdentical: cannot find global correlation coefficient " << tv->GetName() << " in reference" << endl ;
874 ret = false ;
875 }
876 if (ov && deviationCorr(tv->getVal(), ov->getVal())) {
877 isIdenticalErrMsg("global correlation coefficient", tv, ov, verbose);
878 ret = false ;
879 }
880 }
881
882 for (Int_t j=0 ; j<_corrMatrix.GetSize() ; j++) {
883 RooArgList* row = (RooArgList*) _corrMatrix.At(j) ;
884 RooArgList* orow = (RooArgList*) other._corrMatrix.At(j) ;
885 for (Int_t i=0 ; i<row->getSize() ; i++) {
886 auto tv = static_cast<const RooAbsReal*>(row->at(i));
887 auto ov = static_cast<const RooAbsReal*>(orow->find(tv->GetName())) ;
888 if (!ov) {
889 if(verbose) cout << "RooFitResult::isIdentical: cannot find correlation coefficient " << tv->GetName() << " in reference" << endl ;
890 ret = false ;
891 }
892 if (ov && deviationCorr(tv->getVal(), ov->getVal())) {
893 isIdenticalErrMsg("correlation coefficient", tv, ov, verbose);
894 ret = false ;
895 }
896 }
897 }
898 }
899
900 return ret ;
901}
902
903
904
905////////////////////////////////////////////////////////////////////////////////
906/// Import the results of the last fit performed by gMinuit, interpreting
907/// the fit parameters as the given varList of parameters.
908
910{
911 // Verify length of supplied varList
912 if (varList.getSize()>0 && varList.getSize()!=gMinuit->fNu) {
913 oocoutE(nullptr,InputArguments) << "RooFitResult::lastMinuitFit: ERROR: supplied variable list must be either empty " << endl
914 << " or match the number of variables of the last fit (" << gMinuit->fNu << ")" << endl ;
915 return nullptr;
916 }
917
918 // Verify that all members of varList are of type RooRealVar
919 for(RooAbsArg* arg : varList) {
920 if (!dynamic_cast<RooRealVar*>(arg)) {
921 oocoutE(nullptr,InputArguments) << "RooFitResult::lastMinuitFit: ERROR: variable '" << arg->GetName() << "' is not of type RooRealVar" << endl ;
922 return nullptr;
923 }
924 }
925
926 RooFitResult* r = new RooFitResult("lastMinuitFit","Last MINUIT fit") ;
927
928 // Extract names of fit parameters from MINUIT
929 // and construct corresponding RooRealVars
930 RooArgList constPars("constPars") ;
931 RooArgList floatPars("floatPars") ;
932
933 Int_t i ;
934 for (i = 1; i <= gMinuit->fNu; ++i) {
935 if (gMinuit->fNvarl[i-1] < 0) continue;
936 Int_t l = gMinuit->fNiofex[i-1];
937 TString varName(gMinuit->fCpnam[i-1]) ;
938 bool isConst(l==0) ;
939
940 double xlo = gMinuit->fAlim[i-1];
941 double xhi = gMinuit->fBlim[i-1];
942 double xerr = gMinuit->fWerr[l-1];
943 double xval = gMinuit->fU[i-1] ;
944
945 std::unique_ptr<RooRealVar> var;
946 if (varList.empty()) {
947
948 if ((xlo<xhi) && !isConst) {
949 var = std::make_unique<RooRealVar>(varName,varName,xval,xlo,xhi) ;
950 } else {
951 var = std::make_unique<RooRealVar>(varName,varName,xval) ;
952 }
953 var->setConstant(isConst) ;
954 } else {
955
956 var = std::unique_ptr<RooRealVar>{static_cast<RooRealVar*>(varList.at(i-1)->Clone())};
957 var->setConstant(isConst) ;
958 var->setVal(xval) ;
959 if (xlo<xhi) {
960 var->setRange(xlo,xhi) ;
961 }
962 if (varName.CompareTo(var->GetName())) {
963 oocoutI(nullptr,Eval) << "RooFitResult::lastMinuitFit: fit parameter '" << varName
964 << "' stored in variable '" << var->GetName() << "'" << endl ;
965 }
966
967 }
968
969 if (isConst) {
970 constPars.addOwned(std::move(var));
971 } else {
972 var->setError(xerr) ;
973 floatPars.addOwned(std::move(var));
974 }
975 }
976
977 Int_t icode,npari,nparx ;
978 double fmin,edm,errdef ;
979 gMinuit->mnstat(fmin,edm,errdef,npari,nparx,icode) ;
980
981 r->setConstParList(constPars) ;
982 r->setInitParList(floatPars) ;
983 r->setFinalParList(floatPars) ;
984 r->setMinNLL(fmin) ;
985 r->setEDM(edm) ;
986 r->setCovQual(icode) ;
987 r->setStatus(gMinuit->fStatus) ;
988 r->fillCorrMatrix() ;
989
990 return r ;
991}
992
993
994
995////////////////////////////////////////////////////////////////////////////////
996/// Import the results of the last fit performed by gMinuit, interpreting
997/// the fit parameters as the given varList of parameters.
998
1000{
1001 // Verify that all members of varList are of type RooRealVar
1002 for(RooAbsArg * arg : paramList) {
1003 if (!dynamic_cast<RooRealVar *>(arg)) {
1004 oocoutE(nullptr, InputArguments) << "RooFitResult::lastMinuitFit: ERROR: variable '" << arg->GetName()
1005 << "' is not of type RooRealVar" << endl;
1006 return nullptr;
1007 }
1008 }
1009
1010 RooFitResult *r = new RooFitResult("lastMinuitFit", "Last MINUIT fit");
1011
1012 // Extract names of fit parameters from MINUIT
1013 // and construct corresponding RooRealVars
1014 RooArgList constPars("constPars");
1015 RooArgList floatPars("floatPars");
1016
1017 for(RooAbsArg* arg : paramList) {
1018 if (arg->isConstant()) {
1019 constPars.addClone(*arg);
1020 } else {
1021 floatPars.addClone(*arg);
1022 }
1023 }
1024
1025 r->setConstParList(constPars);
1026 r->setInitParList(floatPars);
1027 r->setFinalParList(floatPars);
1028 r->setMinNLL(0);
1029 r->setEDM(0);
1030 r->setCovQual(0);
1031 r->setStatus(0);
1032 r->fillPrefitCorrMatrix();
1033
1034 return r;
1035}
1036
1037////////////////////////////////////////////////////////////////////////////////
1038/// Store externally provided correlation matrix in this RooFitResult ;
1039
1041{
1042 // Delete any previous matrices
1043 if (_VM) {
1044 delete _VM ;
1045 }
1046 if (_CM) {
1047 delete _CM ;
1048 }
1049
1050 // Clone input covariance matrix ;
1051 _VM = (TMatrixDSym*) V.Clone() ;
1052
1053 // Now construct correlation matrix from it
1054 _CM = (TMatrixDSym*) _VM->Clone() ;
1055 for (Int_t i=0 ; i<_CM->GetNrows() ; i++) {
1056 for (Int_t j=0 ; j<_CM->GetNcols() ; j++) {
1057 if (i!=j) {
1058 (*_CM)(i,j) = (*_CM)(i,j) / sqrt((*_CM)(i,i)*(*_CM)(j,j)) ;
1059 }
1060 }
1061 }
1062 for (Int_t i=0 ; i<_CM->GetNrows() ; i++) {
1063 (*_CM)(i,i) = 1.0 ;
1064 }
1065
1066 _covQual = -1 ;
1067}
1068
1069
1070
1071////////////////////////////////////////////////////////////////////////////////
1072/// Return TH2D of correlation matrix
1073
1075{
1076 Int_t n = _CM->GetNcols() ;
1077
1078 TH2D* hh = new TH2D(name,name,n,0,n,n,0,n) ;
1079
1080 for (Int_t i = 0 ; i<n ; i++) {
1081 for (Int_t j = 0 ; j<n; j++) {
1082 hh->Fill(i+0.5,n-j-0.5,(*_CM)(i,j)) ;
1083 }
1084 hh->GetXaxis()->SetBinLabel(i+1,_finalPars->at(i)->GetName()) ;
1085 hh->GetYaxis()->SetBinLabel(n-i,_finalPars->at(i)->GetName()) ;
1086 }
1087 hh->SetMinimum(-1) ;
1088 hh->SetMaximum(+1) ;
1089
1090
1091 return hh ;
1092}
1093
1094
1095
1096
1097////////////////////////////////////////////////////////////////////////////////
1098/// Return covariance matrix
1099
1101{
1102 return *_VM ;
1103}
1104
1105
1106
1107
1108////////////////////////////////////////////////////////////////////////////////
1109/// Return a reduced covariance matrix (Note that Vred _is_ a simple sub-matrix of V,
1110/// row/columns are ordered to matched the convention given in input argument 'params'
1111
1113{
1114 const TMatrixDSym& V = covarianceMatrix() ;
1115
1116
1117 // Make sure that all given params were floating parameters in the represented fit
1118 RooArgList params2 ;
1119 for(RooAbsArg* arg : params) {
1120 if (_finalPars->find(arg->GetName())) {
1121 params2.add(*arg) ;
1122 } else {
1123 coutW(InputArguments) << "RooFitResult::reducedCovarianceMatrix(" << GetName() << ") WARNING input variable "
1124 << arg->GetName() << " was not a floating parameters in fit result and is ignored" << endl ;
1125 }
1126 }
1127
1128 // fix for bug ROOT-8044
1129 // use same order given bby vector params
1130 vector<int> indexMap(params2.getSize());
1131 for (int i=0 ; i<params2.getSize() ; i++) {
1132 indexMap[i] = _finalPars->index(params2[i].GetName());
1133 assert(indexMap[i] < V.GetNrows());
1134 }
1135
1136 TMatrixDSym Vred(indexMap.size());
1137 for (int i = 0; i < Vred.GetNrows(); ++i) {
1138 for (int j = 0; j < Vred.GetNcols(); ++j) {
1139 Vred(i,j) = V( indexMap[i], indexMap[j]);
1140 }
1141 }
1142 return Vred;
1143}
1144
1145
1146
1147////////////////////////////////////////////////////////////////////////////////
1148/// Return a reduced covariance matrix, which is calculated as
1149/// \f[
1150/// V_\mathrm{red} = \bar{V_{22}} = V_{11} - V_{12} \cdot V_{22}^{-1} \cdot V_{21},
1151/// \f]
1152/// where \f$ V_{11},V_{12},V_{21},V_{22} \f$ represent a block decomposition of the covariance matrix into observables that
1153/// are propagated (labeled by index '1') and that are not propagated (labeled by index '2'), and \f$ \bar{V_{22}} \f$
1154/// is the Shur complement of \f$ V_{22} \f$, calculated as shown above.
1155///
1156/// (Note that \f$ V_\mathrm{red} \f$ is *not* a simple sub-matrix of \f$ V \f$)
1157
1159{
1160 const TMatrixDSym& V = covarianceMatrix() ;
1161
1162 // Handle case where V==Vred here
1163 if (V.GetNcols()==params.getSize()) {
1164 return V ;
1165 }
1166
1167 double det = V.Determinant() ;
1168
1169 if (det<=0) {
1170 coutE(Eval) << "RooFitResult::conditionalCovarianceMatrix(" << GetName() << ") ERROR: covariance matrix is not positive definite (|V|="
1171 << det << ") cannot reduce it" << endl ;
1172 throw string("RooFitResult::conditionalCovarianceMatrix() ERROR, input covariance matrix is not positive definite") ;
1173 }
1174
1175 // Make sure that all given params were floating parameters in the represented fit
1176 RooArgList params2 ;
1177 for(RooAbsArg* arg : params) {
1178 if (_finalPars->find(arg->GetName())) {
1179 params2.add(*arg) ;
1180 } else {
1181 coutW(InputArguments) << "RooFitResult::conditionalCovarianceMatrix(" << GetName() << ") WARNING input variable "
1182 << arg->GetName() << " was not a floating parameters in fit result and is ignored" << endl ;
1183 }
1184 }
1185
1186 // Need to order params in vector in same order as in covariance matrix
1187 RooArgList params3 ;
1188 for(RooAbsArg* arg : *_finalPars) {
1189 if (params2.find(arg->GetName())) {
1190 params3.add(*arg) ;
1191 }
1192 }
1193
1194 // Find (subset) of parameters that are stored in the covariance matrix
1195 vector<int> map1, map2 ;
1196 for (int i=0 ; i<_finalPars->getSize() ; i++) {
1197 if (params3.find(_finalPars->at(i)->GetName())) {
1198 map1.push_back(i) ;
1199 } else {
1200 map2.push_back(i) ;
1201 }
1202 }
1203
1204 // Rearrange matrix in block form with 'params' first and 'others' last
1205 // (preserving relative order)
1206 TMatrixDSym S11, S22 ;
1207 TMatrixD S12, S21 ;
1208 RooMultiVarGaussian::blockDecompose(V,map1,map2,S11,S12,S21,S22) ;
1209
1210 // Constructed conditional matrix form -1
1211 // F(X1|X2) --> CovI --> S22bar = S11 - S12 S22 S21
1212
1213 // Do eigenvalue decomposition
1214 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
1215 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
1216
1217 // Convert explicitly to symmetric form
1218 TMatrixDSym Vred(S22bar.GetNcols()) ;
1219 for (int i=0 ; i<Vred.GetNcols() ; i++) {
1220 for (int j=i ; j<Vred.GetNcols() ; j++) {
1221 Vred(i,j) = (S22bar(i,j) + S22bar(j,i))/2 ;
1222 Vred(j,i) = Vred(i,j) ;
1223 }
1224 }
1225
1226 return Vred ;
1227}
1228
1229
1230
1231////////////////////////////////////////////////////////////////////////////////
1232/// Return correlation matrix ;
1233
1235{
1236 return *_CM ;
1237}
1238
1239
1240
1241////////////////////////////////////////////////////////////////////////////////
1242/// Return a p.d.f that represents the fit result as a multi-variate probability densisty
1243/// function on the floating fit parameters, including correlations
1244
1246{
1247 const TMatrixDSym& V = covarianceMatrix() ;
1248 double det = V.Determinant() ;
1249
1250 if (det<=0) {
1251 coutE(Eval) << "RooFitResult::createHessePdf(" << GetName() << ") ERROR: covariance matrix is not positive definite (|V|="
1252 << det << ") cannot construct p.d.f" << endl ;
1253 return nullptr ;
1254 }
1255
1256 // Make sure that all given params were floating parameters in the represented fit
1257 RooArgList params2 ;
1258 for(RooAbsArg* arg : params) {
1259 if (_finalPars->find(arg->GetName())) {
1260 params2.add(*arg) ;
1261 } else {
1262 coutW(InputArguments) << "RooFitResult::createHessePdf(" << GetName() << ") WARNING input variable "
1263 << arg->GetName() << " was not a floating parameters in fit result and is ignored" << endl ;
1264 }
1265 }
1266
1267 // Need to order params in vector in same order as in covariance matrix
1268 RooArgList params3 ;
1269 for(RooAbsArg* arg : *_finalPars) {
1270 if (params2.find(arg->GetName())) {
1271 params3.add(*arg) ;
1272 }
1273 }
1274
1275
1276 // Handle special case of representing full covariance matrix here
1277 if (params3.getSize()==_finalPars->getSize()) {
1278
1279 RooArgList mu ;
1280 for (Int_t i=0 ; i<_finalPars->getSize() ; i++) {
1281 RooRealVar* parclone = (RooRealVar*) _finalPars->at(i)->Clone(Form("%s_centralvalue",_finalPars->at(i)->GetName())) ;
1282 parclone->setConstant(true) ;
1283 mu.add(*parclone) ;
1284 }
1285
1286 string name = Form("pdf_%s",GetName()) ;
1287 string title = Form("P.d.f of %s",GetTitle()) ;
1288
1289 // Create p.d.f.
1290 RooAbsPdf* mvg = new RooMultiVarGaussian(name.c_str(),title.c_str(),params3,mu,V) ;
1291 mvg->addOwnedComponents(mu) ;
1292 return mvg ;
1293 }
1294
1295 // -> ->
1296 // Handle case of conditional p.d.f. MVG(p1|p2) here
1297
1298 // Find (subset) of parameters that are stored in the covariance matrix
1299 vector<int> map1, map2 ;
1300 for (int i=0 ; i<_finalPars->getSize() ; i++) {
1301 if (params3.find(_finalPars->at(i)->GetName())) {
1302 map1.push_back(i) ;
1303 } else {
1304 map2.push_back(i) ;
1305 }
1306 }
1307
1308 // Rearrange matrix in block form with 'params' first and 'others' last
1309 // (preserving relative order)
1310 TMatrixDSym S11, S22 ;
1311 TMatrixD S12, S21 ;
1312 RooMultiVarGaussian::blockDecompose(V,map1,map2,S11,S12,S21,S22) ;
1313
1314 // Calculate offset vectors mu1 and mu2
1315 RooArgList mu1 ;
1316 for (UInt_t i=0 ; i<map1.size() ; i++) {
1317 RooRealVar* parclone = (RooRealVar*) _finalPars->at(map1[i])->Clone(Form("%s_centralvalue",_finalPars->at(map1[i])->GetName())) ;
1318 parclone->setConstant(true) ;
1319 mu1.add(*parclone) ;
1320 }
1321
1322 // Constructed conditional matrix form -1
1323 // F(X1|X2) --> CovI --> S22bar = S11 - S12 S22 S21
1324
1325 // Do eigenvalue decomposition
1326 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
1327 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
1328
1329 // Convert explicitly to symmetric form
1330 TMatrixDSym Vred(S22bar.GetNcols()) ;
1331 for (int i=0 ; i<Vred.GetNcols() ; i++) {
1332 for (int j=i ; j<Vred.GetNcols() ; j++) {
1333 Vred(i,j) = (S22bar(i,j) + S22bar(j,i))/2 ;
1334 Vred(j,i) = Vred(i,j) ;
1335 }
1336 }
1337 string name = Form("pdf_%s",GetName()) ;
1338 string title = Form("P.d.f of %s",GetTitle()) ;
1339
1340 // Create p.d.f.
1341 RooAbsPdf* ret = new RooMultiVarGaussian(name.c_str(),title.c_str(),params3,mu1,Vred) ;
1342 ret->addOwnedComponents(mu1) ;
1343 return ret ;
1344}
1345
1346
1347
1348////////////////////////////////////////////////////////////////////////////////
1349/// Change name of RooFitResult object
1350
1352{
1353 if (_dir) _dir->GetList()->Remove(this);
1355 if (_dir) _dir->GetList()->Add(this);
1356}
1357
1358
1359////////////////////////////////////////////////////////////////////////////////
1360/// Change name and title of RooFitResult object
1361
1362void RooFitResult::SetNameTitle(const char *name, const char* title)
1363{
1364 if (_dir) _dir->GetList()->Remove(this);
1365 TNamed::SetNameTitle(name,title) ;
1366 if (_dir) _dir->GetList()->Add(this);
1367}
1368
1369
1370////////////////////////////////////////////////////////////////////////////////
1371/// Print name of fit result
1372
1373void RooFitResult::printName(ostream& os) const
1374{
1375 os << GetName() ;
1376}
1377
1378
1379////////////////////////////////////////////////////////////////////////////////
1380/// Print title of fit result
1381
1382void RooFitResult::printTitle(ostream& os) const
1383{
1384 os << GetTitle() ;
1385}
1386
1387
1388////////////////////////////////////////////////////////////////////////////////
1389/// Print class name of fit result
1390
1391void RooFitResult::printClassName(ostream& os) const
1392{
1393 os << ClassName() ;
1394}
1395
1396
1397////////////////////////////////////////////////////////////////////////////////
1398/// Print arguments of fit result, i.e. the parameters of the fit
1399
1400void RooFitResult::printArgs(ostream& os) const
1401{
1402 os << "[constPars=" << *_constPars << ",floatPars=" << *_finalPars << "]" ;
1403}
1404
1405
1406
1407////////////////////////////////////////////////////////////////////////////////
1408/// Print the value of the fit result, i.e.g the status, minimized FCN, edm and covariance quality code
1409
1410void RooFitResult::printValue(ostream& os) const
1411{
1412 os << "(status=" << _status << ",FCNmin=" << _minNLL << ",EDM=" << _edm << ",covQual=" << _covQual << ")" ;
1413}
1414
1415
1416////////////////////////////////////////////////////////////////////////////////
1417/// Configure default contents to be printed
1418
1420{
1421 return kName|kClassName|kArgs|kValue ;
1422}
1423
1424
1425////////////////////////////////////////////////////////////////////////////////
1426/// Configure mapping of Print() arguments to RooPrintable print styles
1427
1429{
1430 if (!opt || strlen(opt)==0) {
1431 return kStandard ;
1432 }
1434}
1435
1436
1437////////////////////////////////////////////////////////////////////////////////
1438/// Stream an object of class RooFitResult.
1439
1441{
1442 if (R__b.IsReading()) {
1443 UInt_t R__s, R__c;
1444 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
1445 if (R__v>3) {
1446 R__b.ReadClassBuffer(RooFitResult::Class(),this,R__v,R__s,R__c);
1449 } else {
1450 // backward compatibitily streaming
1451 TNamed::Streamer(R__b);
1454 R__b >> _status;
1455 R__b >> _covQual;
1456 R__b >> _numBadNLL;
1457 R__b >> _minNLL;
1458 R__b >> _edm;
1459 R__b >> _constPars;
1460 R__b >> _initPars;
1461 R__b >> _finalPars;
1462 R__b >> _globalCorr;
1463 _corrMatrix.Streamer(R__b);
1464 R__b.CheckByteCount(R__s, R__c, RooFitResult::IsA());
1465
1466 // Now fill new-style covariance and correlation matrix information
1467 // from legacy form
1468 _CM = new TMatrixDSym(_finalPars->getSize()) ;
1469 _VM = new TMatrixDSym(_CM->GetNcols()) ;
1470 _GC = new TVectorD(_CM->GetNcols()) ;
1471
1472 for (unsigned int i = 0; i < (unsigned int)_CM->GetNcols() ; ++i) {
1473
1474 // Find the next global correlation slot to fill, skipping fixed parameters
1475 auto& gcVal = static_cast<RooRealVar&>((*_globalCorr)[i]);
1476 (*_GC)(i) = gcVal.getVal() ;
1477
1478 // Fill a row of the correlation matrix
1479 auto corrMatrixCol = static_cast<RooArgList const&>(*_corrMatrix.At(i));
1480 for (unsigned int it = 0; it < (unsigned int)_CM->GetNcols() ; ++it) {
1481 auto& cVal = static_cast<RooRealVar&>(corrMatrixCol[it]);
1482 double value = cVal.getVal() ;
1483 (*_CM)(it,i) = value ;
1484 (*_CM)(i,it) = value;
1485 (*_VM)(it,i) = value*((RooRealVar*)_finalPars->at(i))->getError()*((RooRealVar*)_finalPars->at(it))->getError() ;
1486 (*_VM)(i,it) = (*_VM)(it,i) ;
1487 }
1488 }
1489 }
1490
1491 } else {
1493 }
1494}
1495
#define g(i)
Definition RSha256.hxx:105
#define s1(x)
Definition RSha256.hxx:91
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define coutW(a)
#define oocoutE(o, a)
#define oocoutI(o, a)
#define coutE(a)
short Version_t
Definition RtypesCore.h:65
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kBlue
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
static void indent(ostringstream &buf, int indent_level)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
char name[80]
Definition TGX11.cxx:110
TMatrixTSym< Double_t > TMatrixDSym
TMatrixT< Float_t > TMatrix
Definition TMatrix.h:24
R__EXTERN TMinuit * gMinuit
Definition TMinuit.h:271
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
TVectorT< Float_t > TVector
Definition TVectorfwd.h:23
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
static void ioStreamerPass2Finalize()
Method called by workspace container to finalize schema evolution issues that cannot be handled in a ...
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:91
RooAbsCollection * snapshot(bool deepCopy=true) const
Take a snap shot of current collection contents.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
void setConstant(bool value=true)
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooDirItem is a utility base class for RooFit objects that are to be attached to ROOT directories.
Definition RooDirItem.h:22
virtual void Streamer(TBuffer &)
void removeFromDir(TObject *obj)
Remove object from directory it was added to.
TDirectory * _dir
! Associated directory
Definition RooDirItem.h:33
void appendToDir(TObject *obj, bool forceMemoryResident=false)
Append object to directory.
A RooEllipse is a two-dimensional ellipse that can be used to represent an error contour.
Definition RooEllipse.h:22
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
TMatrixDSym conditionalCovarianceMatrix(const RooArgList &params) const
Return a reduced covariance matrix, which is calculated as.
void fillCorrMatrix()
Internal utility method to extract the correlation matrix and the global correlation coefficients fro...
double correlation(const RooAbsArg &par1, const RooAbsArg &par2) const
Return correlation between par1 and par2.
~RooFitResult() override
Destructor.
TList _corrMatrix
! Correlation matrix (list of RooArgLists)
std::vector< std::pair< std::string, int > > _statusHistory
History of status codes.
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
TMatrixDSym * _CM
Correlation matrix.
bool isIdentical(const RooFitResult &other, double tol=1e-6, double tolCorr=1e-4, bool verbose=true) const
Return true if this fit result is identical to other within tolerances.
void setConstParList(const RooArgList &list)
Fill the list of constant parameters.
Int_t statusCodeHistory(UInt_t icycle) const
double _minNLL
NLL at minimum.
Int_t _numBadNLL
Number calls with bad (zero,negative) likelihood.
TMatrixDSym * _VM
Covariance matrix.
void SetNameTitle(const char *name, const char *title) override
Change name and title of RooFitResult object.
Int_t defaultPrintContents(Option_t *opt) const override
Configure default contents to be printed.
void printClassName(std::ostream &os) const override
Print class name of fit result.
RooArgList * _initPars
List of floating parameters with initial values.
TClass * IsA() const override
RooFitResult(const char *name=nullptr, const char *title=nullptr)
Constructor with name and title.
TMatrixDSym reducedCovarianceMatrix(const RooArgList &params) const
Return a reduced covariance matrix (Note that Vred is a simple sub-matrix of V, row/columns are order...
Int_t _covQual
MINUIT quality code of covariance matrix.
void fillPrefitCorrMatrix()
RooArgList * _constPars
List of constant parameters.
RooArgList * _globalCorr
! List of global correlation coefficients
const RooArgList & randomizePars() const
Generate random perturbations of the final parameters using the covariance matrix.
const RooArgList * globalCorr()
Return the list of all global correlations.
static RooFitResult * prefitResult(const RooArgList &paramList)
Import the results of the last fit performed by gMinuit, interpreting the fit parameters as the given...
void setCovarianceMatrix(TMatrixDSym &V)
Store externally provided correlation matrix in this RooFitResult ;.
RooArgList * _finalPars
List of floating parameters with final values.
double edm() const
Return estimated distance to minimum.
const RooArgList & constPars() const
Return list of constant parameters.
Int_t _status
MINUIT status code.
bool isIdenticalNoCov(const RooFitResult &other, double tol=1e-6, double tolErr=1e-3, bool verbose=true) const
Return true if this fit result is identical to other within tolerances, ignoring the correlation matr...
void SetName(const char *name) override
Change name of RooFitResult object.
const char * statusLabelHistory(UInt_t icycle) const
static TClass * Class()
void printValue(std::ostream &os) const override
Print the value of the fit result, i.e.g the status, minimized FCN, edm and covariance quality code.
void printTitle(std::ostream &os) const override
Print title of fit result.
TH2 * correlationHist(const char *name="correlation_matrix") const
Return TH2D of correlation matrix.
void fillLegacyCorrMatrix() const
Sanity check.
void setInitParList(const RooArgList &list)
Fill the list of initial values of the floating parameters.
RooPlot * plotOn(RooPlot *frame, const RooAbsArg &par1, const RooAbsArg &par2, const char *options="ME") const
Add objects to a 2D plot.
void Streamer(TBuffer &) override
Stream an object of class RooFitResult.
StyleOption defaultPrintStyle(Option_t *opt) const override
Configure mapping of Print() arguments to RooPrintable print styles.
double covariance(Int_t row, Int_t col) const
Return the covariance matrix element addressed with numeric indices.
TMatrixF * _Lt
! triangular matrix used for generate random perturbations
const RooArgList & floatParsFinal() const
Return list of floating parameters after fit.
RooArgList * _randomPars
! List of floating parameters with most recent random perturbation applied
static RooFitResult * lastMinuitFit(const RooArgList &varList=RooArgList())
Import the results of the last fit performed by gMinuit, interpreting the fit parameters as the given...
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print fit result to stream 'os'.
TVectorD * _GC
Global correlation coefficients.
RooAbsPdf * createHessePdf(const RooArgSet &params) const
Return a p.d.f that represents the fit result as a multi-variate probability densisty function on the...
double _edm
Estimated distance to minimum.
void setFinalParList(const RooArgList &list)
Fill the list of final values of the floating parameters.
void printArgs(std::ostream &os) const override
Print arguments of fit result, i.e. the parameters of the fit.
void printName(std::ostream &os) const override
Print name of fit result.
const TMatrixDSym & correlationMatrix() const
Return correlation matrix ;.
Multivariate Gaussian p.d.f.
static void blockDecompose(const TMatrixD &input, const std::vector< int > &map1, const std::vector< int > &map2, TMatrixDSym &S11, TMatrixD &S12, TMatrixD &S21, TMatrixDSym &S22)
Block decomposition of covI according to given maps of observables.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
void addObject(TObject *obj, Option_t *drawOptions="", bool invisible=false)
Add a generic object to this plot.
Definition RooPlot.cxx:383
void addPlotable(RooPlotable *plotable, Option_t *drawOptions="", bool invisible=false, bool refreshNorm=false)
Add the specified plotable object to our plot.
Definition RooPlot.cxx:531
RooPlotable is a 'mix-in' base class that define the standard RooFit plotting and printing methods.
virtual StyleOption defaultPrintStyle(Option_t *opt) const
virtual void Streamer(TBuffer &)
virtual void printValue(std::ostream &os) const
Interface to print value of object.
static double gaussian(TRandom *generator=randomGenerator())
Return a Gaussian random variable with mean 0 and variance 1.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
double getError() const
Definition RooRealVar.h:58
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition TAxis.cxx:886
Create a Box.
Definition TBox.h:22
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual void RemoveAll(TCollection *col)
Remove all objects in collection col from this collection.
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
virtual TList * GetList() const
Definition TDirectory.h:222
The axis painter class.
Definition TGaxis.h:24
TAxis * GetXaxis()
Definition TH1.h:322
virtual void SetMaximum(Double_t maximum=-1111)
Definition TH1.h:400
TAxis * GetYaxis()
Definition TH1.h:323
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:401
2-D histogram with a double per channel (see TH1 documentation)}
Definition TH2.h:301
Service class for 2-D histogram classes.
Definition TH2.h:30
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:347
Use the TLine constructor to create a simple line.
Definition TLine.h:22
void Streamer(TBuffer &) override
Stream all objects in the collection to or from the I/O buffer.
Definition TList.cxx:1191
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:822
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:470
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
Manages Markers.
Definition TMarker.h:22
Int_t GetNrows() const
Int_t GetNoElements() const
Int_t GetNcols() const
Double_t Determinant() const override
TMatrixT.
Definition TMatrixT.h:39
@ kTransposed
Definition TMatrixT.h:58
Double_t * fU
Definition TMinuit.h:68
Int_t fNu
Definition TMinuit.h:130
Double_t * fGlobcc
Definition TMinuit.h:74
TString * fCpnam
Character to be plotted at the X,Y contour positions.
Definition TMinuit.h:165
Int_t * fNvarl
Definition TMinuit.h:126
Double_t * fMATUvline
Definition TMinuit.h:107
Double_t * fBlim
Definition TMinuit.h:70
Double_t * fWerr
Definition TMinuit.h:73
Int_t fStatus
Definition TMinuit.h:154
Int_t fNpar
Definition TMinuit.h:41
Double_t * fAlim
Definition TMinuit.h:69
Int_t * fNiofex
Definition TMinuit.h:127
Double_t * fVhmat
Definition TMinuit.h:89
virtual void mnstat(Double_t &fmin, Double_t &fedm, Double_t &errdef, Int_t &npari, Int_t &nparx, Int_t &istat)
Returns concerning the current status of the minimization.
Definition TMinuit.cxx:7637
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
void Streamer(TBuffer &) override
Stream an object of class TObject.
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition TNamed.cxx:154
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TObject.cxx:223
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
Basic string class.
Definition TString.h:139
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
Definition TString.cxx:450
const char * Data() const
Definition TString.h:380
void ToUpper()
Change string to upper case.
Definition TString.cxx:1183
TString & Append(const char *cs)
Definition TString.h:576
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
TLine * line
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4