Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooProdPdf.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\file RooProdPdf.cxx
19\class RooProdPdf
20\ingroup Roofitcore
21
22RooProdPdf is an efficient implementation of a product of PDFs of the form
23\f[ \prod_{i=1}^{N} \mathrm{PDF}_i (x, \ldots) \f]
24
25PDFs may share observables. If that is the case any irreducible subset
26of PDFs that share observables will be normalised with explicit numeric
27integration as any built-in normalisation will no longer be valid.
28
29Alternatively, products using conditional PDFs can be defined, *e.g.*
30
31\f[ F(x|y) \cdot G(y), \f]
32
33meaning a PDF \f$ F(x) \f$ **given** \f$ y \f$ and a PDF \f$ G(y) \f$.
34In this construction, \f$ F \f$ is only
35normalised w.r.t \f$ x\f$, and \f$ G \f$ is normalised w.r.t \f$ y \f$. The product in this construction
36is properly normalised.
37
38If exactly one of the component PDFs supports extended likelihood fits, the
39product will also be usable in extended mode, returning the number of expected
40events from the extendable component PDF. The extendable component does not
41have to appear in any specific place in the list.
42**/
43
44#include "RooProdPdf.h"
45#include "RooBatchCompute.h"
46#include "RooRealProxy.h"
47#include "RooProdGenContext.h"
48#include "RooGenProdProj.h"
49#include "RooProduct.h"
50#include "RooNameReg.h"
51#include "RooMsgService.h"
52#include "RooFormulaVar.h"
53#include "RooRealVar.h"
54#include "RooAddition.h"
55#include "RooGlobalFunc.h"
56#include "RooConstVar.h"
57#include "RooWorkspace.h"
58#include "RooRangeBoolean.h"
59#include "RooCustomizer.h"
60#include "RooRealIntegral.h"
61#include "RooTrace.h"
62#include "RooHelpers.h"
63#include "RooBatchCompute.h"
64#include "strtok.h"
65
66#include <cstring>
67#include <sstream>
68#include <algorithm>
69
70#ifndef _WIN32
71#include <strings.h>
72#endif
73
74using namespace std;
75
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// Default constructor
81
83 _cacheMgr(this,10)
84{
85 // Default constructor
87}
88
89
90////////////////////////////////////////////////////////////////////////////////
91/// Constructor with 2 PDFs (most frequent use case).
92///
93/// The optional cutOff parameter can be used as a speed optimization if
94/// one or more of the PDF have sizable regions with very small values,
95/// which would pull the entire product of PDFs to zero in those regions.
96///
97/// After each PDF multiplication, the running product is compared with
98/// the cutOff parameter. If the running product is smaller than the
99/// cutOff value, the product series is terminated and remaining PDFs
100/// are not evaluated.
101///
102/// There is no magic value of the cutOff, the user should experiment
103/// to find the appropriate balance between speed and precision.
104/// If a cutoff is specified, the PDFs most likely to be small should
105/// be put first in the product. The default cutOff value is zero.
106///
107
108RooProdPdf::RooProdPdf(const char *name, const char *title,
109 RooAbsPdf& pdf1, RooAbsPdf& pdf2, double cutOff) :
110 RooAbsPdf(name,title),
111 _cacheMgr(this,10),
112 _cutOff(cutOff),
113 _pdfList("!pdfs","List of PDFs",this)
114{
115 _pdfList.add(pdf1) ;
116 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>("nset")) ;
117 if (pdf1.canBeExtended()) {
118 _extendedIndex = _pdfList.index(&pdf1) ;
119 }
120
121 _pdfList.add(pdf2) ;
122 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>("nset")) ;
123
124 if (pdf2.canBeExtended()) {
125 if (_extendedIndex>=0) {
126 // Protect against multiple extended terms
127 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName()
128 << ") multiple components with extended terms detected,"
129 << " product will not be extendible." << endl ;
130 _extendedIndex=-1 ;
131 } else {
133 }
134 }
136}
137
138
139
140////////////////////////////////////////////////////////////////////////////////
141/// Constructor from a list of PDFs.
142///
143/// The optional cutOff parameter can be used as a speed optimization if
144/// one or more of the PDF have sizable regions with very small values,
145/// which would pull the entire product of PDFs to zero in those regions.
146///
147/// After each PDF multiplication, the running product is compared with
148/// the cutOff parameter. If the running product is smaller than the
149/// cutOff value, the product series is terminated and remaining PDFs
150/// are not evaluated.
151///
152/// There is no magic value of the cutOff, the user should experiment
153/// to find the appropriate balance between speed and precision.
154/// If a cutoff is specified, the PDFs most likely to be small should
155/// be put first in the product. The default cutOff value is zero.
156
157RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgList& inPdfList, double cutOff) :
158 RooAbsPdf(name,title),
159 _cacheMgr(this,10),
160 _cutOff(cutOff),
161 _pdfList("!pdfs","List of PDFs",this)
162{
163 addPdfs(inPdfList);
165}
166
167
168
169////////////////////////////////////////////////////////////////////////////////
170/// Constructor from named argument list.
171/// \param[in] name Name used by RooFit
172/// \param[in] title Title used for plotting
173/// \param[in] fullPdfSet Set of "regular" PDFs that are normalised over all their observables
174/// \param[in] arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8 Optional arguments according to table below.
175///
176/// <table>
177/// <tr><th> Argument <th> Description
178/// <tr><td> `Conditional(pdfSet,depSet,depsAreCond=false)` <td> Add PDF to product with condition that it
179/// only be normalized over specified observables. Any remaining observables will be conditional observables.
180/// (Setting `depsAreCond` to true inverts this, so the observables in depSet will be the conditional observables.)
181/// </table>
182///
183/// For example, given a PDF \f$ F(x,y) \f$ and \f$ G(y) \f$,
184///
185/// `RooProdPdf("P", "P", G, Conditional(F,x))` will construct a 2-dimensional PDF as follows:
186/// \f[
187/// P(x,y) = \frac{G(y)}{\int_y G(y)} \cdot \frac{F(x,y)}{\int_x F(x,y)},
188/// \f]
189///
190/// which is a well normalised and properly defined PDF, but different from
191/// \f[
192/// P'(x,y) = \frac{F(x,y) \cdot G(y)}{\int_x\int_y F(x,y) \cdot G(y)}.
193/// \f]
194///
195/// In the former case, the \f$ y \f$ distribution of \f$ P \f$ is identical to that of \f$ G \f$, while
196/// \f$ F \f$ only is used to determine the correlation between \f$ X \f$ and \f$ Y \f$. In the latter
197/// case, the \f$ Y \f$ distribution is defined by the product of \f$ F \f$ and \f$ G \f$.
198///
199/// This \f$ P(x,y) \f$ construction is analoguous to generating events from \f$ F(x,y) \f$ with
200/// a prototype dataset sampled from \f$ G(y) \f$.
201
202RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgSet& fullPdfSet,
203 const RooCmdArg& arg1, const RooCmdArg& arg2,
204 const RooCmdArg& arg3, const RooCmdArg& arg4,
205 const RooCmdArg& arg5, const RooCmdArg& arg6,
206 const RooCmdArg& arg7, const RooCmdArg& arg8) :
207 RooAbsPdf(name,title),
208 _cacheMgr(this,10),
209 _pdfList("!pdfs","List of PDFs",this)
210{
212 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
213 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
214 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
215 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
216
217 initializeFromCmdArgList(fullPdfSet,l) ;
219}
220
221
222
223////////////////////////////////////////////////////////////////////////////////
224/// Constructor from named argument list
225
226RooProdPdf::RooProdPdf(const char* name, const char* title,
227 const RooCmdArg& arg1, const RooCmdArg& arg2,
228 const RooCmdArg& arg3, const RooCmdArg& arg4,
229 const RooCmdArg& arg5, const RooCmdArg& arg6,
230 const RooCmdArg& arg7, const RooCmdArg& arg8) :
231 RooAbsPdf(name,title),
232 _cacheMgr(this,10),
233 _pdfList("!pdfList","List of PDFs",this)
234{
236 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
237 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
238 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
239 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
240
243}
244
245
246
247////////////////////////////////////////////////////////////////////////////////
248/// Internal constructor from list of named arguments
249
250RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgSet& fullPdfSet, const RooLinkedList& cmdArgList) :
251 RooAbsPdf(name,title),
252 _cacheMgr(this,10),
253 _pdfList("!pdfs","List of PDFs",this)
254{
255 initializeFromCmdArgList(fullPdfSet, cmdArgList) ;
257}
258
259
260
261////////////////////////////////////////////////////////////////////////////////
262/// Copy constructor
263
264RooProdPdf::RooProdPdf(const RooProdPdf& other, const char* name) :
265 RooAbsPdf(other,name),
266 _cacheMgr(other._cacheMgr,this),
267 _genCode(other._genCode),
268 _cutOff(other._cutOff),
269 _pdfList("!pdfs",this,other._pdfList),
270 _extendedIndex(other._extendedIndex),
271 _useDefaultGen(other._useDefaultGen),
272 _refRangeName(other._refRangeName),
273 _selfNorm(other._selfNorm),
274 _defNormSet(other._defNormSet)
275{
276 // Clone contents of normalizarion set list
277 for(auto const& nset : other._pdfNSetList) {
278 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>(nset->GetName()));
279 nset->snapshot(*_pdfNSetList.back());
280 }
282}
283
284
285
286////////////////////////////////////////////////////////////////////////////////
287/// Initialize RooProdPdf configuration from given list of RooCmdArg configuration arguments
288/// and set of 'regular' p.d.f.s in product
289
291{
292 Int_t numExtended(0) ;
293
294 // Process set of full PDFS
295 for(auto const* pdf : static_range_cast<RooAbsPdf*>(fullPdfSet)) {
296 _pdfList.add(*pdf) ;
297 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>("nset")) ;
298
299 if (pdf->canBeExtended()) {
301 numExtended++ ;
302 }
303
304 }
305
306 // Process list of conditional PDFs
307 for(auto * carg : static_range_cast<RooCmdArg*>(l)) {
308
309 if (0 == strcmp(carg->GetName(), "Conditional")) {
310
311 Int_t argType = carg->getInt(0) ;
312 auto pdfSet = static_cast<RooArgSet const*>(carg->getSet(0));
313 auto normSet = static_cast<RooArgSet const*>(carg->getSet(1));
314
315 for(auto * thePdf : static_range_cast<RooAbsPdf*>(*pdfSet)) {
316 _pdfList.add(*thePdf) ;
317
318 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>(0 == argType ? "nset" : "cset"));
319 normSet->snapshot(*_pdfNSetList.back());
320
321 if (thePdf->canBeExtended()) {
322 _extendedIndex = _pdfList.index(thePdf) ;
323 numExtended++ ;
324 }
325
326 }
327
328 } else if (0 != strlen(carg->GetName())) {
329 coutW(InputArguments) << "Unknown arg: " << carg->GetName() << endl ;
330 }
331 }
332
333 // Protect against multiple extended terms
334 if (numExtended>1) {
335 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName()
336 << ") WARNING: multiple components with extended terms detected,"
337 << " product will not be extendible." << endl ;
338 _extendedIndex = -1 ;
339 }
340
341
342}
343
344
345
346////////////////////////////////////////////////////////////////////////////////
347/// Destructor
348
350{
352}
353
354
356 int code ;
357 auto cache = static_cast<CacheElem*>(_cacheMgr.getObj(nset, 0, &code)) ;
358
359 // If cache doesn't have our configuration, recalculate here
360 if (!cache) {
361 code = getPartIntList(nset, nullptr) ;
362 cache = static_cast<CacheElem*>(_cacheMgr.getObj(nset, 0, &code)) ;
363 }
364 return cache;
365}
366
367
368////////////////////////////////////////////////////////////////////////////////
369/// Calculate current value of object
370
372{
374}
375
376
377
378////////////////////////////////////////////////////////////////////////////////
379/// Calculate running product of pdfs terms, using the supplied
380/// normalization set in 'normSetList' for each component
381
382double RooProdPdf::calculate(const RooProdPdf::CacheElem& cache, bool /*verbose*/) const
383{
384 if (cache._isRearranged) {
385 if (dologD(Eval)) {
386 cxcoutD(Eval) << "RooProdPdf::calculate(" << GetName() << ") rearranged product calculation"
387 << " calculate: num = " << cache._rearrangedNum->GetName() << " = " << cache._rearrangedNum->getVal() << endl ;
388// cache._rearrangedNum->printComponentTree("",0,5) ;
389 cxcoutD(Eval) << "calculate: den = " << cache._rearrangedDen->GetName() << " = " << cache._rearrangedDen->getVal() << endl ;
390// cache._rearrangedDen->printComponentTree("",0,5) ;
391 }
392
393 return cache._rearrangedNum->getVal() / cache._rearrangedDen->getVal();
394 } else {
395
396 double value = 1.0;
397 assert(cache._normList.size() == cache._partList.size());
398 for (std::size_t i = 0; i < cache._partList.size(); ++i) {
399 const auto& partInt = static_cast<const RooAbsReal&>(cache._partList[i]);
400 const auto normSet = cache._normList[i].get();
401
402 const double piVal = partInt.getVal(normSet->getSize() > 0 ? normSet : nullptr);
403 value *= piVal ;
404 if (value <= _cutOff) break;
405 }
406
407 return value ;
408 }
409}
410
411////////////////////////////////////////////////////////////////////////////////
412/// Evaluate product of PDFs in batch mode.
413void RooProdPdf::calculateBatch(const RooProdPdf::CacheElem& cache, cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
414{
416
417 if (cache._isRearranged) {
418 auto numerator = dataMap.at(cache._rearrangedNum.get());
419 auto denominator = dataMap.at(cache._rearrangedDen.get());
420 dispatch->compute(stream, RooBatchCompute::Ratio, output, nEvents, {numerator, denominator});
421 } else {
423 factors.reserve(cache._partList.size());
424 for (const RooAbsArg *i : cache._partList) {
425 auto span = dataMap.at(i);
426 factors.push_back(span);
427 }
428 RooBatchCompute::ArgVector special{static_cast<double>(factors.size())};
429 dispatch->compute(stream, RooBatchCompute::ProdPdf, output, nEvents, factors, special);
430 }
431}
432
433namespace {
434
435template<class T>
436void eraseNullptrs(std::vector<T*>& v) {
437 v.erase(std::remove_if(v.begin(), v.end(), [](T* x){ return x == nullptr; } ), v.end());
438}
439
440void removeCommon(std::vector<RooAbsArg*> &v, std::span<RooAbsArg * const> other) {
441
442 for (auto const& arg : other) {
443 auto namePtrMatch = [&arg](const RooAbsArg* elm) {
444 return elm != nullptr && elm->namePtr() == arg->namePtr();
445 };
446
447 auto found = std::find_if(v.begin(), v.end(), namePtrMatch);
448 if(found != v.end()) {
449 *found = nullptr;
450 }
451 }
452 eraseNullptrs(v);
453}
454
455void addCommon(std::vector<RooAbsArg*> &v, std::vector<RooAbsArg*> const& o1, std::vector<RooAbsArg*> const& o2) {
456
457 for (auto const& arg : o1) {
458 auto namePtrMatch = [&arg](const RooAbsArg* elm) {
459 return elm->namePtr() == arg->namePtr();
460 };
461
462 if(std::find_if(o2.begin(), o2.end(), namePtrMatch) != o2.end()) {
463 v.push_back(arg);
464 }
465 }
466}
467
468}
469
470
471////////////////////////////////////////////////////////////////////////////////
472/// Factorize product in irreducible terms for given choice of integration/normalization
473
474void RooProdPdf::factorizeProduct(const RooArgSet& normSet, const RooArgSet& intSet,
475 RooLinkedList& termList, RooLinkedList& normList,
476 RooLinkedList& impDepList, RooLinkedList& crossDepList,
477 RooLinkedList& intList) const
478{
479 // List of all term dependents: normalization and imported
480 std::vector<RooArgSet> depAllList;
481 std::vector<RooArgSet> depIntNoNormList;
482
483 // Setup lists for factorization terms and their dependents
484 RooArgSet* term(0);
485 RooArgSet* termNormDeps(0);
486 RooArgSet* termIntDeps(0);
487 RooArgSet* termIntNoNormDeps(0);
488
489 std::vector<RooAbsArg*> pdfIntNoNormDeps;
490 std::vector<RooAbsArg*> pdfIntSet;
491 std::vector<RooAbsArg*> pdfNSet;
492 std::vector<RooAbsArg*> pdfCSet;
493 std::vector<RooAbsArg*> pdfNormDeps; // Dependents to be normalized for the PDF
494 std::vector<RooAbsArg*> pdfAllDeps; // All dependents of this PDF
495
496 // Loop over the PDFs
497 for(std::size_t iPdf = 0; iPdf < _pdfList.size(); ++iPdf) {
498 RooAbsPdf& pdf = static_cast<RooAbsPdf&>(_pdfList[iPdf]);
499 RooArgSet& pdfNSetOrig = *_pdfNSetList[iPdf];
500
501 pdfNSet.clear();
502 pdfCSet.clear();
503
504 // Make iterator over tree leaf node list to get the observables.
505 // This code is borrowed from RooAgsPdf::getObservables.
506 // RooAbsArg::treeNodeServer list is relatively expensive, so we only do it
507 // once and use it in a lambda function.
508 RooArgSet pdfLeafList("leafNodeServerList") ;
509 pdf.treeNodeServerList(&pdfLeafList,0,false,true,true) ;
510 auto getObservablesOfCurrentPdf = [&pdfLeafList](
511 std::vector<RooAbsArg*> & out,
512 const RooArgSet& dataList) {
513 for (const auto arg : pdfLeafList) {
514 if (arg->dependsOnValue(dataList) && arg->isLValue()) {
515 out.push_back(arg) ;
516 }
517 }
518 };
519
520 // Reduce pdfNSet to actual dependents
521 if (0 == strcmp("cset", pdfNSetOrig.GetName())) {
522 getObservablesOfCurrentPdf(pdfNSet, normSet);
523 removeCommon(pdfNSet, pdfNSetOrig.get());
524 pdfCSet = pdfNSetOrig.get();
525 } else {
526 // Interpret at NSet
527 getObservablesOfCurrentPdf(pdfNSet, pdfNSetOrig);
528 }
529
530
531 pdfNormDeps.clear();
532 pdfAllDeps.clear();
533
534 // Make list of all dependents of this PDF
535 getObservablesOfCurrentPdf(pdfAllDeps, normSet);
536
537
538// cout << GetName() << ": pdf = " << pdf->GetName() << " pdfAllDeps = " << pdfAllDeps << " pdfNSet = " << *pdfNSet << " pdfCSet = " << *pdfCSet << endl;
539
540 // Make list of normalization dependents for this PDF;
541 if (!pdfNSet.empty()) {
542 // PDF is conditional
543 addCommon(pdfNormDeps, pdfAllDeps, pdfNSet);
544 } else {
545 // PDF is regular
546 pdfNormDeps = pdfAllDeps;
547 }
548
549// cout << GetName() << ": pdfNormDeps for " << pdf->GetName() << " = " << pdfNormDeps << endl;
550
551 pdfIntSet.clear();
552 getObservablesOfCurrentPdf(pdfIntSet, intSet) ;
553
554 // WVE if we have no norm deps, conditional observables should be taken out of pdfIntSet
555 if (pdfNormDeps.empty() && !pdfCSet.empty()) {
556 removeCommon(pdfIntSet, pdfCSet);
557// cout << GetName() << ": have no norm deps, removing conditional observables from intset" << endl;
558 }
559
560 pdfIntNoNormDeps.clear();
561 pdfIntNoNormDeps = pdfIntSet;
562 removeCommon(pdfIntNoNormDeps, pdfNormDeps);
563
564// cout << GetName() << ": pdf = " << pdf->GetName() << " intset = " << *pdfIntSet << " pdfIntNoNormDeps = " << pdfIntNoNormDeps << endl;
565
566 // Check if this PDF has dependents overlapping with one of the existing terms
567 bool done = false;
568 int j = 0;
569 auto lIter = termList.begin();
570 auto ldIter = normList.begin();
571 for(;lIter != termList.end(); (++lIter, ++ldIter, ++j)) {
572 termNormDeps = static_cast<RooArgSet*>(*ldIter);
573 term = static_cast<RooArgSet*>(*lIter);
574 // PDF should be added to existing term if
575 // 1) It has overlapping normalization dependents with any other PDF in existing term
576 // 2) It has overlapping dependents of any class for which integration is requested
577 // 3) If normalization happens over multiple ranges, and those ranges are both defined
578 // in either observable
579
580 bool normOverlap = termNormDeps->overlaps(pdfNormDeps.begin(), pdfNormDeps.end());
581 //bool intOverlap = pdfIntSet->overlaps(*termAllDeps);
582
583 if (normOverlap) {
584// cout << GetName() << ": this term overlaps with term " << (*term) << " in normalization observables" << endl;
585
586 term->add(pdf);
587 termNormDeps->add(pdfNormDeps.begin(), pdfNormDeps.end(), false);
588 depAllList[j].add(pdfAllDeps.begin(), pdfAllDeps.end(), false);
589 if (termIntDeps) {
590 termIntDeps->add(pdfIntSet.begin(), pdfIntSet.end(), false);
591 }
592 if (termIntNoNormDeps) {
593 termIntNoNormDeps->add(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), false);
594 }
595 termIntNoNormDeps->add(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), false);
596 done = true;
597 }
598 }
599
600 // If not, create a new term
601 if (!done) {
602 if (!(pdfNormDeps.empty() && pdfAllDeps.empty() &&
603 pdfIntSet.empty()) || normSet.empty()) {
604 term = new RooArgSet("term");
605 termNormDeps = new RooArgSet("termNormDeps");
606 depAllList.emplace_back(pdfAllDeps.begin(), pdfAllDeps.end(), "termAllDeps");
607 termIntDeps = new RooArgSet(pdfIntSet.begin(), pdfIntSet.end(), "termIntDeps");
608 depIntNoNormList.emplace_back(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), "termIntNoNormDeps");
609 termIntNoNormDeps = &depIntNoNormList.back();
610
611 term->add(pdf);
612 termNormDeps->add(pdfNormDeps.begin(), pdfNormDeps.end(), false);
613
614 termList.Add(term);
615 normList.Add(termNormDeps);
616 intList.Add(termIntDeps);
617 }
618 }
619
620 }
621
622 // Loop over list of terms again to determine 'imported' observables
623 int i = 0;
624 RooArgSet *normDeps;
625 auto lIter = termList.begin();
626 auto ldIter = normList.begin();
627 for(;lIter != termList.end(); (++lIter, ++ldIter, ++i)) {
628 normDeps = static_cast<RooArgSet*>(*ldIter);
629 term = static_cast<RooArgSet*>(*lIter);
630 // Make list of wholly imported dependents
631 RooArgSet impDeps(depAllList[i]);
632 impDeps.remove(*normDeps, true, true);
633 auto snap = new RooArgSet;
634 impDeps.snapshot(*snap);
635 impDepList.Add(snap);
636// cout << GetName() << ": list of imported dependents for term " << (*term) << " set to " << impDeps << endl ;
637
638 // Make list of cross dependents (term is self contained for these dependents,
639 // but components import dependents from other components)
640 auto crossDeps = std::unique_ptr<RooAbsCollection>{depIntNoNormList[i].selectCommon(*normDeps)};
641 snap = new RooArgSet;
642 crossDeps->snapshot(*snap);
643 crossDepList.Add(snap);
644// cout << GetName() << ": list of cross dependents for term " << (*term) << " set to " << *crossDeps << endl ;
645 }
646
647 return;
648}
649
650
651
652
653////////////////////////////////////////////////////////////////////////////////
654/// Return list of (partial) integrals of product terms for integration
655/// of p.d.f over observables iset while normalization over observables nset.
656/// Also return list of normalization sets to be used to evaluate
657/// each component in the list correctly.
658
659Int_t RooProdPdf::getPartIntList(const RooArgSet* nset, const RooArgSet* iset, const char* isetRangeName) const
660{
661 // Check if this configuration was created before
662 Int_t sterileIdx(-1);
663
664 if (static_cast<CacheElem*>(_cacheMgr.getObj(nset,iset,&sterileIdx,isetRangeName))) {
665 return _cacheMgr.lastIndex();
666 }
667
668 std::unique_ptr<CacheElem> cache = createCacheElem(nset, iset, isetRangeName);
669
670 // Store the partial integral list and return the assigned code
671 return _cacheMgr.setObj(nset, iset, cache.release(), RooNameReg::ptr(isetRangeName));
672}
673
674
675
676std::unique_ptr<RooProdPdf::CacheElem> RooProdPdf::createCacheElem(const RooArgSet* nset,
677 const RooArgSet* iset,
678 const char* isetRangeName) const
679{
680// cout << " FOLKERT::RooProdPdf::getPartIntList(" << GetName() <<") nset = " << (nset?*nset:RooArgSet()) << endl
681// << " _normRange = " << _normRange << endl
682// << " iset = " << (iset?*iset:RooArgSet()) << endl
683// << " isetRangeName = " << (isetRangeName?isetRangeName:"<null>") << endl ;
684
685 // Create containers for partial integral components to be generated
686 auto cache = std::make_unique<CacheElem>();
687
688 // Factorize the product in irreducible terms for this nset
689 RooLinkedList terms, norms, imp, ints, cross;
690// cout << "RooProdPdf::getPIL -- now calling factorizeProduct()" << endl ;
691
692
693 // Normalization set used for factorization
694 RooArgSet factNset(nset ? (*nset) : _defNormSet);
695// cout << GetName() << "factNset = " << factNset << endl ;
696
697 factorizeProduct(factNset, iset ? (*iset) : RooArgSet(), terms, norms, imp, cross, ints);
698
699 RooArgSet *norm, *integ, *xdeps, *imps;
700
701 // Group irriducible terms that need to be (partially) integrated together
702 std::list<std::vector<RooArgSet*>> groupedList;
703 RooArgSet outerIntDeps;
704// cout << "RooProdPdf::getPIL -- now calling groupProductTerms()" << endl;
705 groupProductTerms(groupedList, outerIntDeps, terms, norms, imp, ints, cross);
706
707 // Loop over groups
708// cout<<"FK: pdf("<<GetName()<<") Starting selecting F(x|y)!"<<endl;
709 // Find groups of type F(x|y), i.e. termImpSet!=0, construct ratio object
710 std::map<std::string, RooArgSet> ratioTerms;
711 for (auto const& group : groupedList) {
712 if (1 == group.size()) {
713// cout<<"FK: Starting Single Term"<<endl;
714
715 RooArgSet* term = group[0];
716
717 Int_t termIdx = terms.IndexOf(term);
718 norm=(RooArgSet*) norms.At(termIdx);
719 imps=(RooArgSet*)imp.At(termIdx);
720 RooArgSet termNSet(*norm), termImpSet(*imps);
721
722// cout<<"FK: termImpSet.getSize() = "<<termImpSet.getSize()<< " " << termImpSet << endl;
723// cout<<"FK: _refRangeName = "<<_refRangeName<<endl;
724
725 if (!termImpSet.empty() && 0 != _refRangeName) {
726
727// cout << "WVE now here" << endl;
728
729 // WVE we can skip this if the ref range is equal to the normalization range
730 bool rangeIdentical(true);
731// cout << "_normRange = " << _normRange << " _refRangeName = " << RooNameReg::str(_refRangeName) << endl ;
732 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
733 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
734 if (_normRange.Length() > 0) {
735 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
736 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
737 }
738 else{
739 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
740 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
741 }
742 }
743// cout<<"FK: rangeIdentical Single = "<<(rangeIdentical ? 'T':'F')<<endl;
744 // coverity[CONSTANT_EXPRESSION_RESULT]
745 // LM : avoid making integral ratio if range is the same. Why was not included ??? (same at line 857)
746 if (!rangeIdentical ) {
747// cout << "PREPARING RATIO HERE (SINGLE TERM)" << endl ;
748 auto ratio = makeCondPdfRatioCorr(*(RooAbsReal*)term->first(), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
749 std::ostringstream str; termImpSet.printValue(str);
750// cout << GetName() << "inserting ratio term" << endl;
751 ratioTerms[str.str()].addOwned(std::move(ratio));
752 }
753 }
754
755 } else {
756// cout<<"FK: Starting Composite Term"<<endl;
757
758 for (auto const& term : group) {
759
760 Int_t termIdx = terms.IndexOf(term);
761 norm=(RooArgSet*) norms.At(termIdx);
762 imps=(RooArgSet*)imp.At(termIdx);
763 RooArgSet termNSet(*norm), termImpSet(*imps);
764
765 if (!termImpSet.empty() && 0 != _refRangeName) {
766
767 // WVE we can skip this if the ref range is equal to the normalization range
768 bool rangeIdentical(true);
769 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
770 if(_normRange.Length() > 0) {
771 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
772 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
773 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
774 }
775 } else {
776 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
777 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
778 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
779 }
780 }
781// cout<<"FK: rangeIdentical Composite = "<<(rangeIdentical ? 'T':'F') <<endl;
782 if (!rangeIdentical ) {
783// cout << "PREPARING RATIO HERE (COMPOSITE TERM)" << endl ;
784 auto ratio = makeCondPdfRatioCorr(*(RooAbsReal*)term->first(), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
785 std::ostringstream str; termImpSet.printValue(str);
786 ratioTerms[str.str()].addOwned(std::move(ratio));
787 }
788 }
789 }
790 }
791
792 }
793
794 // Find groups with y as termNSet
795 // Replace G(y) with (G(y),ratio)
796 for (auto const& group : groupedList) {
797 for (auto const& term : group) {
798 Int_t termIdx = terms.IndexOf(term);
799 norm = (RooArgSet*) norms.At(termIdx);
800 imps = (RooArgSet*) imp.At(termIdx);
801 RooArgSet termNSet(*norm), termImpSet(*imps);
802
803 // If termNset matches index of ratioTerms, insert ratio here
804 ostringstream str; termNSet.printValue(str);
805 if (!ratioTerms[str.str()].empty()) {
806// cout << "MUST INSERT RATIO OBJECT IN TERM (COMPOSITE)" << *term << endl;
807 term->add(ratioTerms[str.str()]);
808 cache->_ownedList.addOwned(std::move(ratioTerms[str.str()]));
809 }
810 }
811 }
812
813 for (auto const& group : groupedList) {
814// cout << GetName() << ":now processing group" << endl;
815// group->Print("1");
816
817 if (1 == group.size()) {
818// cout << "processing atomic item" << endl;
819 RooArgSet* term = group[0];
820
821 Int_t termIdx = terms.IndexOf(term);
822 norm = (RooArgSet*) norms.At(termIdx);
823 integ = (RooArgSet*) ints.At(termIdx);
824 xdeps = (RooArgSet*) cross.At(termIdx);
825 imps = (RooArgSet*) imp.At(termIdx);
826
827 RooArgSet termNSet, termISet, termXSet, termImpSet;
828
829 // Take list of normalization, integrated dependents from factorization algorithm
830 termISet.add(*integ);
831 termNSet.add(*norm);
832
833 // Cross-imported integrated dependents
834 termXSet.add(*xdeps);
835 termImpSet.add(*imps);
836
837 // Add prefab term to partIntList.
838 bool isOwned(false);
839 vector<RooAbsReal*> func = processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned);
840 if (func[0]) {
841 cache->_partList.add(*func[0]);
842 if (isOwned) cache->_ownedList.addOwned(*func[0]);
843
844 cache->_normList.emplace_back(std::make_unique<RooArgSet>());
845 norm->snapshot(*cache->_normList.back(), false);
846
847 cache->_numList.addOwned(*func[1]);
848 cache->_denList.addOwned(*func[2]);
849// cout << "func[0]=" << func[0]->ClassName() << "::" << func[0]->GetName() << endl;
850// cout << "func[1]=" << func[1]->ClassName() << "::" << func[1]->GetName() << endl;
851// cout << "func[2]=" << func[2]->ClassName() << "::" << func[2]->GetName() << endl;
852 }
853 } else {
854// cout << "processing composite item" << endl;
855 RooArgSet compTermSet, compTermNorm, compTermNum, compTermDen;
856 for (auto const& term : group) {
857// cout << GetName() << ": processing term " << (*term) << " of composite item" << endl ;
858 Int_t termIdx = terms.IndexOf(term);
859 norm = (RooArgSet*) norms.At(termIdx);
860 integ = (RooArgSet*) ints.At(termIdx);
861 xdeps = (RooArgSet*) cross.At(termIdx);
862 imps = (RooArgSet*) imp.At(termIdx);
863
864 RooArgSet termNSet, termISet, termXSet, termImpSet;
865 termISet.add(*integ);
866 termNSet.add(*norm);
867 termXSet.add(*xdeps);
868 termImpSet.add(*imps);
869
870 // Remove outer integration dependents from termISet
871 termISet.remove(outerIntDeps, true, true);
872
873 bool isOwned = false;
874 vector<RooAbsReal*> func = processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned, true);
875// cout << GetName() << ": created composite term component " << func[0]->GetName() << endl;
876 if (func[0]) {
877 compTermSet.add(*func[0]);
878 if (isOwned) cache->_ownedList.addOwned(*func[0]);
879 compTermNorm.add(*norm, false);
880
881 compTermNum.add(*func[1]);
882 compTermDen.add(*func[2]);
883 //cache->_numList.add(*func[1]);
884 //cache->_denList.add(*func[2]);
885
886 }
887 }
888
889// cout << GetName() << ": constructing special composite product" << endl;
890// cout << GetName() << ": compTermSet = " ; compTermSet.Print("1");
891
892 // WVE THIS NEEDS TO BE REARRANGED
893
894 // compTermset is set van partial integrals to be multiplied
895 // prodtmp = product (compTermSet)
896 // inttmp = int ( prodtmp ) d (outerIntDeps) _range_isetRangeName
897
898 const std::string prodname = makeRGPPName("SPECPROD", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
899 auto prodtmp = std::make_unique<RooProduct>(prodname.c_str(), prodname.c_str(), compTermSet);
900
901 const std::string intname = makeRGPPName("SPECINT", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
902 auto inttmp = std::make_unique<RooRealIntegral>(intname.c_str(), intname.c_str(), *prodtmp, outerIntDeps, nullptr, nullptr, isetRangeName);
903 inttmp->setStringAttribute("PROD_TERM_TYPE", "SPECINT");
904
905 cache->_partList.add(*inttmp);
906
907 // Product of numerator terms
908 const string prodname_num = makeRGPPName("SPECPROD_NUM", compTermNum, RooArgSet(), RooArgSet(), 0);
909 auto prodtmp_num = std::make_unique<RooProduct>(prodname_num.c_str(), prodname_num.c_str(), compTermNum);
910 prodtmp_num->addOwnedComponents(compTermNum);
911
912 // Product of denominator terms
913 const string prodname_den = makeRGPPName("SPECPROD_DEN", compTermDen, RooArgSet(), RooArgSet(), 0);
914 auto prodtmp_den = std::make_unique<RooProduct>(prodname_den.c_str(), prodname_den.c_str(), compTermDen);
915 prodtmp_den->addOwnedComponents(compTermDen);
916
917 // Ratio
918 std::string name = Form("SPEC_RATIO(%s,%s)", prodname_num.c_str(), prodname_den.c_str());
919 auto ndr = std::make_unique<RooFormulaVar>(name.c_str(), "@0/@1", RooArgList(*prodtmp_num, *prodtmp_den));
920
921 // Integral of ratio
922 std::unique_ptr<RooAbsReal> numtmp{ndr->createIntegral(outerIntDeps,isetRangeName)};
923 numtmp->addOwnedComponents(std::move(ndr));
924
925 cache->_ownedList.addOwned(std::move(prodtmp));
926 cache->_ownedList.addOwned(std::move(inttmp));
927 cache->_ownedList.addOwned(std::move(prodtmp_num));
928 cache->_ownedList.addOwned(std::move(prodtmp_den));
929 cache->_numList.addOwned(std::move(numtmp));
930 cache->_denList.addOwned(*(RooAbsArg*)RooFit::RooConst(1).clone("1"));
931 cache->_normList.emplace_back(std::make_unique<RooArgSet>());
932 compTermNorm.snapshot(*cache->_normList.back(), false);
933 }
934 }
935
936 // Need to rearrange product in case of multiple ranges
937 if (_normRange.Contains(",")) {
938 rearrangeProduct(*cache);
939 }
940
941 // We own contents of all lists filled by factorizeProduct()
942 terms.Delete();
943 ints.Delete();
944 imp.Delete();
945 norms.Delete();
946 cross.Delete();
947
948 return cache;
949}
950
951
952
953////////////////////////////////////////////////////////////////////////////////
954/// For single normalization ranges
955
956std::unique_ptr<RooAbsReal> RooProdPdf::makeCondPdfRatioCorr(RooAbsReal& pdf, const RooArgSet& termNset, const RooArgSet& /*termImpSet*/, const char* normRangeTmp, const char* refRange) const
957{
958 RooAbsReal* ratio_num = pdf.createIntegral(termNset,normRangeTmp) ;
959 RooAbsReal* ratio_den = pdf.createIntegral(termNset,refRange) ;
960 auto ratio = std::make_unique<RooFormulaVar>(Form("ratio(%s,%s)",ratio_num->GetName(),ratio_den->GetName()),"@0/@1",
961 RooArgList(*ratio_num,*ratio_den)) ;
962
963 ratio->addOwnedComponents(RooArgSet(*ratio_num,*ratio_den)) ;
964 ratio->setAttribute("RATIO_TERM") ;
965 return ratio ;
966}
967
968
969
970
971////////////////////////////////////////////////////////////////////////////////
972
974{
975 RooAbsReal *part, *num, *den ;
976 RooArgSet nomList ;
977
978 list<string> rangeComps ;
979 {
980 std::vector<char> buf(strlen(_normRange.Data()) + 1);
981 strcpy(buf.data(),_normRange.Data()) ;
982 char* save(0) ;
983 char* token = R__STRTOK_R(buf.data(),",",&save) ;
984 while(token) {
985 rangeComps.push_back(token) ;
986 token = R__STRTOK_R(0,",",&save) ;
987 }
988 }
989
990
991 map<string,RooArgSet> denListList ;
992 RooArgSet specIntDeps ;
993 string specIntRange ;
994
995// cout << "THIS IS REARRANGEPRODUCT" << endl ;
996
997 for (std::size_t i = 0; i < cache._partList.size(); i++) {
998
999 part = static_cast<RooAbsReal*>(cache._partList.at(i));
1000 num = static_cast<RooAbsReal*>(cache._numList.at(i));
1001 den = static_cast<RooAbsReal*>(cache._denList.at(i));
1002 i++;
1003
1004// cout << "now processing part " << part->GetName() << " of type " << part->getStringAttribute("PROD_TERM_TYPE") << endl ;
1005// cout << "corresponding numerator = " << num->GetName() << endl ;
1006// cout << "corresponding denominator = " << den->GetName() << endl ;
1007
1008
1009 RooFormulaVar* ratio(0) ;
1010 RooArgSet origNumTerm ;
1011
1012 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1013
1014 RooRealIntegral* orig = (RooRealIntegral*) num;
1015 RooFormulaVar* specratio = (RooFormulaVar*) &orig->integrand() ;
1016 RooProduct* func = (RooProduct*) specratio->getParameter(0) ;
1017
1018 std::unique_ptr<RooArgSet> components{orig->getComponents()};
1019 for(RooAbsArg * carg : *components) {
1020 if (carg->getAttribute("RATIO_TERM")) {
1021 ratio = (RooFormulaVar*)carg ;
1022 break ;
1023 }
1024 }
1025
1026 if (ratio) {
1027 RooCustomizer cust(*func,"blah") ;
1028 cust.replaceArg(*ratio,RooFit::RooConst(1)) ;
1029 RooAbsArg* funcCust = cust.build() ;
1030// cout << "customized function = " << endl ;
1031// funcCust->printComponentTree() ;
1032 nomList.add(*funcCust) ;
1033 } else {
1034 nomList.add(*func) ;
1035 }
1036
1037
1038 } else {
1039
1040 // Find the ratio term
1041 RooAbsReal* func = num;
1042 // If top level object is integral, navigate to integrand
1043 if (func->InheritsFrom(RooRealIntegral::Class())) {
1044 func = (RooAbsReal*) &((RooRealIntegral*)(func))->integrand();
1045 }
1046 if (func->InheritsFrom(RooProduct::Class())) {
1047// cout << "product term found: " ; func->Print() ;
1048 for(RooAbsArg * arg : static_cast<RooProduct*>(func)->components()) {
1049 if (arg->getAttribute("RATIO_TERM")) {
1050 ratio = (RooFormulaVar*)(arg) ;
1051 } else {
1052 origNumTerm.add(*arg) ;
1053 }
1054 }
1055 }
1056
1057 if (ratio) {
1058// cout << "Found ratio term in numerator: " << ratio->GetName() << endl ;
1059// cout << "Adding only original term to numerator: " << origNumTerm << endl ;
1060 nomList.add(origNumTerm) ;
1061 } else {
1062 nomList.add(*num) ;
1063 }
1064
1065 }
1066
1067 for (list<string>::iterator iter = rangeComps.begin() ; iter != rangeComps.end() ; ++iter) {
1068 // If denominator is an integral, make a clone with the integration range adjusted to
1069 // the selected component of the normalization integral
1070// cout << "NOW PROCESSING DENOMINATOR " << den->ClassName() << "::" << den->GetName() << endl ;
1071
1072 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1073
1074// cout << "create integral: SPECINT case" << endl ;
1075 RooRealIntegral* orig = (RooRealIntegral*) num;
1076 RooFormulaVar* specRatio = (RooFormulaVar*) &orig->integrand() ;
1077 specIntDeps.add(orig->intVars()) ;
1078 if (orig->intRange()) {
1079 specIntRange = orig->intRange() ;
1080 }
1081 //RooProduct* numtmp = (RooProduct*) specRatio->getParameter(0) ;
1082 RooProduct* dentmp = (RooProduct*) specRatio->getParameter(1) ;
1083
1084// cout << "numtmp = " << numtmp->ClassName() << "::" << numtmp->GetName() << endl ;
1085// cout << "dentmp = " << dentmp->ClassName() << "::" << dentmp->GetName() << endl ;
1086
1087// cout << "denominator components are " << dentmp->components() << endl ;
1088 for (auto* parg : static_range_cast<RooAbsReal*>(dentmp->components())) {
1089// cout << "now processing denominator component " << parg->ClassName() << "::" << parg->GetName() << endl ;
1090
1091 if (ratio && parg->dependsOn(*ratio)) {
1092// cout << "depends in value of ratio" << endl ;
1093
1094 // Make specialize ratio instance
1095 RooAbsReal* specializedRatio = specializeRatio(*(RooFormulaVar*)ratio,iter->c_str()) ;
1096
1097// cout << "specRatio = " << endl ;
1098// specializedRatio->printComponentTree() ;
1099
1100 // Replace generic ratio with specialized ratio
1101 RooAbsArg *partCust(0) ;
1102 if (parg->InheritsFrom(RooAddition::Class())) {
1103
1104
1105
1106 RooAddition* tmpadd = (RooAddition*)(parg) ;
1107
1108 RooCustomizer cust(*tmpadd->list1().first(),Form("blah_%s",iter->c_str())) ;
1109 cust.replaceArg(*ratio,*specializedRatio) ;
1110 partCust = cust.build() ;
1111
1112 } else {
1113 RooCustomizer cust(*parg,Form("blah_%s",iter->c_str())) ;
1114 cust.replaceArg(*ratio,*specializedRatio) ;
1115 partCust = cust.build() ;
1116 }
1117
1118 // Print customized denominator
1119// cout << "customized function = " << endl ;
1120// partCust->printComponentTree() ;
1121
1122 RooAbsReal* specializedPartCust = specializeIntegral(*(RooAbsReal*)partCust,iter->c_str()) ;
1123
1124 // Finally divide again by ratio
1125 string name = Form("%s_divided_by_ratio",specializedPartCust->GetName()) ;
1126 RooFormulaVar* specIntFinal = new RooFormulaVar(name.c_str(),"@0/@1",RooArgList(*specializedPartCust,*specializedRatio)) ;
1127
1128 denListList[*iter].add(*specIntFinal) ;
1129 } else {
1130
1131// cout << "does NOT depend on value of ratio" << endl ;
1132// parg->Print("t") ;
1133
1134 denListList[*iter].add(*specializeIntegral(*parg,iter->c_str())) ;
1135
1136 }
1137 }
1138// cout << "end iteration over denominator components" << endl ;
1139 } else {
1140
1141 if (ratio) {
1142
1143 RooAbsReal* specRatio = specializeRatio(*(RooFormulaVar*)ratio,iter->c_str()) ;
1144
1145 // If integral is 'Int r(y)*g(y) dy ' then divide a posteriori by r(y)
1146// cout << "have ratio, orig den = " << den->GetName() << endl ;
1147
1148 RooArgSet tmp(origNumTerm) ;
1149 tmp.add(*specRatio) ;
1150 const string pname = makeRGPPName("PROD",tmp,RooArgSet(),RooArgSet(),0) ;
1151 RooProduct* specDenProd = new RooProduct(pname.c_str(),pname.c_str(),tmp) ;
1152 RooAbsReal* specInt(0) ;
1153
1155 specInt = specDenProd->createIntegral(((RooRealIntegral*)den)->intVars(),iter->c_str()) ;
1156 } else if (den->InheritsFrom(RooAddition::Class())) {
1157 RooAddition* orig = (RooAddition*)den ;
1158 RooRealIntegral* origInt = (RooRealIntegral*) orig->list1().first() ;
1159 specInt = specDenProd->createIntegral(origInt->intVars(),iter->c_str()) ;
1160 } else {
1161 throw string("this should not happen") ;
1162 }
1163
1164 //RooAbsReal* specInt = specializeIntegral(*den,iter->c_str()) ;
1165 string name = Form("%s_divided_by_ratio",specInt->GetName()) ;
1166 RooFormulaVar* specIntFinal = new RooFormulaVar(name.c_str(),"@0/@1",RooArgList(*specInt,*specRatio)) ;
1167 denListList[*iter].add(*specIntFinal) ;
1168 } else {
1169 denListList[*iter].add(*specializeIntegral(*den,iter->c_str())) ;
1170 }
1171
1172 }
1173 }
1174
1175 }
1176
1177 // Do not rearrage terms if numerator and denominator are effectively empty
1178 if (nomList.empty()) {
1179 return ;
1180 }
1181
1182 string name = Form("%s_numerator",GetName()) ;
1183 // WVE FIX THIS (2)
1184
1185 RooAbsReal* numerator = new RooProduct(name.c_str(),name.c_str(),nomList) ;
1186
1187 RooArgSet products ;
1188// cout << "nomList = " << nomList << endl ;
1189 for (map<string,RooArgSet>::iterator iter = denListList.begin() ; iter != denListList.end() ; ++iter) {
1190// cout << "denList[" << iter->first << "] = " << iter->second << endl ;
1191 name = Form("%s_denominator_comp_%s",GetName(),iter->first.c_str()) ;
1192 // WVE FIX THIS (2)
1193 RooProduct* prod_comp = new RooProduct(name.c_str(),name.c_str(),iter->second) ;
1194 products.add(*prod_comp) ;
1195 }
1196 name = Form("%s_denominator_sum",GetName()) ;
1197 RooAbsReal* norm = new RooAddition(name.c_str(),name.c_str(),products) ;
1198 norm->addOwnedComponents(products) ;
1199
1200 if (specIntDeps.getSize()>0) {
1201 // Apply posterior integration required for SPECINT case
1202
1203 string namesr = Form("SPEC_RATIO(%s,%s)",numerator->GetName(),norm->GetName()) ;
1204 RooFormulaVar* ndr = new RooFormulaVar(namesr.c_str(),"@0/@1",RooArgList(*numerator,*norm)) ;
1205
1206 // Integral of ratio
1207 RooAbsReal* numtmp = ndr->createIntegral(specIntDeps,specIntRange.c_str()) ;
1208
1209 numerator = numtmp ;
1210 norm = (RooAbsReal*) RooFit::RooConst(1).Clone() ;
1211 }
1212
1213
1214// cout << "numerator" << endl ;
1215// numerator->printComponentTree("",0,5) ;
1216// cout << "denominator" << endl ;
1217// norm->printComponentTree("",0,5) ;
1218
1219
1220 // WVE DEBUG
1221 //RooMsgService::instance().debugWorkspace()->import(RooArgSet(*numerator,*norm)) ;
1222
1223 cache._rearrangedNum.reset(numerator);
1224 cache._rearrangedDen.reset(norm);
1225 cache._isRearranged = true ;
1226
1227}
1228
1229
1230////////////////////////////////////////////////////////////////////////////////
1231
1232RooAbsReal* RooProdPdf::specializeRatio(RooFormulaVar& input, const char* targetRangeName) const
1233{
1234 RooRealIntegral* numint = (RooRealIntegral*) input.getParameter(0) ;
1235 RooRealIntegral* denint = (RooRealIntegral*) input.getParameter(1) ;
1236
1237 RooAbsReal* numint_spec = specializeIntegral(*numint,targetRangeName) ;
1238
1239 RooAbsReal* ret = new RooFormulaVar(Form("ratio(%s,%s)",numint_spec->GetName(),denint->GetName()),"@0/@1",RooArgList(*numint_spec,*denint)) ;
1240 ret->addOwnedComponents(*numint_spec) ;
1241
1242 return ret ;
1243}
1244
1245
1246
1247////////////////////////////////////////////////////////////////////////////////
1248
1249RooAbsReal* RooProdPdf::specializeIntegral(RooAbsReal& input, const char* targetRangeName) const
1250{
1251 if (input.InheritsFrom(RooRealIntegral::Class())) {
1252
1253 // If input is integral, recreate integral but override integration range to be targetRangeName
1255// cout << "creating integral: integrand = " << orig->integrand().GetName() << " vars = " << orig->intVars() << " range = " << targetRangeName << endl ;
1256 return orig->integrand().createIntegral(orig->intVars(),targetRangeName) ;
1257
1258 } else if (input.InheritsFrom(RooAddition::Class())) {
1259
1260 // If input is sum of integrals, recreate integral from first component of set, but override integration range to be targetRangeName
1261 RooAddition* orig = (RooAddition*)&input ;
1262 RooRealIntegral* origInt = (RooRealIntegral*) orig->list1().first() ;
1263// cout << "creating integral from addition: integrand = " << origInt->integrand().GetName() << " vars = " << origInt->intVars() << " range = " << targetRangeName << endl ;
1264 return origInt->integrand().createIntegral(origInt->intVars(),targetRangeName) ;
1265
1266 } else {
1267
1268// cout << "specializeIntegral: unknown input type " << input.ClassName() << "::" << input.GetName() << endl ;
1269 }
1270
1271 return &input ;
1272}
1273
1274
1275////////////////////////////////////////////////////////////////////////////////
1276/// Group product into terms that can be calculated independently
1277
1278void RooProdPdf::groupProductTerms(std::list<std::vector<RooArgSet*>>& groupedTerms, RooArgSet& outerIntDeps,
1279 const RooLinkedList& terms, const RooLinkedList& norms,
1280 const RooLinkedList& imps, const RooLinkedList& ints, const RooLinkedList& /*cross*/) const
1281{
1282 // Start out with each term in its own group
1283 for(auto * term : static_range_cast<RooArgSet*>(terms)) {
1284 groupedTerms.emplace_back();
1285 groupedTerms.back().emplace_back(term) ;
1286 }
1287
1288 // Make list of imported dependents that occur in any term
1289 RooArgSet allImpDeps ;
1290 for(auto * impDeps : static_range_cast<RooArgSet*>(imps)) {
1291 allImpDeps.add(*impDeps,false) ;
1292 }
1293
1294 // Make list of integrated dependents that occur in any term
1295 RooArgSet allIntDeps ;
1296 for(auto * intDeps : static_range_cast<RooArgSet*>(ints)) {
1297 allIntDeps.add(*intDeps,false) ;
1298 }
1299
1300 outerIntDeps.removeAll() ;
1301 outerIntDeps.add(*std::unique_ptr<RooArgSet>{static_cast<RooArgSet*>(allIntDeps.selectCommon(allImpDeps))});
1302
1303 // Now iteratively merge groups that should be (partially) integrated together
1304 for(RooAbsArg * outerIntDep : outerIntDeps) {
1305
1306 // Collect groups that feature this dependent
1307 std::vector<RooArgSet*>* newGroup = nullptr ;
1308
1309 // Loop over groups
1310 bool needMerge = false ;
1311 auto group = groupedTerms.begin();
1312 auto nGroups = groupedTerms.size();
1313 for (size_t iGroup = 0; iGroup < nGroups; ++iGroup) {
1314
1315 // See if any term in this group depends in any ay on outerDepInt
1316 for (auto const& term2 : *group) {
1317
1318 Int_t termIdx = terms.IndexOf(term2) ;
1319 RooArgSet* termNormDeps = (RooArgSet*) norms.At(termIdx) ;
1320 RooArgSet* termIntDeps = (RooArgSet*) ints.At(termIdx) ;
1321 RooArgSet* termImpDeps = (RooArgSet*) imps.At(termIdx) ;
1322
1323 if (termNormDeps->contains(*outerIntDep) ||
1324 termIntDeps->contains(*outerIntDep) ||
1325 termImpDeps->contains(*outerIntDep)) {
1326 needMerge = true ;
1327 }
1328
1329 }
1330
1331 if (needMerge) {
1332 // Create composite group if not yet existing
1333 if (newGroup==nullptr) {
1334 groupedTerms.emplace_back() ;
1335 newGroup = &groupedTerms.back() ;
1336 }
1337
1338 // Add terms of this group to new term
1339 for (auto& term2 : *group) {
1340 newGroup->emplace_back(term2) ;
1341 }
1342
1343 // Remove this non-owning group from list
1344 group = groupedTerms.erase(group);
1345 } else {
1346 ++group;
1347 }
1348 }
1349
1350 }
1351}
1352
1353
1354
1355////////////////////////////////////////////////////////////////////////////////
1356/// Calculate integrals of factorized product terms over observables iset while normalized
1357/// to observables in nset.
1358
1359std::vector<RooAbsReal*> RooProdPdf::processProductTerm(const RooArgSet* nset, const RooArgSet* iset, const char* isetRangeName,
1360 const RooArgSet* term,const RooArgSet& termNSet, const RooArgSet& termISet,
1361 bool& isOwned, bool forceWrap) const
1362{
1363 vector<RooAbsReal*> ret(3) ; ret[0] = 0 ; ret[1] = 0 ; ret[2] = 0 ;
1364
1365 // CASE I: factorizing term: term is integrated over all normalizing observables
1366 // -----------------------------------------------------------------------------
1367 // Check if all observbales of this term are integrated. If so the term cancels
1368 if (termNSet.getSize()>0 && termNSet.getSize()==termISet.getSize() && isetRangeName==0) {
1369
1370
1371 //cout << "processProductTerm(" << GetName() << ") case I " << endl ;
1372
1373 // Term factorizes
1374 return ret ;
1375 }
1376
1377 // CASE II: Dropped terms: if term is entirely unnormalized, it should be dropped
1378 // ------------------------------------------------------------------------------
1379 if (nset && termNSet.empty()) {
1380
1381 //cout << "processProductTerm(" << GetName() << ") case II " << endl ;
1382
1383 // Drop terms that are not asked to be normalized
1384 return ret ;
1385 }
1386
1387 if (iset && termISet.getSize()>0) {
1388 if (term->getSize()==1) {
1389
1390 // CASE IIIa: Normalized and partially integrated single PDF term
1391 //---------------------------------------------------------------
1392
1393 RooAbsPdf* pdf = (RooAbsPdf*) term->first() ;
1394
1395 RooAbsReal* partInt = pdf->createIntegral(termISet,termNSet,isetRangeName) ;
1396 partInt->setOperMode(operMode()) ;
1397 partInt->setStringAttribute("PROD_TERM_TYPE","IIIa") ;
1398
1399 isOwned=true ;
1400
1401 //cout << "processProductTerm(" << GetName() << ") case IIIa func = " << partInt->GetName() << endl ;
1402
1403 ret[0] = partInt ;
1404
1405 // Split mode results
1406 ret[1] = pdf->createIntegral(termISet,isetRangeName) ;
1407 ret[2] = pdf->createIntegral(termNSet,normRange()) ;
1408
1409 return ret ;
1410
1411
1412 } else {
1413
1414 // CASE IIIb: Normalized and partially integrated composite PDF term
1415 //---------------------------------------------------------------
1416
1417 // Use auxiliary class RooGenProdProj to calculate this term
1418 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1419 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName) ;
1420 partInt->setStringAttribute("PROD_TERM_TYPE","IIIb") ;
1421 partInt->setOperMode(operMode()) ;
1422
1423 //cout << "processProductTerm(" << GetName() << ") case IIIb func = " << partInt->GetName() << endl ;
1424
1425 isOwned=true ;
1426 ret[0] = partInt ;
1427
1428 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),0) ;
1429
1430 // WVE FIX THIS
1431 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1432
1433 ret[1] = tmp_prod->createIntegral(termISet,isetRangeName) ;
1434 ret[2] = tmp_prod->createIntegral(termNSet,normRange()) ;
1435
1436 return ret ;
1437 }
1438 }
1439
1440 // CASE IVa: Normalized non-integrated composite PDF term
1441 // -------------------------------------------------------
1442 if (nset && nset->getSize()>0 && term->getSize()>1) {
1443 // Composite term needs normalized integration
1444
1445 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1446 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName,normRange()) ;
1448
1449 partInt->setStringAttribute("PROD_TERM_TYPE","IVa") ;
1450 partInt->setOperMode(operMode()) ;
1451
1452 //cout << "processProductTerm(" << GetName() << ") case IVa func = " << partInt->GetName() << endl ;
1453
1454 isOwned=true ;
1455 ret[0] = partInt ;
1456
1457 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),0) ;
1458
1459 // WVE FIX THIS
1460 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1461
1462 ret[1] = tmp_prod->createIntegral(termISet,isetRangeName) ;
1463 ret[2] = tmp_prod->createIntegral(termNSet,normRange()) ;
1464
1465 return ret ;
1466 }
1467
1468 // CASE IVb: Normalized, non-integrated single PDF term
1469 // -----------------------------------------------------
1470 for (auto* pdf : static_range_cast<RooAbsPdf*>(*term)) {
1471
1472 if (forceWrap) {
1473
1474 // Construct representative name of normalization wrapper
1475 TString name(pdf->GetName()) ;
1476 name.Append("_NORM[") ;
1477 bool first(true) ;
1478 for (auto const* arg : termNSet) {
1479 if (!first) {
1480 name.Append(",") ;
1481 } else {
1482 first=false ;
1483 }
1484 name.Append(arg->GetName()) ;
1485 }
1486 if (normRange()) {
1487 name.Append("|") ;
1488 name.Append(normRange()) ;
1489 }
1490 name.Append("]") ;
1491
1492 RooAbsReal* partInt = new RooRealIntegral(name.Data(),name.Data(),*pdf,RooArgSet(),&termNSet) ;
1493 partInt->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1494 isOwned=true ;
1495
1496 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << partInt->GetName() << endl ;
1497
1498 ret[0] = partInt ;
1499
1500 ret[1] = pdf->createIntegral(RooArgSet()) ;
1501 ret[2] = pdf->createIntegral(termNSet,normRange()) ;
1502
1503 return ret ;
1504
1505
1506 } else {
1507 isOwned=false ;
1508
1509 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << pdf->GetName() << endl ;
1510
1511
1512 pdf->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1513 ret[0] = pdf ;
1514
1515 ret[1] = pdf->createIntegral(RooArgSet()) ;
1516 ret[2] = termNSet.getSize()>0 ? pdf->createIntegral(termNSet,normRange()) : ((RooAbsReal*)RooFit::RooConst(1).clone("1")) ;
1517 return ret ;
1518 }
1519 }
1520
1521 coutE(Eval) << "RooProdPdf::processProductTerm(" << GetName() << ") unidentified term!!!" << endl ;
1522 return ret ;
1523}
1524
1525
1526
1527
1528////////////////////////////////////////////////////////////////////////////////
1529/// Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1530
1531std::string RooProdPdf::makeRGPPName(const char* pfx, const RooArgSet& term, const RooArgSet& iset,
1532 const RooArgSet& nset, const char* isetRangeName) const
1533{
1534 // Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1535
1536 std::ostringstream os(pfx);
1537 os << "[";
1538
1539 // Encode component names
1540 bool first(true) ;
1541 for (auto const* pdf : static_range_cast<RooAbsPdf*>(term)) {
1542 if (!first) os << "_X_";
1543 first = false;
1544 os << pdf->GetName();
1545 }
1546 os << "]" << integralNameSuffix(iset,&nset,isetRangeName,true);
1547
1548 return os.str();
1549}
1550
1551
1552
1553////////////////////////////////////////////////////////////////////////////////
1554/// Force RooRealIntegral to offer all observables for internal integration
1555
1557{
1558 return true ;
1559}
1560
1561
1562
1563////////////////////////////////////////////////////////////////////////////////
1564/// Determine which part (if any) of given integral can be performed analytically.
1565/// If any analytical integration is possible, return integration scenario code.
1566///
1567/// RooProdPdf implements two strategies in implementing analytical integrals
1568///
1569/// First, PDF components whose entire set of dependents are requested to be integrated
1570/// can be dropped from the product, as they will integrate out to 1 by construction
1571///
1572/// Second, RooProdPdf queries each remaining component PDF for its analytical integration
1573/// capability of the requested set ('allVars'). It finds the largest common set of variables
1574/// that can be integrated by all remaining components. If such a set exists, it reconfirms that
1575/// each component is capable of analytically integrating the common set, and combines the components
1576/// individual integration codes into a single integration code valid for RooProdPdf.
1577
1579 const RooArgSet* normSet, const char* rangeName) const
1580{
1581 if (_forceNumInt) return 0 ;
1582
1583 // Declare that we can analytically integrate all requested observables
1584 analVars.add(allVars) ;
1585
1586 // Retrieve (or create) the required partial integral list
1587 Int_t code = getPartIntList(normSet,&allVars,rangeName);
1588
1589 return code+1 ;
1590}
1591
1592
1593
1594
1595////////////////////////////////////////////////////////////////////////////////
1596/// Return analytical integral defined by given scenario code
1597
1598double RooProdPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
1599{
1600 // No integration scenario
1601 if (code==0) {
1602 return getVal(normSet) ;
1603 }
1604
1605
1606 // WVE needs adaptation for rangename feature
1607
1608 // Partial integration scenarios
1609 CacheElem* cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
1610
1611 // If cache has been sterilized, revive this slot
1612 if (cache==0) {
1613 std::unique_ptr<RooArgSet> vars{getParameters(RooArgSet())} ;
1614 RooArgSet nset = _cacheMgr.selectFromSet1(*vars, code-1) ;
1615 RooArgSet iset = _cacheMgr.selectFromSet2(*vars, code-1) ;
1616
1617 Int_t code2 = getPartIntList(&nset, &iset, rangeName) ;
1618
1619 // preceding call to getPartIntList guarantees non-null return
1620 // coverity[NULL_RETURNS]
1621 cache = (CacheElem*) _cacheMgr.getObj(&nset,&iset,&code2,rangeName) ;
1622 }
1623
1624 double val = calculate(*cache,true) ;
1625// cout << "RPP::aIWN(" << GetName() << ") ,code = " << code << ", value = " << val << endl ;
1626
1627 return val ;
1628}
1629
1630
1631
1632////////////////////////////////////////////////////////////////////////////////
1633/// If this product contains exactly one extendable p.d.f return the extension abilities of
1634/// that p.d.f, otherwise return CanNotBeExtended
1635
1637{
1639}
1640
1641
1642
1643////////////////////////////////////////////////////////////////////////////////
1644/// Return the expected number of events associated with the extendable input PDF
1645/// in the product. If there is no extendable term, abort.
1646
1647double RooProdPdf::expectedEvents(const RooArgSet* nset) const
1648{
1649 if (_extendedIndex<0) {
1650 coutF(Generation) << "Requesting expected number of events from a RooProdPdf that does not contain an extended p.d.f" << endl ;
1651 throw std::logic_error(std::string("RooProdPdf ") + GetName() + " could not be extended.");
1652 }
1653
1654 return ((RooAbsPdf*)_pdfList.at(_extendedIndex))->expectedEvents(nset) ;
1655}
1656
1657
1658
1659////////////////////////////////////////////////////////////////////////////////
1660/// Return generator context optimized for generating events from product p.d.f.s
1661
1663 const RooArgSet* auxProto, bool verbose) const
1664{
1665 if (_useDefaultGen) return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ;
1666 return new RooProdGenContext(*this,vars,prototype,auxProto,verbose) ;
1667}
1668
1669
1670
1671////////////////////////////////////////////////////////////////////////////////
1672/// Query internal generation capabilities of component p.d.f.s and aggregate capabilities
1673/// into master configuration passed to the generator context
1674
1675Int_t RooProdPdf::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
1676{
1677 if (!_useDefaultGen) return 0 ;
1678
1679 // Find the subset directVars that only depend on a single PDF in the product
1680 RooArgSet directSafe ;
1681 for (auto const* arg : directVars) {
1682 if (isDirectGenSafe(*arg)) directSafe.add(*arg) ;
1683 }
1684
1685
1686 // Now find direct integrator for relevant components ;
1687 std::vector<Int_t> code;
1688 code.reserve(64);
1689 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1690 RooArgSet pdfDirect ;
1691 Int_t pdfCode = pdf->getGenerator(directSafe,pdfDirect,staticInitOK);
1692 code.push_back(pdfCode);
1693 if (pdfCode != 0) {
1694 generateVars.add(pdfDirect) ;
1695 }
1696 }
1697
1698
1699 if (generateVars.getSize()>0) {
1700 Int_t masterCode = _genCode.store(code) ;
1701 return masterCode+1 ;
1702 } else {
1703 return 0 ;
1704 }
1705}
1706
1707
1708
1709////////////////////////////////////////////////////////////////////////////////
1710/// Forward one-time initialization call to component generation initialization
1711/// methods.
1712
1714{
1715 if (!_useDefaultGen) return ;
1716
1717 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1718 Int_t i(0) ;
1719 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1720 if (codeList[i]!=0) {
1721 pdf->initGenerator(codeList[i]) ;
1722 }
1723 i++ ;
1724 }
1725}
1726
1727
1728
1729////////////////////////////////////////////////////////////////////////////////
1730/// Generate a single event with configuration specified by 'code'
1731/// Defer internal generation to components as encoded in the _genCode
1732/// registry for given generator code.
1733
1735{
1736 if (!_useDefaultGen) return ;
1737
1738 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1739 Int_t i(0) ;
1740 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1741 if (codeList[i]!=0) {
1742 pdf->generateEvent(codeList[i]) ;
1743 }
1744 i++ ;
1745 }
1746}
1747
1748
1749
1750////////////////////////////////////////////////////////////////////////////////
1751/// Return RooAbsArg components contained in the cache
1752
1754{
1755 RooArgList ret ;
1756 ret.add(_partList) ;
1757 ret.add(_numList) ;
1758 ret.add(_denList) ;
1759 if (_rearrangedNum) ret.add(*_rearrangedNum) ;
1760 if (_rearrangedDen) ret.add(*_rearrangedDen) ;
1761 return ret ;
1762
1763}
1764
1765
1766
1767////////////////////////////////////////////////////////////////////////////////
1768/// Hook function to print cache contents in tree printing of RooProdPdf
1769
1770void RooProdPdf::CacheElem::printCompactTreeHook(ostream& os, const char* indent, Int_t curElem, Int_t maxElem)
1771{
1772 if (curElem==0) {
1773 os << indent << "RooProdPdf begin partial integral cache" << endl ;
1774 }
1775
1776 auto indent2 = std::string(indent) + "[" + std::to_string(curElem) + "]";
1777 for(auto const& arg : _partList) {
1778 arg->printCompactTree(os,indent2.c_str()) ;
1779 }
1780
1781 if (curElem==maxElem) {
1782 os << indent << "RooProdPdf end partial integral cache" << endl ;
1783 }
1784}
1785
1786
1787
1788////////////////////////////////////////////////////////////////////////////////
1789/// Forward determination of safety of internal generator code to
1790/// component p.d.f that would generate the given observable
1791
1793{
1794 // Only override base class behaviour if default generator method is enabled
1795 if (!_useDefaultGen) return RooAbsPdf::isDirectGenSafe(arg) ;
1796
1797 // Argument may appear in only one PDF component
1798 RooAbsPdf* thePdf(0) ;
1799 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1800
1801 if (pdf->dependsOn(arg)) {
1802 // Found PDF depending on arg
1803
1804 // If multiple PDFs depend on arg directGen is not safe
1805 if (thePdf) return false ;
1806
1807 thePdf = pdf ;
1808 }
1809 }
1810 // Forward call to relevant component PDF
1811 return thePdf?(thePdf->isDirectGenSafe(arg)):false ;
1812}
1813
1814
1815
1816////////////////////////////////////////////////////////////////////////////////
1817/// Look up user specified normalization set for given input PDF component
1818
1820{
1821 Int_t idx = _pdfList.index(&pdf) ;
1822 if (idx<0) return nullptr;
1823 return _pdfNSetList[idx].get() ;
1824}
1825
1826
1827
1828/// Add some full PDFs to the factors of this RooProdPdf.
1830{
1831 size_t numExtended = (_extendedIndex==-1) ? 0 : 1;
1832
1833 for(auto arg : pdfs) {
1834 RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg);
1835 if (!pdf) {
1836 coutW(InputArguments) << "RooProdPdf::addPdfs(" << GetName() << ") list arg "
1837 << arg->GetName() << " is not a PDF, ignored" << endl ;
1838 continue;
1839 }
1840 if(pdf->canBeExtended()) {
1841 if (_extendedIndex == -1) {
1843 } else {
1844 numExtended++;
1845 }
1846 }
1847 _pdfList.add(*pdf);
1848 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>("nset"));
1849 }
1850
1851 // Protect against multiple extended terms
1852 if (numExtended>1) {
1853 coutW(InputArguments) << "RooProdPdf::addPdfs(" << GetName()
1854 << ") WARNING: multiple components with extended terms detected,"
1855 << " product will not be extendible." << endl ;
1856 _extendedIndex = -1 ;
1857 }
1858
1859 // Reset cache
1860 _cacheMgr.reset() ;
1861
1862}
1863
1864/// Remove some PDFs from the factors of this RooProdPdf.
1866{
1867 // Remember what the extended PDF is
1868 RooAbsArg const* extPdf = _extendedIndex >= 0 ? &_pdfList[_extendedIndex] : nullptr;
1869
1870 // Actually remove the PDFs and associated nsets
1871 for(size_t i=0;i < _pdfList.size(); i++) {
1872 if(pdfs.contains(_pdfList[i])) {
1874 _pdfNSetList.erase(_pdfNSetList.begin()+i);
1875 i--;
1876 }
1877 }
1878
1879 // Since we may have removed PDFs from the list, the index of the extended
1880 // PDF in the list needs to be updated. The new index might also be -1 if the
1881 // extended PDF got removed.
1882 if(extPdf) {
1883 _extendedIndex = _pdfList.index(*extPdf);
1884 }
1885
1886 // Reset cache
1887 _cacheMgr.reset() ;
1888}
1889
1890
1891namespace {
1892
1893std::vector<TNamed const*> sortedNamePtrs(RooAbsCollection const& col)
1894{
1895 std::vector<TNamed const*> ptrs;
1896 ptrs.reserve(col.size());
1897 for(RooAbsArg* arg : col) {
1898 ptrs.push_back(arg->namePtr());
1899 }
1900 std::sort(ptrs.begin(), ptrs.end());
1901 return ptrs;
1902}
1903
1904bool sortedNamePtrsOverlap(std::vector<TNamed const*> const& ptrsA, std::vector<TNamed const*> const& ptrsB)
1905{
1906 auto pA = ptrsA.begin();
1907 auto pB = ptrsB.begin();
1908 while (pA != ptrsA.end() && pB != ptrsB.end()) {
1909 if (*pA < *pB) {
1910 ++pA;
1911 } else if (*pB < *pA) {
1912 ++pB;
1913 } else {
1914 return true;
1915 }
1916 }
1917 return false;
1918}
1919
1920} // namespace
1921
1922
1923////////////////////////////////////////////////////////////////////////////////
1924/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
1925/// The observables set is required to distinguish unambiguously p.d.f in terms
1926/// of observables and parameters, which are not constraints, and p.d.fs in terms
1927/// of parameters only, which can serve as constraints p.d.f.s
1928
1929RooArgSet* RooProdPdf::getConstraints(const RooArgSet& observables, RooArgSet& constrainedParams,
1930 bool stripDisconnected, bool removeConstraintsFromPdf) const
1931{
1932 RooArgSet constraints ;
1933 RooArgSet pdfParams, conParams ;
1934
1935 // For the optimized implementation of checking if two collections overlap by name.
1936 auto observablesNamePtrs = sortedNamePtrs(observables);
1937 auto constrainedParamsNamePtrs = sortedNamePtrs(constrainedParams);
1938
1939 // Loop over PDF components
1940 for(auto * pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1941 // A constraint term is a p.d.f that does not depend on any of the listed observables
1942 // but does depends on any of the parameters that should be constrained
1943 RooArgSet tmp;
1944 pdf->getParameters(nullptr, tmp);
1945 auto tmpNamePtrs = sortedNamePtrs(tmp);
1946 // Before, there were calls to `pdf->dependsOn()` here, but they were very
1947 // expensive for large computation graphs! Given that we have to traverse
1948 // the computation graph with a call to `pdf->getParameters()` anyway, we
1949 // can just check if the set of all variables operlaps with the observables
1950 // or constraind parameters.
1951 //
1952 // We are using an optimized implementation of overlap checking. Because
1953 // the overlap is checked by name, we can check overlap of the
1954 // corresponding name pointers. The optimization can't be in
1955 // RooAbsCollection itself, because it is crucial that the memory for the
1956 // non-tmp name pointers is not reallocated for each pdf.
1957 if (!sortedNamePtrsOverlap(tmpNamePtrs, observablesNamePtrs) &&
1958 sortedNamePtrsOverlap(tmpNamePtrs, constrainedParamsNamePtrs)) {
1959 constraints.add(*pdf) ;
1960 conParams.add(tmp,true) ;
1961 } else {
1962 // We only want to add parameter, not observables. Since a call like
1963 // `pdf->getParameters(&observables)` would be expensive, we take the set
1964 // of all variables and remove the ovservables, which is much cheaper. In
1965 // a call to `pdf->getParameters(&observables)`, the observables are
1966 // matched by name, so we have to pass the `matchByNameOnly` here.
1967 tmp.remove(observables, /*silent=*/false, /*matchByNameOnly=*/true);
1968 pdfParams.add(tmp,true) ;
1969 }
1970 }
1971
1972 // Remove the constraints now from the PDF if the caller requested it
1973 if(removeConstraintsFromPdf) {
1974 const_cast<RooProdPdf*>(this)->removePdfs(constraints);
1975 }
1976
1977 // Strip any constraints that are completely decoupled from the other product terms
1978 RooArgSet* finalConstraints = new RooArgSet("constraints") ;
1979 for(auto * pdf : static_range_cast<RooAbsPdf*>(constraints)) {
1980 if (pdf->dependsOnValue(pdfParams) || !stripDisconnected) {
1981 finalConstraints->add(*pdf) ;
1982 } else {
1983 coutI(Minimization) << "RooProdPdf::getConstraints(" << GetName() << ") omitting term " << pdf->GetName()
1984 << " as constraint term as it does not share any parameters with the other pdfs in product. "
1985 << "To force inclusion in likelihood, add an explicit Constrain() argument for the target parameter" << endl ;
1986 }
1987 }
1988
1989 // Now remove from constrainedParams all parameters that occur exclusively in constraint term and not in regular pdf term
1990
1991 RooArgSet cexl;
1992 conParams.selectCommon(constrainedParams, cexl);
1993 cexl.remove(pdfParams,true,true) ;
1994 constrainedParams.remove(cexl,true,true) ;
1995
1996 return finalConstraints ;
1997}
1998
1999
2000
2001
2002////////////////////////////////////////////////////////////////////////////////
2003/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
2004/// The observables set is required to distinguish unambiguously p.d.f in terms
2005/// of observables and parameters, which are not constraints, and p.d.fs in terms
2006/// of parameters only, which can serve as constraints p.d.f.s
2007
2009{
2010 RooArgSet* connectedPars = new RooArgSet("connectedPars") ;
2011 for (const auto arg : _pdfList) {
2012 // Check if term is relevant
2013 if (arg->dependsOn(observables)) {
2014 RooArgSet tmp;
2015 arg->getParameters(&observables, tmp);
2016 connectedPars->add(tmp) ;
2017 }
2018 }
2019 return connectedPars ;
2020}
2021
2022
2023
2024
2025////////////////////////////////////////////////////////////////////////////////
2026
2027void RooProdPdf::getParametersHook(const RooArgSet* nset, RooArgSet* params, bool stripDisconnected) const
2028{
2029 if (!stripDisconnected) return ;
2030 if (!nset || nset->empty()) return ;
2031
2032 // Get/create appropriate term list for this normalization set
2033 Int_t code = getPartIntList(nset, nullptr);
2034 RooArgList & plist = static_cast<CacheElem*>(_cacheMgr.getObj(nset, &code))->_partList;
2035
2036 // Strip any terms from params that do not depend on any term
2037 RooArgSet tostrip ;
2038 for (auto param : *params) {
2039 bool anyDep(false) ;
2040 for (auto term : plist) {
2041 if (term->dependsOnValue(*param)) {
2042 anyDep=true ;
2043 }
2044 }
2045 if (!anyDep) {
2046 tostrip.add(*param) ;
2047 }
2048 }
2049
2050 if (!tostrip.empty()) {
2051 params->remove(tostrip,true,true);
2052 }
2053
2054}
2055
2056
2057
2058////////////////////////////////////////////////////////////////////////////////
2059/// Interface function used by test statistics to freeze choice of range
2060/// for interpretation of conditional product terms
2061
2062void RooProdPdf::selectNormalizationRange(const char* rangeName, bool force)
2063{
2064 if (!force && _refRangeName) {
2065 return ;
2066 }
2067
2068 fixRefRange(rangeName) ;
2069}
2070
2071
2072
2073
2074////////////////////////////////////////////////////////////////////////////////
2075
2076void RooProdPdf::fixRefRange(const char* rangeName)
2077{
2078 _refRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
2079}
2080
2081
2082
2083////////////////////////////////////////////////////////////////////////////////
2084/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2085
2086std::list<double>* RooProdPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
2087{
2088 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2089 if (std::list<double>* hint = pdf->plotSamplingHint(obs,xlo,xhi)) {
2090 return hint ;
2091 }
2092 }
2093
2094 return nullptr;
2095}
2096
2097
2098
2099////////////////////////////////////////////////////////////////////////////////
2100/// If all components that depend on obs are binned that so is the product
2101
2103{
2104 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2105 if (pdf->dependsOn(obs) && !pdf->isBinnedDistribution(obs)) {
2106 return false ;
2107 }
2108 }
2109
2110 return true ;
2111}
2112
2113
2114
2115
2116
2117
2118////////////////////////////////////////////////////////////////////////////////
2119/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2120
2121std::list<double>* RooProdPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
2122{
2123 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2124 if (std::list<double>* hint = pdf->binBoundaries(obs,xlo,xhi)) {
2125 return hint ;
2126 }
2127 }
2128
2129 return nullptr;
2130}
2131
2132
2133////////////////////////////////////////////////////////////////////////////////
2134/// Label OK'ed components of a RooProdPdf with cache-and-track, _and_ label all RooProdPdf
2135/// descendants with extra information about (conditional) normalization, needed to be able
2136/// to Cache-And-Track them outside the RooProdPdf context.
2137
2139{
2140 for (const auto parg : _pdfList) {
2141
2142 if (parg->canNodeBeCached()==Always) {
2143 trackNodes.add(*parg) ;
2144// cout << "tracking node RooProdPdf component " << parg << " " << parg->ClassName() << "::" << parg->GetName() << endl ;
2145
2146 // Additional processing to fix normalization sets in case product defines conditional observables
2147 if (RooArgSet* pdf_nset = findPdfNSet(static_cast<RooAbsPdf&>(*parg))) {
2148 // Check if conditional normalization is specified
2150 if (string("nset")==pdf_nset->GetName() && !pdf_nset->empty()) {
2151 parg->setStringAttribute("CATNormSet",getColonSeparatedNameString(*pdf_nset).c_str()) ;
2152 }
2153 if (string("cset")==pdf_nset->GetName()) {
2154 parg->setStringAttribute("CATCondSet",getColonSeparatedNameString(*pdf_nset).c_str()) ;
2155 }
2156 } else {
2157 coutW(Optimization) << "RooProdPdf::setCacheAndTrackHints(" << GetName() << ") WARNING product pdf does not specify a normalization set for component " << parg->GetName() << endl ;
2158 }
2159 }
2160 }
2161}
2162
2163
2164
2165////////////////////////////////////////////////////////////////////////////////
2166/// Customized printing of arguments of a RooProdPdf to more intuitively reflect the contents of the
2167/// product operator construction
2168
2169void RooProdPdf::printMetaArgs(ostream& os) const
2170{
2171 for (int i=0 ; i<_pdfList.getSize() ; i++) {
2172 if (i>0) os << " * " ;
2173 RooArgSet* ncset = _pdfNSetList[i].get() ;
2174 os << _pdfList.at(i)->GetName() ;
2175 if (ncset->getSize()>0) {
2176 if (string("nset")==ncset->GetName()) {
2177 os << *ncset ;
2178 } else {
2179 os << "|" ;
2180 bool first(true) ;
2181 for (auto const* arg : *ncset) {
2182 if (!first) {
2183 os << "," ;
2184 } else {
2185 first = false ;
2186 }
2187 os << arg->GetName() ;
2188 }
2189 }
2190 }
2191 }
2192 os << " " ;
2193}
2194
2195
2196
2197////////////////////////////////////////////////////////////////////////////////
2198/// Implement support for node removal
2199
2200bool RooProdPdf::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
2201{
2202 if (nameChange && _pdfList.find("REMOVAL_DUMMY")) {
2203
2204 cxcoutD(LinkStateMgmt) << "RooProdPdf::redirectServersHook(" << GetName() << "): removing REMOVAL_DUMMY" << endl ;
2205
2206 // Remove node from _pdfList proxy and remove corresponding entry from normset list
2207 RooAbsArg* pdfDel = _pdfList.find("REMOVAL_DUMMY") ;
2208
2209 _pdfNSetList.erase(_pdfNSetList.begin() + _pdfList.index("REMOVAL_DUMMY")) ;
2210 _pdfList.remove(*pdfDel) ;
2211
2212 // Clear caches
2213 _cacheMgr.reset() ;
2214 }
2215
2216 // If the replaced server is an observable that is used in any of the
2217 // normalization sets for conditional fits, replace the element in the
2218 // normalization set too.
2219 for(std::unique_ptr<RooArgSet> const& normSet : _pdfNSetList) {
2220 for(RooAbsArg * arg : *normSet) {
2221 if(RooAbsArg * newArg = arg->findNewServer(newServerList, nameChange)) {
2222 // Need to do some tricks here because it's not possible to replace in
2223 // an owning RooAbsCollection.
2224 normSet->releaseOwnership();
2225 normSet->replace(*std::unique_ptr<RooAbsArg>{arg}, *newArg->cloneTree());
2226 normSet->takeOwnership();
2227 }
2228 }
2229 }
2230
2231 return RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
2232}
2233
2234void RooProdPdf::CacheElem::writeToStream(std::ostream& os) const {
2235 using namespace RooHelpers;
2236 os << "_partList\n";
2237 os << getColonSeparatedNameString(_partList) << "\n";
2238 os << "_numList\n";
2239 os << getColonSeparatedNameString(_numList) << "\n";
2240 os << "_denList\n";
2241 os << getColonSeparatedNameString(_denList) << "\n";
2242 os << "_ownedList\n";
2243 os << getColonSeparatedNameString(_ownedList) << "\n";
2244 os << "_normList\n";
2245 for(auto const& set : _normList) {
2246 os << getColonSeparatedNameString(*set) << "\n";
2247 }
2248 os << "_isRearranged" << "\n";
2249 os << _isRearranged << "\n";
2250 os << "_rearrangedNum" << "\n";
2251 if(_rearrangedNum) {
2252 os << getColonSeparatedNameString(*_rearrangedNum) << "\n";
2253 } else {
2254 os << "nullptr" << "\n";
2255 }
2256 os << "_rearrangedDen" << "\n";
2257 if(_rearrangedDen) {
2258 os << getColonSeparatedNameString(*_rearrangedDen) << "\n";
2259 } else {
2260 os << "nullptr" << "\n";
2261 }
2262}
2263
2264std::unique_ptr<RooArgSet> RooProdPdf::fillNormSetForServer(RooArgSet const &normSet, RooAbsArg const &server) const
2265{
2266 if (normSet.empty())
2267 return nullptr;
2268 auto *pdfNset = findPdfNSet(static_cast<RooAbsPdf const &>(server));
2269 if (pdfNset && !pdfNset->empty()) {
2270 std::unique_ptr<RooArgSet> out;
2271 if (0 == strcmp("cset", pdfNset->GetName())) {
2272 // If the name of the normalization set is "cset", it doesn't contain the
2273 // normalization set but the conditional observables that should *not* be
2274 // normalized over.
2275 out = std::make_unique<RooArgSet>(normSet);
2276 RooArgSet common;
2277 out->selectCommon(*pdfNset, common);
2278 out->remove(common);
2279 } else {
2280 out = std::make_unique<RooArgSet>(*pdfNset);
2281 }
2282 // prefix also the arguments in the normSets if they have not already been
2283 if (auto prefix = getStringAttribute("__prefix__")) {
2284 for (RooAbsArg *arg : *out) {
2285 if (!arg->getStringAttribute("__prefix__")) {
2286 arg->SetName((std::string(prefix) + arg->GetName()).c_str());
2287 arg->setStringAttribute("__prefix__", prefix);
2288 }
2289 }
2290 }
2291 return out;
2292 } else {
2293 return nullptr;
2294 }
2295}
2296
2297/// A RooProdPdf with a fixed normalization set can be replaced by this class.
2298/// Its purpose is to provide the right client-server interface for the
2299/// evaluation of RooProdPdf cache elements that were created for a given
2300/// normalization set.
2302public:
2303 RooFixedProdPdf(std::unique_ptr<RooProdPdf> &&prodPdf, RooArgSet const &normSet)
2304 : RooAbsPdf(prodPdf->GetName(), prodPdf->GetTitle()), _normSet{normSet},
2305 _servers("!servers", "List of servers", this), _prodPdf{std::move(prodPdf)}
2306 {
2307 initialize();
2308 }
2309 RooFixedProdPdf(const RooFixedProdPdf &other, const char *name = nullptr)
2310 : RooAbsPdf(other, name), _normSet{other._normSet},
2311 _servers("!servers", "List of servers", this), _prodPdf{static_cast<RooProdPdf *>(other._prodPdf->Clone())}
2312 {
2313 initialize();
2314 }
2315 TObject *clone(const char *newname) const override { return new RooFixedProdPdf(*this, newname); }
2316
2317 bool selfNormalized() const override { return true; }
2318
2319 inline bool canComputeBatchWithCuda() const override { return true; }
2320
2321 void computeBatch(cudaStream_t *stream, double *output, size_t nEvents,
2322 RooFit::Detail::DataMap const &dataMap) const override
2323 {
2324 _prodPdf->calculateBatch(*_cache, stream, output, nEvents, dataMap);
2325 }
2326
2327 ExtendMode extendMode() const override { return _prodPdf->extendMode(); }
2328 double expectedEvents(const RooArgSet * /*nset*/) const override { return _prodPdf->expectedEvents(&_normSet); }
2329
2330 // Analytical Integration handling
2331 bool forceAnalyticalInt(const RooAbsArg &dep) const override { return _prodPdf->forceAnalyticalInt(dep); }
2332 Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet,
2333 const char *rangeName = nullptr) const override
2334 {
2335 return _prodPdf->getAnalyticalIntegralWN(allVars, analVars, normSet, rangeName);
2336 }
2337 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &numVars, const char *rangeName = nullptr) const override
2338 {
2339 return _prodPdf->getAnalyticalIntegral(allVars, numVars, rangeName);
2340 }
2341 double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName) const override
2342 {
2343 return _prodPdf->analyticalIntegralWN(code, normSet, rangeName);
2344 }
2345 double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override
2346 {
2347 return _prodPdf->analyticalIntegral(code, rangeName);
2348 }
2349
2350private:
2352 {
2353 _cache = _prodPdf->createCacheElem(&_normSet, nullptr);
2354 auto &cache = *_cache;
2355
2356 // The actual servers for a given normalization set depend on whether the
2357 // cache is rearranged or not. See RooProdPdf::calculateBatch to see
2358 // which args in the cache are used directly.
2359 if (cache._isRearranged) {
2360 _servers.add(*cache._rearrangedNum);
2361 _servers.add(*cache._rearrangedDen);
2362 } else {
2363 for (std::size_t i = 0; i < cache._partList.size(); ++i) {
2364 _servers.add(cache._partList[i]);
2365 }
2366 }
2367 }
2368
2369 double evaluate() const override { return _prodPdf->calculate(*_cache); }
2370
2372 std::unique_ptr<RooProdPdf::CacheElem> _cache;
2374 std::unique_ptr<RooProdPdf> _prodPdf;
2375};
2376
2377std::unique_ptr<RooAbsArg>
2379{
2380 std::unique_ptr<RooProdPdf> prodPdfClone{static_cast<RooProdPdf *>(this->Clone())};
2381 prodPdfClone->setAttribute("_COMPILED");
2382
2383 RooArgList serverClones;
2384 for (const auto server : prodPdfClone->servers()) {
2385 auto nsetForServer = fillNormSetForServer(normSet, *server);
2386 RooArgSet const &nset = nsetForServer ? *nsetForServer : normSet;
2387
2388 auto depList = new RooArgSet; // INTENTIONAL LEAK FOR NOW!
2389 server->getObservables(&nset, *depList);
2390
2391 if (auto serverClone = ctx.compile(*server, *prodPdfClone, *depList)) {
2392 serverClones.add(*serverClone);
2393 }
2394 }
2395 prodPdfClone->redirectServers(serverClones, false, true);
2396
2397 auto fixedProdPdf = std::make_unique<RooFixedProdPdf>(std::move(prodPdfClone), normSet);
2398 fixedProdPdf->setAttribute("_COMPILED");
2399
2400 return fixedProdPdf;
2401}
#define coutI(a)
#define cxcoutD(a)
#define coutW(a)
#define dologD(a)
#define coutF(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
#define ClassImp(name)
Definition Rtypes.h:375
static void indent(ostringstream &buf, int indent_level)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void input
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2468
const std::vector< Int_t > & retrieve(Int_t masterCode) const
Retrieve the array of integer codes associated with the given master code.
Int_t store(const std::vector< Int_t > &codeList, RooArgSet *set1=nullptr, RooArgSet *set2=nullptr, RooArgSet *set3=nullptr, RooArgSet *set4=nullptr)
Store given arrays of integer codes, and up to four RooArgSets in the registry (each setX pointer may...
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:78
RooExpensiveObjectCache & expensiveObjectCache() const
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
void setOperMode(OperMode mode, bool recurseADirty=true)
Set the operation mode of this node.
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition RooAbsArg.h:512
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
friend class RooRealIntegral
Definition RooAbsArg.h:642
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:90
RooArgSet * getComponents() const
Create a RooArgSet with all components (branch nodes) of the expression tree headed by this object.
void treeNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=nullptr, bool doBranch=true, bool doLeaf=true, bool valueOnly=false, bool recurseNonDerived=false) const
Fill supplied list with nodes of the arg tree, following all server links, starting with ourself as t...
RooAbsArg * findNewServer(const RooAbsCollection &newSet, bool nameChange) const
Find the new server in the specified set that matches the old server.
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:493
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
Storage_t const & get() const
Const access to the underlying stl container.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
void printValue(std::ostream &os) const override
Print value of collection, i.e.
Int_t getSize() const
Return the number of elements in the collection.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
bool overlaps(Iterator_t otherCollBegin, Iterator_t otherCollEnd) const
Storage_t::size_type size() const
RooAbsArg * first() const
bool selectCommon(const RooAbsCollection &refColl, RooAbsCollection &outColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
TString _normRange
Normalization range.
Definition RooAbsPdf.h:391
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:368
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:271
@ CanNotBeExtended
Definition RooAbsPdf.h:265
const char * normRange() const
Definition RooAbsPdf.h:301
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const
Interface function to create a generator context from a p.d.f.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:105
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:490
TString integralNameSuffix(const RooArgSet &iset, const RooArgSet *nset=nullptr, const char *rangeName=nullptr, bool omitEmpty=false) const
Construct std::string with unique suffix name to give to integral object that encodes integrated obse...
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition RooAddition.h:27
const RooArgList & list1() const
Definition RooAddition.h:41
static TClass * Class()
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
RooArgSet selectFromSet1(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 1.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 2.
void reset()
Clear the cache.
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition RooCmdArg.h:26
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false) override
Remove object 'var' from set and deregister 'var' as server to owner.
TObject * clone(const char *newname) const override
Definition RooConstVar.h:32
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, bool verbose=false)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
T * compile(T &arg, RooAbsArg &owner, RooArgSet const &normSet)
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
A RooProdPdf with a fixed normalization set can be replaced by this class.
std::unique_ptr< RooProdPdf::CacheElem > _cache
void computeBatch(cudaStream_t *stream, double *output, size_t nEvents, RooFit::Detail::DataMap const &dataMap) const override
Base function for computing multiple values of a RooAbsReal.
RooSetProxy _servers
RooFixedProdPdf(const RooFixedProdPdf &other, const char *name=nullptr)
RooFixedProdPdf(std::unique_ptr< RooProdPdf > &&prodPdf, RooArgSet const &normSet)
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &numVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
TObject * clone(const char *newname) const override
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
std::unique_ptr< RooProdPdf > _prodPdf
double expectedEvents(const RooArgSet *) const override
Return expected number of events to be used in calculation of extended likelihood.
bool canComputeBatchWithCuda() const override
bool forceAnalyticalInt(const RooAbsArg &dep) const override
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName) const override
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
ExtendMode extendMode() const override
Returns ability of PDF to provide extended likelihood terms.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooAbsArg * getParameter(const char *name) const
Return pointer to parameter with given name.
General form of projected integral of product of PDFs, utility class for RooProdPdf.
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
TObject * At(int index) const
Return object stored in sequential position given by index.
RooLinkedListIterImpl end() const
void Delete(Option_t *o=nullptr) override
Remove all elements in collection and delete all elements NB: Collection does not own elements,...
virtual void Add(TObject *arg)
RooLinkedListIterImpl begin() const
Int_t IndexOf(const char *name) const
Return position of given object in list.
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:37
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
RooArgList containedArgs(Action) override
Return RooAbsArg components contained in the cache.
std::unique_ptr< RooAbsReal > _rearrangedNum
Definition RooProdPdf.h:146
void printCompactTreeHook(std::ostream &, const char *, Int_t, Int_t) override
Hook function to print cache contents in tree printing of RooProdPdf.
std::vector< std::unique_ptr< RooArgSet > > _normList
Definition RooProdPdf.h:144
std::unique_ptr< RooAbsReal > _rearrangedDen
Definition RooProdPdf.h:147
void writeToStream(std::ostream &os) const
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components of a RooProdPdf with cache-and-track, and label all RooProdPdf descendants wit...
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Query internal generation capabilities of component p.d.f.s and aggregate capabilities into master co...
void rearrangeProduct(CacheElem &) const
void factorizeProduct(const RooArgSet &normSet, const RooArgSet &intSet, RooLinkedList &termList, RooLinkedList &normList, RooLinkedList &impDepList, RooLinkedList &crossDepList, RooLinkedList &intList) const
Factorize product in irreducible terms for given choice of integration/normalization.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral defined by given scenario code.
RooArgSet * getConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, bool stripDisconnected, bool removeConstraintsFromPdf=false) const override
Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
~RooProdPdf() override
Destructor.
Int_t _extendedIndex
Index of extended PDF (if any)
Definition RooProdPdf.h:177
RooProdPdf()
Default constructor.
void removePdfs(RooAbsCollection const &pdfs)
Remove some PDFs from the factors of this RooProdPdf.
bool _useDefaultGen
Use default or distributed event generator.
Definition RooProdPdf.h:180
std::vector< std::unique_ptr< RooArgSet > > _pdfNSetList
List of PDF component normalization sets.
Definition RooProdPdf.h:176
std::unique_ptr< RooArgSet > fillNormSetForServer(RooArgSet const &normSet, RooAbsArg const &server) const
bool forceAnalyticalInt(const RooAbsArg &dep) const override
Force RooRealIntegral to offer all observables for internal integration.
void calculateBatch(const RooProdPdf::CacheElem &cache, cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const
Evaluate product of PDFs in batch mode.
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
void getParametersHook(const RooArgSet *, RooArgSet *, bool stripDisconnected) const override
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Return generator context optimized for generating events from product p.d.f.s.
TNamed * _refRangeName
Reference range name for interpretation of conditional products.
Definition RooProdPdf.h:182
RooAICRegistry _genCode
! Registry of composite direct generator codes
Definition RooProdPdf.h:172
void addPdfs(RooAbsCollection const &pdfs)
Add some full PDFs to the factors of this RooProdPdf.
RooListProxy _pdfList
List of PDF components.
Definition RooProdPdf.h:175
Int_t getPartIntList(const RooArgSet *nset, const RooArgSet *iset, const char *isetRangeName=nullptr) const
Return list of (partial) integrals of product terms for integration of p.d.f over observables iset wh...
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooProdPdf to more intuitively reflect the contents of the prod...
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
RooAbsReal * specializeIntegral(RooAbsReal &orig, const char *targetRangeName) const
RooObjCacheManager _cacheMgr
Definition RooProdPdf.h:156
std::string makeRGPPName(const char *pfx, const RooArgSet &term, const RooArgSet &iset, const RooArgSet &nset, const char *isetRangeName) const
Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
bool isDirectGenSafe(const RooAbsArg &arg) const override
Forward determination of safety of internal generator code to component p.d.f that would generate the...
std::unique_ptr< RooAbsReal > makeCondPdfRatioCorr(RooAbsReal &term, const RooArgSet &termNset, const RooArgSet &termImpSet, const char *normRange, const char *refRange) const
For single normalization ranges.
RooArgSet * findPdfNSet(RooAbsPdf const &pdf) const
Look up user specified normalization set for given input PDF component.
ExtendMode extendMode() const override
If this product contains exactly one extendable p.d.f return the extension abilities of that p....
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
double expectedEvents(const RooArgSet *nset) const override
Return the expected number of events associated with the extendable input PDF in the product.
std::vector< RooAbsReal * > processProductTerm(const RooArgSet *nset, const RooArgSet *iset, const char *isetRangeName, const RooArgSet *term, const RooArgSet &termNSet, const RooArgSet &termISet, bool &isOwned, bool forceWrap=false) const
Calculate integrals of factorized product terms over observables iset while normalized to observables...
TObject * clone(const char *newname) const override
Definition RooProdPdf.h:55
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Determine which part (if any) of given integral can be performed analytically.
friend class RooFixedProdPdf
Definition RooProdPdf.h:167
bool isBinnedDistribution(const RooArgSet &obs) const override
If all components that depend on obs are binned that so is the product.
friend class RooProdGenContext
Definition RooProdPdf.h:166
RooArgSet * getConnectedParameters(const RooArgSet &observables) const
Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
double calculate(const RooProdPdf::CacheElem &cache, bool verbose=false) const
Calculate running product of pdfs terms, using the supplied normalization set in 'normSetList' for ea...
RooArgSet _defNormSet
Default normalization set.
Definition RooProdPdf.h:185
RooAbsReal * specializeRatio(RooFormulaVar &input, const char *targetRangeName) const
CacheElem * getCacheElem(RooArgSet const *nset) const
The cache manager.
bool redirectServersHook(const RooAbsCollection &, bool, bool, bool) override
Implement support for node removal.
void initGenerator(Int_t code) override
Forward one-time initialization call to component generation initialization methods.
void generateEvent(Int_t code) override
Generate a single event with configuration specified by 'code' Defer internal generation to component...
void groupProductTerms(std::list< std::vector< RooArgSet * > > &groupedTerms, RooArgSet &outerIntDeps, const RooLinkedList &terms, const RooLinkedList &norms, const RooLinkedList &imps, const RooLinkedList &ints, const RooLinkedList &cross) const
Group product into terms that can be calculated independently.
void fixRefRange(const char *rangeName)
std::unique_ptr< CacheElem > createCacheElem(const RooArgSet *nset, const RooArgSet *iset, const char *isetRangeName=nullptr) const
double evaluate() const override
Calculate current value of object.
void initializeFromCmdArgList(const RooArgSet &fullPdfSet, const RooLinkedList &l)
Initialize RooProdPdf configuration from given list of RooCmdArg configuration arguments and set of '...
void selectNormalizationRange(const char *rangeName=nullptr, bool force=false) override
Interface function used by test statistics to freeze choice of range for interpretation of conditiona...
double _cutOff
Cutoff parameter for running product.
Definition RooProdPdf.h:174
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
static TClass * Class()
RooArgList components()
Definition RooProduct.h:48
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooArgSet intVars() const
const RooAbsReal & integrand() const
static TClass * Class()
const char * intRange() const
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Mother of all ROOT objects.
Definition TObject.h:41
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:421
const char * Data() const
Definition TString.h:380
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
RooConstVar & RooConst(double val)
Double_t x[n]
Definition legend1.C:17
std::vector< RooSpan< const double > > VarVector
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
std::vector< double > ArgVector
std::string getColonSeparatedNameString(RooArgSet const &argSet)
Create a string with all sorted names of RooArgSet elements separated by colons.
Definition first.py:1
TLine l
Definition textangle.C:4
static void output()