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