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