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