Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooHistPdf.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofit:$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 \file RooHistPdf.cxx
19 \class RooHistPdf
20 \ingroup Roofitcore
21 
22 RooHistPdf implements a probablity density function sampled from a
23 multidimensional histogram. The histogram distribution is explicitly
24 normalized by RooHistPdf and can have an arbitrary number of real or
25 discrete dimensions.
26 **/
27 
28 #include "RooFit.h"
29 #include "Riostream.h"
30 
31 #include "RooHistPdf.h"
32 #include "RooDataHist.h"
33 #include "RooMsgService.h"
34 #include "RooRealVar.h"
35 #include "RooCategory.h"
36 #include "RooWorkspace.h"
37 #include "RooGlobalFunc.h"
38 
39 #include "TError.h"
40 
41 using namespace std;
42 
43 ClassImp(RooHistPdf);
44 ;
45 
46 
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 /// Default constructor
50 /// coverity[UNINIT_CTOR]
51 
52 RooHistPdf::RooHistPdf() : _dataHist(0), _totVolume(0), _unitNorm(kFALSE)
53 {
54 
55 }
56 
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 /// Constructor from a RooDataHist. RooDataHist dimensions
60 /// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
61 /// RooHistPdf neither owns or clone 'dhist' and the user must ensure the input histogram exists
62 /// for the entire life span of this PDF.
63 
64 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet& vars,
65  const RooDataHist& dhist, Int_t intOrder) :
66  RooAbsPdf(name,title),
67  _pdfObsList("pdfObs","List of p.d.f. observables",this),
68  _dataHist((RooDataHist*)&dhist),
69  _codeReg(10),
70  _intOrder(intOrder),
71  _cdfBoundaries(kFALSE),
72  _totVolume(0),
73  _unitNorm(kFALSE)
74 {
75  _histObsList.addClone(vars) ;
76  _pdfObsList.add(vars) ;
77 
78  // Verify that vars and dhist.get() have identical contents
79  const RooArgSet* dvars = dhist.get() ;
80  if (vars.getSize()!=dvars->getSize()) {
81  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
82  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
83  assert(0) ;
84  }
85  for (const auto arg : vars) {
86  if (!dvars->find(arg->GetName())) {
87  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
88  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
89  assert(0) ;
90  }
91  }
92 
93 
94  // Adjust ranges of _histObsList to those of _dataHist
95  for (const auto hobs : _histObsList) {
96  // Guaranteed to succeed, since checked above in ctor
97  RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
98  RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
99  if (dhreal){
100  ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
101  }
102  }
103 
104 }
105 
106 
107 
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Constructor from a RooDataHist. The first list of observables are the p.d.f.
111 /// observables, which may any RooAbsReal (function or variable). The second list
112 /// are the corresponding observables in the RooDataHist which must be of type
113 /// RooRealVar or RooCategory This constructor thus allows to apply a coordinate transformation
114 /// on the histogram data to be applied.
115 
116 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs,
117  const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder) :
118  RooAbsPdf(name,title),
119  _pdfObsList("pdfObs","List of p.d.f. observables",this),
120  _dataHist((RooDataHist*)&dhist),
121  _codeReg(10),
122  _intOrder(intOrder),
123  _cdfBoundaries(kFALSE),
124  _totVolume(0),
125  _unitNorm(kFALSE)
126 {
127  _histObsList.addClone(histObs) ;
128  _pdfObsList.add(pdfObs) ;
129 
130  // Verify that vars and dhist.get() have identical contents
131  const RooArgSet* dvars = dhist.get() ;
132  if (histObs.getSize()!=dvars->getSize()) {
133  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
134  << ") ERROR histogram variable list and RooDataHist must contain the same variables." << endl ;
135  throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
136  }
137 
138  for (const auto arg : histObs) {
139  if (!dvars->find(arg->GetName())) {
140  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
141  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
142  throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
143  }
144  if (!arg->isFundamental()) {
145  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
146  << ") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << endl ;
147  throw(string("RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
148  }
149  }
150 
151 
152  // Adjust ranges of _histObsList to those of _dataHist
153  for (const auto hobs : _histObsList) {
154  // Guaranteed to succeed, since checked above in ctor
155  RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
156  RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
157  if (dhreal){
158  ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
159  }
160  }
161 }
162 
163 
164 
165 ////////////////////////////////////////////////////////////////////////////////
166 /// Copy constructor
167 
168 RooHistPdf::RooHistPdf(const RooHistPdf& other, const char* name) :
169  RooAbsPdf(other,name),
170  _pdfObsList("pdfObs",this,other._pdfObsList),
171  _dataHist(other._dataHist),
172  _codeReg(other._codeReg),
173  _intOrder(other._intOrder),
174  _cdfBoundaries(other._cdfBoundaries),
175  _totVolume(other._totVolume),
176  _unitNorm(other._unitNorm)
177 {
178  _histObsList.addClone(other._histObsList) ;
179 
180 }
181 
182 
183 
184 
185 ////////////////////////////////////////////////////////////////////////////////
186 /// Destructor
187 
188 RooHistPdf::~RooHistPdf()
189 {
190 
191 }
192 
193 
194 
195 
196 
197 ////////////////////////////////////////////////////////////////////////////////
198 /// Return the current value: The value of the bin enclosing the current coordinates
199 /// of the observables, normalized by the histograms contents. Interpolation
200 /// is applied if the RooHistPdf is configured to do that.
201 
202 Double_t RooHistPdf::evaluate() const
203 {
204  // Transfer values from
205  for (unsigned int i=0; i < _pdfObsList.size(); ++i) {
206  RooAbsArg* harg = _histObsList[i];
207  RooAbsArg* parg = _pdfObsList[i];
208 
209  if (harg != parg) {
210  parg->syncCache() ;
211  harg->copyCache(parg,kTRUE) ;
212  if (!harg->inRange(0)) {
213  return 0 ;
214  }
215  }
216  }
217 
218  Double_t ret = _dataHist->weight(_histObsList, _intOrder, !_unitNorm, _cdfBoundaries);
219 // cout << "RooHistPdf::evaluate(" << GetName() << ") ret = " << ret << " ";
220 // cout << _histObsList[0] << " ";
221 // _histObsList[0]->Print("");
222 // _dataHist->Print("V");
223 // _dataHist->dump2();
224 
225  if (ret<0) {
226  ret=0 ;
227  }
228  return ret ;
229 }
230 
231 
232 ////////////////////////////////////////////////////////////////////////////////
233 /// Return the total volume spanned by the observables of the RooHistPdf
234 
235 Double_t RooHistPdf::totVolume() const
236 {
237  // Return previously calculated value, if any
238  if (_totVolume>0) {
239  return _totVolume ;
240  }
241  _totVolume = 1. ;
242 
243  for (const auto arg : _histObsList) {
244  RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
245  if (real) {
246  _totVolume *= (real->getMax()-real->getMin()) ;
247  } else {
248  RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
249  if (cat) {
250  _totVolume *= cat->numTypes() ;
251  }
252  }
253  }
254 
255  return _totVolume ;
256 }
257 
258 namespace {
259 bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
260 {
261  const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
262  const RooAbsRealLValue *_y = dynamic_cast<const RooAbsRealLValue*>(&y);
263  if (!_x || !_y) return false;
264  if (!range || !strlen(range) || !_x->hasRange(range) ||
265  _x->getBinningPtr(range)->isParameterized()) {
266  // parameterized ranges may be full range now, but that might change,
267  // so return false
268  if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
269  return false;
270  return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
271  }
272  return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
273 }
274 }
275 
276 
277 ////////////////////////////////////////////////////////////////////////////////
278 /// Determine integration scenario. If no interpolation is used,
279 /// RooHistPdf can perform all integrals over its dependents
280 /// analytically via partial or complete summation of the input
281 /// histogram. If interpolation is used on the integral over
282 /// all histogram observables is supported
283 
284 Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
285 {
286  // First make list of pdf observables to histogram observables
287  // and select only those for which the integral is over the full range
288 
289  Int_t code = 0, frcode = 0;
290  for (unsigned int n=0; n < _pdfObsList.size() && n < _histObsList.size(); ++n) {
291  const auto pa = _pdfObsList[n];
292  const auto ha = _histObsList[n];
293 
294  if (allVars.find(*pa)) {
295  code |= 2 << n;
296  analVars.add(*pa);
297  if (fullRange(*pa, *ha, rangeName)) {
298  frcode |= 2 << n;
299  }
300  }
301  }
302 
303  if (code == frcode) {
304  // integrate over full range of all observables - use bit 0 to indicate
305  // full range integration over all observables
306  code |= 1;
307  }
308  // Disable partial analytical integrals if interpolation is used, and we
309  // integrate over sub-ranges, but leave them enabled when we integrate over
310  // the full range of one or several variables
311  if (_intOrder > 1 && !(code & 1)) {
312  analVars.removeAll();
313  return 0;
314  }
315  return (code >= 2) ? code : 0;
316 }
317 
318 
319 
320 ////////////////////////////////////////////////////////////////////////////////
321 /// Return integral identified by 'code'. The actual integration
322 /// is deferred to RooDataHist::sum() which implements partial
323 /// or complete summation over the histograms contents.
324 
325 Double_t RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
326 {
327  // Simplest scenario, full-range integration over all dependents
328  if (((2 << _histObsList.getSize()) - 1) == code) {
329  return _dataHist->sum(kFALSE);
330  }
331 
332  // Partial integration scenario, retrieve set of variables, calculate partial
333  // sum, figure out integration ranges (if needed)
334  RooArgSet intSet;
335  std::map<const RooAbsArg*, std::pair<Double_t, Double_t> > ranges;
336  for (unsigned int n=0; n < _pdfObsList.size() && n < _histObsList.size(); ++n) {
337  const auto pa = _pdfObsList[n];
338  const auto ha = _histObsList[n];
339 
340  if (code & (2 << n)) {
341  intSet.add(*ha);
342  }
343  if (!(code & 1)) {
344  RooAbsRealLValue* rlv = dynamic_cast<RooAbsRealLValue*>(pa);
345  if (rlv) {
346  const RooAbsBinning* binning = rlv->getBinningPtr(rangeName);
347  if (rangeName && rlv->hasRange(rangeName)) {
348  ranges[ha] = std::make_pair(
349  rlv->getMin(rangeName), rlv->getMax(rangeName));
350  } else if (binning) {
351  if (!binning->isParameterized()) {
352  ranges[ha] = std::make_pair(
353  binning->lowBound(), binning->highBound());
354  } else {
355  ranges[ha] = std::make_pair(
356  binning->lowBoundFunc()->getVal(), binning->highBoundFunc()->getVal());
357  }
358  }
359  }
360  }
361  // WVE must sync hist slice list values to pdf slice list
362  // Transfer values from
363  if (ha != pa) {
364  pa->syncCache();
365  ha->copyCache(pa,kTRUE);
366  }
367  }
368 
369  Double_t ret = (code & 1) ?
370  _dataHist->sum(intSet,_histObsList,kTRUE,kTRUE) :
371  _dataHist->sum(intSet,_histObsList,kFALSE,kTRUE, ranges);
372 
373  // cout << "intSet = " << intSet << endl ;
374  // cout << "slice position = " << endl ;
375  // _histObsList.Print("v") ;
376  // cout << "RooHistPdf::ai(" << GetName() << ") code = " << code << " ret = " << ret << endl ;
377 
378  return ret ;
379 }
380 
381 
382 
383 ////////////////////////////////////////////////////////////////////////////////
384 /// Return sampling hint for making curves of (projections) of this function
385 /// as the recursive division strategy of RooCurve cannot deal efficiently
386 /// with the vertical lines that occur in a non-interpolated histogram
387 
388 list<Double_t>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
389 {
390  // No hints are required when interpolation is used
391  if (_intOrder>0) {
392  return 0 ;
393  }
394 
395  // Check that observable is in dataset, if not no hint is generated
396  RooAbsArg* dhObs = nullptr;
397  for (unsigned int i=0; i < _pdfObsList.size(); ++i) {
398  RooAbsArg* histObs = _histObsList[i];
399  RooAbsArg* pdfObs = _pdfObsList[i];
400  if (TString(obs.GetName())==pdfObs->GetName()) {
401  dhObs = _dataHist->get()->find(histObs->GetName()) ;
402  break;
403  }
404  }
405 
406  if (!dhObs) {
407  return 0 ;
408  }
409  RooAbsLValue* lval = dynamic_cast<RooAbsLValue*>(dhObs) ;
410  if (!lval) {
411  return 0 ;
412  }
413 
414  // Retrieve position of all bin boundaries
415 
416  const RooAbsBinning* binning = lval->getBinningPtr(0) ;
417  Double_t* boundaries = binning->array() ;
418 
419  list<Double_t>* hint = new list<Double_t> ;
420 
421  // Widen range slighty
422  xlo = xlo - 0.01*(xhi-xlo) ;
423  xhi = xhi + 0.01*(xhi-xlo) ;
424 
425  Double_t delta = (xhi-xlo)*1e-8 ;
426 
427  // Construct array with pairs of points positioned epsilon to the left and
428  // right of the bin boundaries
429  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
430  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
431  hint->push_back(boundaries[i]-delta) ;
432  hint->push_back(boundaries[i]+delta) ;
433  }
434  }
435 
436  return hint ;
437 }
438 
439 
440 
441 ////////////////////////////////////////////////////////////////////////////////
442 /// Return sampling hint for making curves of (projections) of this function
443 /// as the recursive division strategy of RooCurve cannot deal efficiently
444 /// with the vertical lines that occur in a non-interpolated histogram
445 
446 std::list<Double_t>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
447 {
448  // No hints are required when interpolation is used
449  if (_intOrder>0) {
450  return 0 ;
451  }
452 
453  // Check that observable is in dataset, if not no hint is generated
454  RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
455  if (!lvarg) {
456  return 0 ;
457  }
458 
459  // Retrieve position of all bin boundaries
460  const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
461  Double_t* boundaries = binning->array() ;
462 
463  list<Double_t>* hint = new list<Double_t> ;
464 
465  // Construct array with pairs of points positioned epsilon to the left and
466  // right of the bin boundaries
467  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
468  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
469  hint->push_back(boundaries[i]) ;
470  }
471  }
472 
473  return hint ;
474 }
475 
476 
477 
478 
479 ////////////////////////////////////////////////////////////////////////////////
480 /// Only handle case of maximum in all variables
481 
482 Int_t RooHistPdf::getMaxVal(const RooArgSet& vars) const
483 {
484  RooAbsCollection* common = _pdfObsList.selectCommon(vars) ;
485  if (common->getSize()==_pdfObsList.getSize()) {
486  delete common ;
487  return 1;
488  }
489  delete common ;
490  return 0 ;
491 }
492 
493 
494 ////////////////////////////////////////////////////////////////////////////////
495 
496 Double_t RooHistPdf::maxVal(Int_t code) const
497 {
498  R__ASSERT(code==1) ;
499 
500  Double_t max(-1) ;
501  for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
502  _dataHist->get(i) ;
503  Double_t wgt = _dataHist->weight() ;
504  if (wgt>max) max=wgt ;
505  }
506 
507  return max*1.05 ;
508 }
509 
510 
511 
512 
513 ////////////////////////////////////////////////////////////////////////////////
514 
515 Bool_t RooHistPdf::areIdentical(const RooDataHist& dh1, const RooDataHist& dh2)
516 {
517  if (fabs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return kFALSE ;
518  if (dh1.numEntries() != dh2.numEntries()) return kFALSE ;
519  for (int i=0 ; i < dh1.numEntries() ; i++) {
520  dh1.get(i) ;
521  dh2.get(i) ;
522  if (fabs(dh1.weight()-dh2.weight())>1e-8) return kFALSE ;
523  }
524  return kTRUE ;
525 }
526 
527 
528 
529 ////////////////////////////////////////////////////////////////////////////////
530 /// Check if our datahist is already in the workspace
531 
532 Bool_t RooHistPdf::importWorkspaceHook(RooWorkspace& ws)
533 {
534  std::list<RooAbsData*> allData = ws.allData() ;
535  std::list<RooAbsData*>::const_iterator iter ;
536  for (iter = allData.begin() ; iter != allData.end() ; ++iter) {
537  // If your dataset is already in this workspace nothing needs to be done
538  if (*iter == _dataHist) {
539  return kFALSE ;
540  }
541  }
542 
543  // Check if dataset with given name already exists
544  RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName()) ;
545 
546  if (wsdata) {
547 
548  // Yes it exists - now check if it is identical to our internal histogram
549  if (wsdata->InheritsFrom(RooDataHist::Class())) {
550 
551  // Check if histograms are identical
552  if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
553 
554  // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
555  _dataHist = (RooDataHist*) wsdata ;
556  } else {
557 
558  // not identical, clone rename and import
559  TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
560  Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
561  if (flag) {
562  coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
563  return kTRUE ;
564  }
565  _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
566  }
567 
568  } else {
569 
570  // Exists and is NOT of correct type: clone rename and import
571  TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
572  Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
573  if (flag) {
574  coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
575  return kTRUE ;
576  }
577  _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
578 
579  }
580  return kFALSE ;
581  }
582 
583  // We need to import our datahist into the workspace
584  ws.import(*_dataHist,RooFit::Embedded()) ;
585 
586  // Redirect our internal pointer to the copy in the workspace
587  _dataHist = (RooDataHist*) ws.embeddedData(_dataHist->GetName()) ;
588  return kFALSE ;
589 }
590 
591 
592 ////////////////////////////////////////////////////////////////////////////////
593 /// Stream an object of class RooHistPdf.
594 
595 void RooHistPdf::Streamer(TBuffer &R__b)
596 {
597  if (R__b.IsReading()) {
598  R__b.ReadClassBuffer(RooHistPdf::Class(),this);
599  // WVE - interim solution - fix proxies here
600  //_proxyList.Clear() ;
601  //registerProxy(_pdfObsList) ;
602  } else {
603  R__b.WriteClassBuffer(RooHistPdf::Class(),this);
604  }
605 }
606