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