52 RooHistPdf::RooHistPdf() : _dataHist(0), _totVolume(0), _unitNorm(kFALSE)
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),
71 _cdfBoundaries(kFALSE),
75 _histObsList.addClone(vars) ;
76 _pdfObsList.add(vars) ;
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 ;
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 ;
95 for (
const auto hobs : _histObsList) {
97 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
98 RooRealVar* dhreal =
dynamic_cast<RooRealVar*
>(dhobs) ;
100 ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
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),
123 _cdfBoundaries(kFALSE),
127 _histObsList.addClone(histObs) ;
128 _pdfObsList.add(pdfObs) ;
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")) ;
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")) ;
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.")) ;
153 for (
const auto hobs : _histObsList) {
155 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
156 RooRealVar* dhreal =
dynamic_cast<RooRealVar*
>(dhobs) ;
158 ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
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)
178 _histObsList.addClone(other._histObsList) ;
188 RooHistPdf::~RooHistPdf()
202 Double_t RooHistPdf::evaluate()
const
205 for (
unsigned int i=0; i < _pdfObsList.size(); ++i) {
206 RooAbsArg* harg = _histObsList[i];
207 RooAbsArg* parg = _pdfObsList[i];
211 harg->copyCache(parg,kTRUE) ;
212 if (!harg->inRange(0)) {
218 Double_t ret = _dataHist->weight(_histObsList, _intOrder, !_unitNorm, _cdfBoundaries);
235 Double_t RooHistPdf::totVolume()
const
243 for (
const auto arg : _histObsList) {
244 RooRealVar* real =
dynamic_cast<RooRealVar*
>(arg) ;
246 _totVolume *= (real->getMax()-real->getMin()) ;
248 RooCategory* cat =
dynamic_cast<RooCategory*
>(arg) ;
250 _totVolume *= cat->numTypes() ;
259 bool fullRange(
const RooAbsArg& x,
const RooAbsArg& y ,
const char* range)
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()) {
268 if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
270 return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
272 return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
284 Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars,
const char* rangeName)
const
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];
294 if (allVars.find(*pa)) {
297 if (fullRange(*pa, *ha, rangeName)) {
303 if (code == frcode) {
311 if (_intOrder > 1 && !(code & 1)) {
312 analVars.removeAll();
315 return (code >= 2) ? code : 0;
325 Double_t RooHistPdf::analyticalIntegral(Int_t code,
const char* rangeName)
const
328 if (((2 << _histObsList.getSize()) - 1) == code) {
329 return _dataHist->sum(kFALSE);
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];
340 if (code & (2 << n)) {
344 RooAbsRealLValue* rlv =
dynamic_cast<RooAbsRealLValue*
>(pa);
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());
355 ranges[ha] = std::make_pair(
356 binning->lowBoundFunc()->getVal(), binning->highBoundFunc()->getVal());
365 ha->copyCache(pa,kTRUE);
369 Double_t ret = (code & 1) ?
370 _dataHist->sum(intSet,_histObsList,kTRUE,kTRUE) :
371 _dataHist->sum(intSet,_histObsList,kFALSE,kTRUE, ranges);
388 list<Double_t>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi)
const
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()) ;
409 RooAbsLValue* lval =
dynamic_cast<RooAbsLValue*
>(dhObs) ;
416 const RooAbsBinning* binning = lval->getBinningPtr(0) ;
417 Double_t* boundaries = binning->array() ;
419 list<Double_t>* hint =
new list<Double_t> ;
422 xlo = xlo - 0.01*(xhi-xlo) ;
423 xhi = xhi + 0.01*(xhi-xlo) ;
425 Double_t delta = (xhi-xlo)*1e-8 ;
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) ;
446 std::list<Double_t>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi)
const
454 RooAbsLValue* lvarg =
dynamic_cast<RooAbsLValue*
>(_dataHist->get()->find(obs.GetName())) ;
460 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
461 Double_t* boundaries = binning->array() ;
463 list<Double_t>* hint =
new list<Double_t> ;
467 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
468 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
469 hint->push_back(boundaries[i]) ;
482 Int_t RooHistPdf::getMaxVal(
const RooArgSet& vars)
const
484 RooAbsCollection* common = _pdfObsList.selectCommon(vars) ;
485 if (common->getSize()==_pdfObsList.getSize()) {
496 Double_t RooHistPdf::maxVal(Int_t code)
const
501 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
503 Double_t wgt = _dataHist->weight() ;
504 if (wgt>max) max=wgt ;
515 Bool_t RooHistPdf::areIdentical(
const RooDataHist& dh1,
const RooDataHist& dh2)
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++) {
522 if (fabs(dh1.weight()-dh2.weight())>1e-8)
return kFALSE ;
532 Bool_t RooHistPdf::importWorkspaceHook(RooWorkspace& ws)
534 std::list<RooAbsData*> allData = ws.allData() ;
535 std::list<RooAbsData*>::const_iterator iter ;
536 for (iter = allData.begin() ; iter != allData.end() ; ++iter) {
538 if (*iter == _dataHist) {
544 RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName()) ;
549 if (wsdata->InheritsFrom(RooDataHist::Class())) {
552 if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
555 _dataHist = (RooDataHist*) wsdata ;
559 TString uniqueName = Form(
"%s_%s",_dataHist->GetName(),GetName()) ;
560 Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
562 coutE(ObjectHandling) <<
" RooHistPdf::importWorkspaceHook(" << GetName() <<
") unable to import clone of underlying RooDataHist with unique name " << uniqueName <<
", abort" << endl ;
565 _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
571 TString uniqueName = Form(
"%s_%s",_dataHist->GetName(),GetName()) ;
572 Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
574 coutE(ObjectHandling) <<
" RooHistPdf::importWorkspaceHook(" << GetName() <<
") unable to import clone of underlying RooDataHist with unique name " << uniqueName <<
", abort" << endl ;
577 _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
584 ws.import(*_dataHist,RooFit::Embedded()) ;
587 _dataHist = (RooDataHist*) ws.embeddedData(_dataHist->GetName()) ;
595 void RooHistPdf::Streamer(TBuffer &R__b)
597 if (R__b.IsReading()) {
598 R__b.ReadClassBuffer(RooHistPdf::Class(),
this);
603 R__b.WriteClassBuffer(RooHistPdf::Class(),
this);