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