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
22Efficient 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 "RooFitImplHelpers.h"
63#include "strtok.h"
64
65#include <algorithm>
66#include <array>
67#include <cstring>
68#include <sstream>
69
70#ifndef _WIN32
71#include <strings.h>
72#endif
73
74using std::endl, std::string, std::vector, std::list, std::ostream, std::map, std::ostringstream;
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 analogous 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, nullptr, &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, nullptr, &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->empty() ? 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.
414{
415 if (cache._isRearranged) {
416 auto numerator = ctx.at(cache._rearrangedNum.get());
417 auto denominator = ctx.at(cache._rearrangedDen.get());
418 RooBatchCompute::compute(ctx.config(caller), RooBatchCompute::Ratio, ctx.output(), {numerator, denominator});
419 } else {
420 std::vector<std::span<const double>> factors;
421 factors.reserve(cache._partList.size());
422 for (const RooAbsArg *i : cache._partList) {
423 auto span = ctx.at(i);
424 factors.push_back(span);
425 }
426 std::array<double, 1> special{static_cast<double>(factors.size())};
427 RooBatchCompute::compute(ctx.config(caller), RooBatchCompute::ProdPdf, ctx.output(), factors, special);
428 }
429}
430
431namespace {
432
433template<class T>
434void eraseNullptrs(std::vector<T*>& v) {
435 v.erase(std::remove_if(v.begin(), v.end(), [](T* x){ return x == nullptr; } ), v.end());
436}
437
438void removeCommon(std::vector<RooAbsArg*> &v, std::span<RooAbsArg * const> other) {
439
440 for (auto const& arg : other) {
441 auto namePtrMatch = [&arg](const RooAbsArg* elm) {
442 return elm != nullptr && elm->namePtr() == arg->namePtr();
443 };
444
445 auto found = std::find_if(v.begin(), v.end(), namePtrMatch);
446 if(found != v.end()) {
447 *found = nullptr;
448 }
449 }
450 eraseNullptrs(v);
451}
452
453void addCommon(std::vector<RooAbsArg*> &v, std::vector<RooAbsArg*> const& o1, std::vector<RooAbsArg*> const& o2) {
454
455 for (auto const& arg : o1) {
456 auto namePtrMatch = [&arg](const RooAbsArg* elm) {
457 return elm->namePtr() == arg->namePtr();
458 };
459
460 if(std::find_if(o2.begin(), o2.end(), namePtrMatch) != o2.end()) {
461 v.push_back(arg);
462 }
463 }
464}
465
466}
467
468
469////////////////////////////////////////////////////////////////////////////////
470/// Factorize product in irreducible terms for given choice of integration/normalization
471
472void RooProdPdf::factorizeProduct(const RooArgSet& normSet, const RooArgSet& intSet,
473 RooLinkedList& termList, RooLinkedList& normList,
474 RooLinkedList& impDepList, RooLinkedList& crossDepList,
475 RooLinkedList& intList) const
476{
477 // List of all term dependents: normalization and imported
478 std::vector<RooArgSet> depAllList;
479 std::vector<RooArgSet> depIntNoNormList;
480
481 // Setup lists for factorization terms and their dependents
482 RooArgSet* term(nullptr);
483 RooArgSet* termNormDeps(nullptr);
484 RooArgSet* termIntDeps(nullptr);
485 RooArgSet* termIntNoNormDeps(nullptr);
486
487 std::vector<RooAbsArg*> pdfIntNoNormDeps;
488 std::vector<RooAbsArg*> pdfIntSet;
489 std::vector<RooAbsArg*> pdfNSet;
490 std::vector<RooAbsArg*> pdfCSet;
491 std::vector<RooAbsArg*> pdfNormDeps; // Dependents to be normalized for the PDF
492 std::vector<RooAbsArg*> pdfAllDeps; // All dependents of this PDF
493
494 // Loop over the PDFs
495 for(std::size_t iPdf = 0; iPdf < _pdfList.size(); ++iPdf) {
496 RooAbsPdf& pdf = static_cast<RooAbsPdf&>(_pdfList[iPdf]);
497 RooArgSet& pdfNSetOrig = *_pdfNSetList[iPdf];
498
499 pdfNSet.clear();
500 pdfCSet.clear();
501
502 // Make iterator over tree leaf node list to get the observables.
503 // This code is borrowed from RooAbsPdf::getObservables().
504 // RooAbsArg::treeNodeServer list is relatively expensive, so we only do it
505 // once and use it in a lambda function.
506 RooArgSet pdfLeafList("leafNodeServerList") ;
507 pdf.treeNodeServerList(&pdfLeafList,nullptr,false,true,true) ;
508 auto getObservablesOfCurrentPdf = [&pdfLeafList](
509 std::vector<RooAbsArg*> & out,
510 const RooArgSet& dataList) {
511 for (const auto arg : pdfLeafList) {
512 if (arg->dependsOnValue(dataList) && arg->isLValue()) {
513 out.push_back(arg) ;
514 }
515 }
516 };
517
518 // Reduce pdfNSet to actual dependents
519 if (0 == strcmp("cset", pdfNSetOrig.GetName())) {
520 getObservablesOfCurrentPdf(pdfNSet, normSet);
521 removeCommon(pdfNSet, pdfNSetOrig.get());
522 pdfCSet = pdfNSetOrig.get();
523 } else {
524 // Interpret at NSet
525 getObservablesOfCurrentPdf(pdfNSet, pdfNSetOrig);
526 }
527
528
529 pdfNormDeps.clear();
530 pdfAllDeps.clear();
531
532 // Make list of all dependents of this PDF
533 getObservablesOfCurrentPdf(pdfAllDeps, normSet);
534
535
536// cout << GetName() << ": pdf = " << pdf->GetName() << " pdfAllDeps = " << pdfAllDeps << " pdfNSet = " << *pdfNSet << " pdfCSet = " << *pdfCSet << endl;
537
538 // Make list of normalization dependents for this PDF;
539 if (!pdfNSet.empty()) {
540 // PDF is conditional
541 addCommon(pdfNormDeps, pdfAllDeps, pdfNSet);
542 } else {
543 // PDF is regular
544 pdfNormDeps = pdfAllDeps;
545 }
546
547// cout << GetName() << ": pdfNormDeps for " << pdf->GetName() << " = " << pdfNormDeps << endl;
548
549 pdfIntSet.clear();
550 getObservablesOfCurrentPdf(pdfIntSet, intSet) ;
551
552 // WVE if we have no norm deps, conditional observables should be taken out of pdfIntSet
553 if (pdfNormDeps.empty() && !pdfCSet.empty()) {
554 removeCommon(pdfIntSet, pdfCSet);
555// cout << GetName() << ": have no norm deps, removing conditional observables from intset" << endl;
556 }
557
558 pdfIntNoNormDeps.clear();
559 pdfIntNoNormDeps = pdfIntSet;
560 removeCommon(pdfIntNoNormDeps, pdfNormDeps);
561
562// cout << GetName() << ": pdf = " << pdf->GetName() << " intset = " << *pdfIntSet << " pdfIntNoNormDeps = " << pdfIntNoNormDeps << endl;
563
564 // Check if this PDF has dependents overlapping with one of the existing terms
565 bool done = false;
566 int j = 0;
567 auto lIter = termList.begin();
568 auto ldIter = normList.begin();
569 for(;lIter != termList.end(); (++lIter, ++ldIter, ++j)) {
570 termNormDeps = static_cast<RooArgSet*>(*ldIter);
571 term = static_cast<RooArgSet*>(*lIter);
572 // PDF should be added to existing term if
573 // 1) It has overlapping normalization dependents with any other PDF in existing term
574 // 2) It has overlapping dependents of any class for which integration is requested
575 // 3) If normalization happens over multiple ranges, and those ranges are both defined
576 // in either observable
577
578 bool normOverlap = termNormDeps->overlaps(pdfNormDeps.begin(), pdfNormDeps.end());
579 //bool intOverlap = pdfIntSet->overlaps(*termAllDeps);
580
581 if (normOverlap) {
582// cout << GetName() << ": this term overlaps with term " << (*term) << " in normalization observables" << endl;
583
584 term->add(pdf);
585 termNormDeps->add(pdfNormDeps.begin(), pdfNormDeps.end(), false);
586 depAllList[j].add(pdfAllDeps.begin(), pdfAllDeps.end(), false);
587 if (termIntDeps) {
588 termIntDeps->add(pdfIntSet.begin(), pdfIntSet.end(), false);
589 }
590 if (termIntNoNormDeps) {
591 termIntNoNormDeps->add(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), false);
592 }
593 termIntNoNormDeps->add(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), false);
594 done = true;
595 }
596 }
597
598 // If not, create a new term
599 if (!done) {
600 if (!(pdfNormDeps.empty() && pdfAllDeps.empty() &&
601 pdfIntSet.empty()) || normSet.empty()) {
602 term = new RooArgSet("term");
603 termNormDeps = new RooArgSet("termNormDeps");
604 depAllList.emplace_back(pdfAllDeps.begin(), pdfAllDeps.end(), "termAllDeps");
605 termIntDeps = new RooArgSet(pdfIntSet.begin(), pdfIntSet.end(), "termIntDeps");
606 depIntNoNormList.emplace_back(pdfIntNoNormDeps.begin(), pdfIntNoNormDeps.end(), "termIntNoNormDeps");
607 termIntNoNormDeps = &depIntNoNormList.back();
608
609 term->add(pdf);
610 termNormDeps->add(pdfNormDeps.begin(), pdfNormDeps.end(), false);
611
612 termList.Add(term);
613 normList.Add(termNormDeps);
614 intList.Add(termIntDeps);
615 }
616 }
617
618 }
619
620 // Loop over list of terms again to determine 'imported' observables
621 int i = 0;
622 RooArgSet *normDeps;
623 auto lIter = termList.begin();
624 auto ldIter = normList.begin();
625 for(;lIter != termList.end(); (++lIter, ++ldIter, ++i)) {
626 normDeps = static_cast<RooArgSet*>(*ldIter);
627 term = static_cast<RooArgSet*>(*lIter);
628 // Make list of wholly imported dependents
629 RooArgSet impDeps(depAllList[i]);
630 impDeps.remove(*normDeps, true, true);
631 auto snap = new RooArgSet;
632 impDeps.snapshot(*snap);
633 impDepList.Add(snap);
634// cout << GetName() << ": list of imported dependents for term " << (*term) << " set to " << impDeps << endl ;
635
636 // Make list of cross dependents (term is self contained for these dependents,
637 // but components import dependents from other components)
638 auto crossDeps = std::unique_ptr<RooAbsCollection>{depIntNoNormList[i].selectCommon(*normDeps)};
639 snap = new RooArgSet;
640 crossDeps->snapshot(*snap);
641 crossDepList.Add(snap);
642// cout << GetName() << ": list of cross dependents for term " << (*term) << " set to " << *crossDeps << endl ;
643 }
644
645 return;
646}
647
648
649
650
651////////////////////////////////////////////////////////////////////////////////
652/// Return list of (partial) integrals of product terms for integration
653/// of p.d.f over observables iset while normalization over observables nset.
654/// Also return list of normalization sets to be used to evaluate
655/// each component in the list correctly.
656
657Int_t RooProdPdf::getPartIntList(const RooArgSet* nset, const RooArgSet* iset, const char* isetRangeName) const
658{
659 // Check if this configuration was created before
660 Int_t sterileIdx(-1);
661
662 if (static_cast<CacheElem*>(_cacheMgr.getObj(nset,iset,&sterileIdx,isetRangeName))) {
663 return _cacheMgr.lastIndex();
664 }
665
666 std::unique_ptr<CacheElem> cache = createCacheElem(nset, iset, isetRangeName);
667
668 // Store the partial integral list and return the assigned code
669 return _cacheMgr.setObj(nset, iset, cache.release(), RooNameReg::ptr(isetRangeName));
670}
671
672
673
674std::unique_ptr<RooProdPdf::CacheElem> RooProdPdf::createCacheElem(const RooArgSet* nset,
675 const RooArgSet* iset,
676 const char* isetRangeName) const
677{
678// cout << " FOLKERT::RooProdPdf::getPartIntList(" << GetName() <<") nset = " << (nset?*nset:RooArgSet()) << endl
679// << " _normRange = " << _normRange << endl
680// << " iset = " << (iset?*iset:RooArgSet()) << endl
681// << " isetRangeName = " << (isetRangeName?isetRangeName:"<null>") << endl ;
682
683 // Create containers for partial integral components to be generated
684 auto cache = std::make_unique<CacheElem>();
685
686 // Factorize the product in irreducible terms for this nset
687 RooLinkedList terms;
688 RooLinkedList norms;
689 RooLinkedList imp;
690 RooLinkedList ints;
691 RooLinkedList cross;
692 // cout << "RooProdPdf::getPIL -- now calling factorizeProduct()" << endl ;
693
694
695 // Normalization set used for factorization
696 RooArgSet factNset(nset ? (*nset) : _defNormSet);
697// cout << GetName() << "factNset = " << factNset << endl ;
698
699 factorizeProduct(factNset, iset ? (*iset) : RooArgSet(), terms, norms, imp, cross, ints);
700
701 RooArgSet *norm;
702 RooArgSet *integ;
703 RooArgSet *xdeps;
704 RooArgSet *imps;
705
706 // Group irriducible terms that need to be (partially) integrated together
707 std::list<std::vector<RooArgSet*>> groupedList;
708 RooArgSet outerIntDeps;
709// cout << "RooProdPdf::getPIL -- now calling groupProductTerms()" << endl;
710 groupProductTerms(groupedList, outerIntDeps, terms, norms, imp, ints, cross);
711
712 // Loop over groups
713// cout<<"FK: pdf("<<GetName()<<") Starting selecting F(x|y)!"<<endl;
714 // Find groups of type F(x|y), i.e. termImpSet!=0, construct ratio object
715 std::map<std::string, RooArgSet> ratioTerms;
716 for (auto const& group : groupedList) {
717 if (1 == group.size()) {
718// cout<<"FK: Starting Single Term"<<endl;
719
720 RooArgSet* term = group[0];
721
722 Int_t termIdx = terms.IndexOf(term);
723 norm=static_cast<RooArgSet*>(norms.At(termIdx));
724 imps=static_cast<RooArgSet*>(imp.At(termIdx));
725 RooArgSet termNSet(*norm);
726 RooArgSet termImpSet(*imps);
727
728 // cout<<"FK: termImpSet.size() = "<<termImpSet.size()<< " " << termImpSet << endl;
729 // cout<<"FK: _refRangeName = "<<_refRangeName<<endl;
730
731 if (!termImpSet.empty() && nullptr != _refRangeName) {
732
733// cout << "WVE now here" << endl;
734
735 // WVE we can skip this if the ref range is equal to the normalization range
736 bool rangeIdentical(true);
737// cout << "_normRange = " << _normRange << " _refRangeName = " << RooNameReg::str(_refRangeName) << endl ;
738 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
739 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
740 if (_normRange.Length() > 0) {
741 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
742 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
743 }
744 else{
745 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
746 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
747 }
748 }
749// cout<<"FK: rangeIdentical Single = "<<(rangeIdentical ? 'T':'F')<<endl;
750 // coverity[CONSTANT_EXPRESSION_RESULT]
751 // LM : avoid making integral ratio if range is the same. Why was not included ??? (same at line 857)
752 if (!rangeIdentical ) {
753// cout << "PREPARING RATIO HERE (SINGLE TERM)" << endl ;
754 auto ratio = makeCondPdfRatioCorr(*static_cast<RooAbsReal*>(term->first()), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
755 std::ostringstream str; termImpSet.printValue(str);
756// cout << GetName() << "inserting ratio term" << endl;
757 ratioTerms[str.str()].addOwned(std::move(ratio));
758 }
759 }
760
761 } else {
762// cout<<"FK: Starting Composite Term"<<endl;
763
764 for (auto const& term : group) {
765
766 Int_t termIdx = terms.IndexOf(term);
767 norm=static_cast<RooArgSet*>(norms.At(termIdx));
768 imps=static_cast<RooArgSet*>(imp.At(termIdx));
769 RooArgSet termNSet(*norm);
770 RooArgSet termImpSet(*imps);
771
772 if (!termImpSet.empty() && nullptr != _refRangeName) {
773
774 // WVE we can skip this if the ref range is equal to the normalization range
775 bool rangeIdentical(true);
776 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
777 if(_normRange.Length() > 0) {
778 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
779 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
780 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
781 }
782 } else {
783 for (auto const* normObs : static_range_cast<RooRealVar*>(termNSet)) {
784 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = false;
785 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = false;
786 }
787 }
788// cout<<"FK: rangeIdentical Composite = "<<(rangeIdentical ? 'T':'F') <<endl;
789 if (!rangeIdentical ) {
790// cout << "PREPARING RATIO HERE (COMPOSITE TERM)" << endl ;
791 auto ratio = makeCondPdfRatioCorr(*static_cast<RooAbsReal*>(term->first()), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
792 std::ostringstream str; termImpSet.printValue(str);
793 ratioTerms[str.str()].addOwned(std::move(ratio));
794 }
795 }
796 }
797 }
798
799 }
800
801 // Find groups with y as termNSet
802 // Replace G(y) with (G(y),ratio)
803 for (auto const& group : groupedList) {
804 for (auto const& term : group) {
805 Int_t termIdx = terms.IndexOf(term);
806 norm = static_cast<RooArgSet*>(norms.At(termIdx));
807 imps = static_cast<RooArgSet*>(imp.At(termIdx));
808 RooArgSet termNSet(*norm);
809 RooArgSet termImpSet(*imps);
810
811 // If termNset matches index of ratioTerms, insert ratio here
812 ostringstream str; termNSet.printValue(str);
813 if (!ratioTerms[str.str()].empty()) {
814// cout << "MUST INSERT RATIO OBJECT IN TERM (COMPOSITE)" << *term << endl;
815 term->add(ratioTerms[str.str()]);
816 cache->_ownedList.addOwned(std::move(ratioTerms[str.str()]));
817 }
818 }
819 }
820
821 for (auto const& group : groupedList) {
822// cout << GetName() << ":now processing group" << endl;
823// group->Print("1");
824
825 if (1 == group.size()) {
826// cout << "processing atomic item" << endl;
827 RooArgSet* term = group[0];
828
829 Int_t termIdx = terms.IndexOf(term);
830 norm = static_cast<RooArgSet*>(norms.At(termIdx));
831 integ = static_cast<RooArgSet*>(ints.At(termIdx));
832 xdeps = static_cast<RooArgSet*>(cross.At(termIdx));
833 imps = static_cast<RooArgSet*>(imp.At(termIdx));
834
835 RooArgSet termNSet;
836 RooArgSet termISet;
837 RooArgSet termXSet;
838 RooArgSet termImpSet;
839
840 // Take list of normalization, integrated dependents from factorization algorithm
841 termISet.add(*integ);
842 termNSet.add(*norm);
843
844 // Cross-imported integrated dependents
845 termXSet.add(*xdeps);
846 termImpSet.add(*imps);
847
848 // Add prefab term to partIntList.
849 bool isOwned(false);
850 vector<RooAbsReal*> func = processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned);
851 if (func[0]) {
852 cache->_partList.add(*func[0]);
853 if (isOwned) cache->_ownedList.addOwned(std::unique_ptr<RooAbsArg>{func[0]});
854
855 cache->_normList.emplace_back(std::make_unique<RooArgSet>());
856 norm->snapshot(*cache->_normList.back(), false);
857
858 cache->_numList.addOwned(std::unique_ptr<RooAbsArg>{func[1]});
859 cache->_denList.addOwned(std::unique_ptr<RooAbsArg>{func[2]});
860 }
861 } else {
862// cout << "processing composite item" << endl;
863 RooArgSet compTermSet;
864 RooArgSet compTermNorm;
865 RooArgSet compTermNum;
866 RooArgSet compTermDen;
867 for (auto const &term : group) {
868 // cout << GetName() << ": processing term " << (*term) << " of composite item" << endl ;
869 Int_t termIdx = terms.IndexOf(term);
870 norm = static_cast<RooArgSet *>(norms.At(termIdx));
871 integ = static_cast<RooArgSet *>(ints.At(termIdx));
872 xdeps = static_cast<RooArgSet *>(cross.At(termIdx));
873 imps = static_cast<RooArgSet *>(imp.At(termIdx));
874
875 RooArgSet termNSet;
876 RooArgSet termISet;
877 RooArgSet termXSet;
878 RooArgSet termImpSet;
879 termISet.add(*integ);
880 termNSet.add(*norm);
881 termXSet.add(*xdeps);
882 termImpSet.add(*imps);
883
884 // Remove outer integration dependents from termISet
885 termISet.remove(outerIntDeps, true, true);
886
887 bool isOwned = false;
888 vector<RooAbsReal *> func =
889 processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned, true);
890 // cout << GetName() << ": created composite term component " << func[0]->GetName() << endl;
891 if (func[0]) {
892 compTermSet.add(*func[0]);
893 if (isOwned) cache->_ownedList.addOwned(std::unique_ptr<RooAbsArg>{func[0]});
894 compTermNorm.add(*norm, false);
895
896 compTermNum.add(*func[1]);
897 compTermDen.add(*func[2]);
898 //cache->_numList.add(*func[1]);
899 //cache->_denList.add(*func[2]);
900
901 }
902 }
903
904// cout << GetName() << ": constructing special composite product" << endl;
905// cout << GetName() << ": compTermSet = " ; compTermSet.Print("1");
906
907 // WVE THIS NEEDS TO BE REARRANGED
908
909 // compTermset is set van partial integrals to be multiplied
910 // prodtmp = product (compTermSet)
911 // inttmp = int ( prodtmp ) d (outerIntDeps) _range_isetRangeName
912
913 const std::string prodname = makeRGPPName("SPECPROD", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
914 auto prodtmp = std::make_unique<RooProduct>(prodname.c_str(), prodname.c_str(), compTermSet);
915
916 const std::string intname = makeRGPPName("SPECINT", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
917 auto inttmp = std::make_unique<RooRealIntegral>(intname.c_str(), intname.c_str(), *prodtmp, outerIntDeps, nullptr, nullptr, isetRangeName);
918 inttmp->setStringAttribute("PROD_TERM_TYPE", "SPECINT");
919
920 cache->_partList.add(*inttmp);
921
922 // Product of numerator terms
923 const string prodname_num = makeRGPPName("SPECPROD_NUM", compTermNum, RooArgSet(), RooArgSet(), nullptr);
924 auto prodtmp_num = std::make_unique<RooProduct>(prodname_num.c_str(), prodname_num.c_str(), compTermNum);
925 prodtmp_num->addOwnedComponents(compTermNum);
926
927 // Product of denominator terms
928 const string prodname_den = makeRGPPName("SPECPROD_DEN", compTermDen, RooArgSet(), RooArgSet(), nullptr);
929 auto prodtmp_den = std::make_unique<RooProduct>(prodname_den.c_str(), prodname_den.c_str(), compTermDen);
930 prodtmp_den->addOwnedComponents(compTermDen);
931
932 // Ratio
933 std::string name = Form("SPEC_RATIO(%s,%s)", prodname_num.c_str(), prodname_den.c_str());
934 auto ndr = std::make_unique<RooFormulaVar>(name.c_str(), "@0/@1", RooArgList(*prodtmp_num, *prodtmp_den));
935
936 // Integral of ratio
937 std::unique_ptr<RooAbsReal> numtmp{ndr->createIntegral(outerIntDeps,isetRangeName)};
938 numtmp->addOwnedComponents(std::move(ndr));
939
940 cache->_ownedList.addOwned(std::move(prodtmp));
941 cache->_ownedList.addOwned(std::move(inttmp));
942 cache->_ownedList.addOwned(std::move(prodtmp_num));
943 cache->_ownedList.addOwned(std::move(prodtmp_den));
944 cache->_numList.addOwned(std::move(numtmp));
945 cache->_denList.addOwned(std::unique_ptr<RooAbsArg>{static_cast<RooAbsArg*>(RooFit::RooConst(1).clone("1"))});
946 cache->_normList.emplace_back(std::make_unique<RooArgSet>());
947 compTermNorm.snapshot(*cache->_normList.back(), false);
948 }
949 }
950
951 // Need to rearrange product in case of multiple ranges
952 if (_normRange.Contains(",")) {
953 rearrangeProduct(*cache);
954 }
955
956 // We own contents of all lists filled by factorizeProduct()
957 terms.Delete();
958 ints.Delete();
959 imp.Delete();
960 norms.Delete();
961 cross.Delete();
962
963 return cache;
964}
965
966
967
968////////////////////////////////////////////////////////////////////////////////
969/// For single normalization ranges
970
971std::unique_ptr<RooAbsReal> RooProdPdf::makeCondPdfRatioCorr(RooAbsReal& pdf, const RooArgSet& termNset, const RooArgSet& /*termImpSet*/, const char* normRangeTmp, const char* refRange) const
972{
973 std::unique_ptr<RooAbsReal> ratio_num{pdf.createIntegral(termNset,normRangeTmp)};
974 std::unique_ptr<RooAbsReal> ratio_den{pdf.createIntegral(termNset,refRange)};
975 auto ratio = std::make_unique<RooFormulaVar>(Form("ratio(%s,%s)",ratio_num->GetName(),ratio_den->GetName()),"@0/@1",
976 RooArgList(*ratio_num,*ratio_den)) ;
977
978 ratio->addOwnedComponents(std::move(ratio_num));
979 ratio->addOwnedComponents(std::move(ratio_den));
980 ratio->setAttribute("RATIO_TERM") ;
981 return ratio ;
982}
983
984
985
986
987////////////////////////////////////////////////////////////////////////////////
988
990{
991 RooAbsReal *part;
992 RooAbsReal *num;
993 RooAbsReal *den;
994 RooArgSet nomList ;
995
996 list<string> rangeComps ;
997 {
998 std::vector<char> buf(strlen(_normRange.Data()) + 1);
999 strcpy(buf.data(),_normRange.Data()) ;
1000 char* save(nullptr) ;
1001 char* token = R__STRTOK_R(buf.data(),",",&save) ;
1002 while(token) {
1003 rangeComps.push_back(token) ;
1004 token = R__STRTOK_R(nullptr,",",&save) ;
1005 }
1006 }
1007
1008
1009 std::map<std::string,RooArgSet> denListList ;
1010 RooArgSet specIntDeps ;
1011 string specIntRange ;
1012
1013// cout << "THIS IS REARRANGEPRODUCT" << endl ;
1014
1015 for (std::size_t i = 0; i < cache._partList.size(); i++) {
1016
1017 part = static_cast<RooAbsReal*>(cache._partList.at(i));
1018 num = static_cast<RooAbsReal*>(cache._numList.at(i));
1019 den = static_cast<RooAbsReal*>(cache._denList.at(i));
1020 i++;
1021
1022// cout << "now processing part " << part->GetName() << " of type " << part->getStringAttribute("PROD_TERM_TYPE") << endl ;
1023// cout << "corresponding numerator = " << num->GetName() << endl ;
1024// cout << "corresponding denominator = " << den->GetName() << endl ;
1025
1026
1027 RooFormulaVar* ratio(nullptr) ;
1028 RooArgSet origNumTerm ;
1029
1030 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1031
1032 RooRealIntegral* orig = static_cast<RooRealIntegral*>(num);
1033 auto specratio = static_cast<RooFormulaVar const*>(&orig->integrand()) ;
1034 RooProduct* func = static_cast<RooProduct*>(specratio->getParameter(0)) ;
1035
1036 std::unique_ptr<RooArgSet> components{orig->getComponents()};
1037 for(RooAbsArg * carg : *components) {
1038 if (carg->getAttribute("RATIO_TERM")) {
1039 ratio = static_cast<RooFormulaVar*>(carg) ;
1040 break ;
1041 }
1042 }
1043
1044 if (ratio) {
1045 RooCustomizer cust(*func,"blah") ;
1046 cust.replaceArg(*ratio,RooFit::RooConst(1)) ;
1047 RooAbsArg* funcCust = cust.build() ;
1048// cout << "customized function = " << endl ;
1049// funcCust->printComponentTree() ;
1050 nomList.add(*funcCust) ;
1051 } else {
1052 nomList.add(*func) ;
1053 }
1054
1055
1056 } else {
1057
1058 // Find the ratio term
1059 RooAbsReal* func = num;
1060 // If top level object is integral, navigate to integrand
1061 if (func->InheritsFrom(RooRealIntegral::Class())) {
1062 func = const_cast<RooAbsReal*>(&static_cast<RooRealIntegral*>(func)->integrand());
1063 }
1064 if (func->InheritsFrom(RooProduct::Class())) {
1065// cout << "product term found: " ; func->Print() ;
1066 for(RooAbsArg * arg : static_cast<RooProduct*>(func)->components()) {
1067 if (arg->getAttribute("RATIO_TERM")) {
1068 ratio = static_cast<RooFormulaVar*>(arg) ;
1069 } else {
1070 origNumTerm.add(*arg) ;
1071 }
1072 }
1073 }
1074
1075 if (ratio) {
1076// cout << "Found ratio term in numerator: " << ratio->GetName() << endl ;
1077// cout << "Adding only original term to numerator: " << origNumTerm << endl ;
1078 nomList.add(origNumTerm) ;
1079 } else {
1080 nomList.add(*num) ;
1081 }
1082
1083 }
1084
1085 for (list<string>::iterator iter = rangeComps.begin() ; iter != rangeComps.end() ; ++iter) {
1086 // If denominator is an integral, make a clone with the integration range adjusted to
1087 // the selected component of the normalization integral
1088// cout << "NOW PROCESSING DENOMINATOR " << den->ClassName() << "::" << den->GetName() << endl ;
1089
1090 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1091
1092// cout << "create integral: SPECINT case" << endl ;
1093 RooRealIntegral* orig = static_cast<RooRealIntegral*>(num);
1094 auto specRatio = static_cast<RooFormulaVar const*>(&orig->integrand()) ;
1095 specIntDeps.add(orig->intVars()) ;
1096 if (orig->intRange()) {
1097 specIntRange = orig->intRange() ;
1098 }
1099 //RooProduct* numtmp = (RooProduct*) specRatio->getParameter(0) ;
1100 RooProduct* dentmp = static_cast<RooProduct*>(specRatio->getParameter(1)) ;
1101
1102// cout << "numtmp = " << numtmp->ClassName() << "::" << numtmp->GetName() << endl ;
1103// cout << "dentmp = " << dentmp->ClassName() << "::" << dentmp->GetName() << endl ;
1104
1105// cout << "denominator components are " << dentmp->components() << endl ;
1106 for (auto* parg : static_range_cast<RooAbsReal*>(dentmp->components())) {
1107// cout << "now processing denominator component " << parg->ClassName() << "::" << parg->GetName() << endl ;
1108
1109 if (ratio && parg->dependsOn(*ratio)) {
1110// cout << "depends in value of ratio" << endl ;
1111
1112 // Make specialize ratio instance
1113 std::unique_ptr<RooAbsReal> specializedRatio{specializeRatio(*(RooFormulaVar*)ratio,iter->c_str())};
1114
1115// cout << "specRatio = " << endl ;
1116// specializedRatio->printComponentTree() ;
1117
1118 // Replace generic ratio with specialized ratio
1119 RooAbsArg *partCust(nullptr) ;
1120 if (parg->InheritsFrom(RooAddition::Class())) {
1121
1122
1123
1124 RooAddition* tmpadd = static_cast<RooAddition*>(parg) ;
1125
1126 RooCustomizer cust(*tmpadd->list1().first(),Form("blah_%s",iter->c_str())) ;
1127 cust.replaceArg(*ratio,*specializedRatio) ;
1128 partCust = cust.build() ;
1129
1130 } else {
1131 RooCustomizer cust(*parg,Form("blah_%s",iter->c_str())) ;
1132 cust.replaceArg(*ratio,*specializedRatio) ;
1133 partCust = cust.build() ;
1134 }
1135
1136 // Print customized denominator
1137// cout << "customized function = " << endl ;
1138// partCust->printComponentTree() ;
1139
1140 std::unique_ptr<RooAbsReal> specializedPartCust{specializeIntegral(*static_cast<RooAbsReal*>(partCust),iter->c_str())};
1141
1142 // Finally divide again by ratio
1143 string name = Form("%s_divided_by_ratio",specializedPartCust->GetName()) ;
1144 auto specIntFinal = std::make_unique<RooFormulaVar>(name.c_str(),"@0/@1",RooArgList(*specializedPartCust,*specializedRatio)) ;
1145 specIntFinal->addOwnedComponents(std::move(specializedPartCust));
1146 specIntFinal->addOwnedComponents(std::move(specializedRatio));
1147
1148 denListList[*iter].addOwned(std::move(specIntFinal));
1149 } else {
1150
1151// cout << "does NOT depend on value of ratio" << endl ;
1152// parg->Print("t") ;
1153
1154 denListList[*iter].addOwned(specializeIntegral(*parg,iter->c_str()));
1155
1156 }
1157 }
1158// cout << "end iteration over denominator components" << endl ;
1159 } else {
1160
1161 if (ratio) {
1162
1163 std::unique_ptr<RooAbsReal> specRatio{specializeRatio(*(RooFormulaVar*)ratio,iter->c_str())};
1164
1165 // If integral is 'Int r(y)*g(y) dy ' then divide a posteriori by r(y)
1166// cout << "have ratio, orig den = " << den->GetName() << endl ;
1167
1168 RooArgSet tmp(origNumTerm) ;
1169 tmp.add(*specRatio) ;
1170 const string pname = makeRGPPName("PROD",tmp,RooArgSet(),RooArgSet(),nullptr) ;
1171 auto specDenProd = std::make_unique<RooProduct>(pname.c_str(),pname.c_str(),tmp) ;
1172 std::unique_ptr<RooAbsReal> specInt;
1173
1175 specInt = std::unique_ptr<RooAbsReal>{specDenProd->createIntegral((static_cast<RooRealIntegral*>(den))->intVars(),iter->c_str())};
1176 specInt->addOwnedComponents(std::move(specDenProd));
1177 } else if (den->InheritsFrom(RooAddition::Class())) {
1178 RooAddition* orig = static_cast<RooAddition*>(den) ;
1179 RooRealIntegral* origInt = static_cast<RooRealIntegral*>(orig->list1().first()) ;
1180 specInt = std::unique_ptr<RooAbsReal>{specDenProd->createIntegral(origInt->intVars(),iter->c_str())};
1181 specInt->addOwnedComponents(std::move(specDenProd));
1182 } else {
1183 throw string("this should not happen") ;
1184 }
1185
1186 //RooAbsReal* specInt = specializeIntegral(*den,iter->c_str()) ;
1187 string name = Form("%s_divided_by_ratio",specInt->GetName()) ;
1188 auto specIntFinal = std::make_unique<RooFormulaVar>(name.c_str(),"@0/@1",RooArgList(*specInt,*specRatio)) ;
1189 specIntFinal->addOwnedComponents(std::move(specInt));
1190 specIntFinal->addOwnedComponents(std::move(specRatio));
1191 denListList[*iter].addOwned(std::move(specIntFinal));
1192 } else {
1193 denListList[*iter].addOwned(specializeIntegral(*den,iter->c_str()));
1194 }
1195
1196 }
1197 }
1198
1199 }
1200
1201 // Do not rearrage terms if numerator and denominator are effectively empty
1202 if (nomList.empty()) {
1203 return ;
1204 }
1205
1206 string name = Form("%s_numerator",GetName()) ;
1207 // WVE FIX THIS (2)
1208
1209 std::unique_ptr<RooAbsReal> numerator = std::make_unique<RooProduct>(name.c_str(),name.c_str(),nomList) ;
1210
1211 RooArgSet products ;
1212// cout << "nomList = " << nomList << endl ;
1213 for (map<string,RooArgSet>::iterator iter = denListList.begin() ; iter != denListList.end() ; ++iter) {
1214// cout << "denList[" << iter->first << "] = " << iter->second << endl ;
1215 name = Form("%s_denominator_comp_%s",GetName(),iter->first.c_str()) ;
1216 // WVE FIX THIS (2)
1217 RooProduct* prod_comp = new RooProduct(name.c_str(),name.c_str(),iter->second) ;
1218 prod_comp->addOwnedComponents(std::move(iter->second));
1219 products.add(*prod_comp) ;
1220 }
1221 name = Form("%s_denominator_sum",GetName()) ;
1222 RooAbsReal* norm = new RooAddition(name.c_str(),name.c_str(),products) ;
1223 norm->addOwnedComponents(products) ;
1224
1225 if (!specIntDeps.empty()) {
1226 // Apply posterior integration required for SPECINT case
1227
1228 string namesr = Form("SPEC_RATIO(%s,%s)",numerator->GetName(),norm->GetName()) ;
1229 RooFormulaVar* ndr = new RooFormulaVar(namesr.c_str(),"@0/@1",RooArgList(*numerator,*norm)) ;
1230 ndr->addOwnedComponents(std::move(numerator));
1231
1232 // Integral of ratio
1233 numerator = std::unique_ptr<RooAbsReal>{ndr->createIntegral(specIntDeps,specIntRange.c_str())};
1234
1235 norm = static_cast<RooAbsReal*>(RooFit::RooConst(1).Clone()) ;
1236 }
1237
1238
1239// cout << "numerator" << endl ;
1240// numerator->printComponentTree("",0,5) ;
1241// cout << "denominator" << endl ;
1242// norm->printComponentTree("",0,5) ;
1243
1244
1245 // WVE DEBUG
1246 //RooMsgService::instance().debugWorkspace()->import(RooArgSet(*numerator,*norm)) ;
1247
1248 cache._rearrangedNum = std::move(numerator);
1249 cache._rearrangedDen.reset(norm);
1250 cache._isRearranged = true ;
1251
1252}
1253
1254
1255////////////////////////////////////////////////////////////////////////////////
1256
1257std::unique_ptr<RooAbsReal> RooProdPdf::specializeRatio(RooFormulaVar& input, const char* targetRangeName) const
1258{
1259 RooRealIntegral* numint = static_cast<RooRealIntegral*>(input.getParameter(0)) ;
1260 RooRealIntegral* denint = static_cast<RooRealIntegral*>(input.getParameter(1)) ;
1261
1262 std::unique_ptr<RooAbsReal> numint_spec{specializeIntegral(*numint,targetRangeName)};
1263
1264 std::unique_ptr<RooAbsReal> ret = std::make_unique<RooFormulaVar>(Form("ratio(%s,%s)",numint_spec->GetName(),denint->GetName()),"@0/@1",RooArgList(*numint_spec,*denint)) ;
1265 ret->addOwnedComponents(std::move(numint_spec));
1266
1267 return ret;
1268}
1269
1270
1271
1272////////////////////////////////////////////////////////////////////////////////
1273
1274std::unique_ptr<RooAbsReal> RooProdPdf::specializeIntegral(RooAbsReal& input, const char* targetRangeName) const
1275{
1276 if (input.InheritsFrom(RooRealIntegral::Class())) {
1277
1278 // If input is integral, recreate integral but override integration range to be targetRangeName
1279 RooRealIntegral* orig = static_cast<RooRealIntegral*>(&input) ;
1280// cout << "creating integral: integrand = " << orig->integrand().GetName() << " vars = " << orig->intVars() << " range = " << targetRangeName << endl ;
1281 return std::unique_ptr<RooAbsReal>{orig->integrand().createIntegral(orig->intVars(),targetRangeName)};
1282
1283 } else if (input.InheritsFrom(RooAddition::Class())) {
1284
1285 // If input is sum of integrals, recreate integral from first component of set, but override integration range to be targetRangeName
1286 RooAddition* orig = static_cast<RooAddition*>(&input) ;
1287 RooRealIntegral* origInt = static_cast<RooRealIntegral*>(orig->list1().first()) ;
1288// cout << "creating integral from addition: integrand = " << origInt->integrand().GetName() << " vars = " << origInt->intVars() << " range = " << targetRangeName << endl ;
1289 return std::unique_ptr<RooAbsReal>{origInt->integrand().createIntegral(origInt->intVars(),targetRangeName)};
1290 }
1291
1292 std::stringstream errMsg;
1293 errMsg << "specializeIntegral: unknown input type " << input.ClassName() << "::" << input.GetName();
1294 throw std::runtime_error(errMsg.str());
1295}
1296
1297
1298////////////////////////////////////////////////////////////////////////////////
1299/// Group product into terms that can be calculated independently
1300
1301void RooProdPdf::groupProductTerms(std::list<std::vector<RooArgSet*>>& groupedTerms, RooArgSet& outerIntDeps,
1302 const RooLinkedList& terms, const RooLinkedList& norms,
1303 const RooLinkedList& imps, const RooLinkedList& ints, const RooLinkedList& /*cross*/) const
1304{
1305 // Start out with each term in its own group
1306 for(auto * term : static_range_cast<RooArgSet*>(terms)) {
1307 groupedTerms.emplace_back();
1308 groupedTerms.back().emplace_back(term) ;
1309 }
1310
1311 // Make list of imported dependents that occur in any term
1312 RooArgSet allImpDeps ;
1313 for(auto * impDeps : static_range_cast<RooArgSet*>(imps)) {
1314 allImpDeps.add(*impDeps,false) ;
1315 }
1316
1317 // Make list of integrated dependents that occur in any term
1318 RooArgSet allIntDeps ;
1319 for(auto * intDeps : static_range_cast<RooArgSet*>(ints)) {
1320 allIntDeps.add(*intDeps,false) ;
1321 }
1322
1323 outerIntDeps.removeAll() ;
1324 outerIntDeps.add(*std::unique_ptr<RooArgSet>{allIntDeps.selectCommon(allImpDeps)});
1325
1326 // Now iteratively merge groups that should be (partially) integrated together
1327 for(RooAbsArg * outerIntDep : outerIntDeps) {
1328
1329 // Collect groups that feature this dependent
1330 std::vector<RooArgSet*>* newGroup = nullptr ;
1331
1332 // Loop over groups
1333 bool needMerge = false ;
1334 auto group = groupedTerms.begin();
1335 auto nGroups = groupedTerms.size();
1336 for (size_t iGroup = 0; iGroup < nGroups; ++iGroup) {
1337
1338 // See if any term in this group depends in any ay on outerDepInt
1339 for (auto const& term2 : *group) {
1340
1341 Int_t termIdx = terms.IndexOf(term2) ;
1342 RooArgSet* termNormDeps = static_cast<RooArgSet*>(norms.At(termIdx)) ;
1343 RooArgSet* termIntDeps = static_cast<RooArgSet*>(ints.At(termIdx)) ;
1344 RooArgSet* termImpDeps = static_cast<RooArgSet*>(imps.At(termIdx)) ;
1345
1346 if (termNormDeps->contains(*outerIntDep) ||
1347 termIntDeps->contains(*outerIntDep) ||
1348 termImpDeps->contains(*outerIntDep)) {
1349 needMerge = true ;
1350 }
1351
1352 }
1353
1354 if (needMerge) {
1355 // Create composite group if not yet existing
1356 if (newGroup==nullptr) {
1357 groupedTerms.emplace_back() ;
1358 newGroup = &groupedTerms.back() ;
1359 }
1360
1361 // Add terms of this group to new term
1362 for (auto& term2 : *group) {
1363 newGroup->emplace_back(term2) ;
1364 }
1365
1366 // Remove this non-owning group from list
1367 group = groupedTerms.erase(group);
1368 } else {
1369 ++group;
1370 }
1371 }
1372
1373 }
1374}
1375
1376
1377
1378////////////////////////////////////////////////////////////////////////////////
1379/// Calculate integrals of factorized product terms over observables iset while normalized
1380/// to observables in nset.
1381
1382std::vector<RooAbsReal*> RooProdPdf::processProductTerm(const RooArgSet* nset, const RooArgSet* iset, const char* isetRangeName,
1383 const RooArgSet* term,const RooArgSet& termNSet, const RooArgSet& termISet,
1384 bool& isOwned, bool forceWrap) const
1385{
1386 vector<RooAbsReal*> ret(3) ; ret[0] = nullptr ; ret[1] = nullptr ; ret[2] = nullptr ;
1387
1388 // CASE I: factorizing term: term is integrated over all normalizing observables
1389 // -----------------------------------------------------------------------------
1390 // Check if all observbales of this term are integrated. If so the term cancels
1391 if (!termNSet.empty() && termNSet.size()==termISet.size() && isetRangeName==nullptr) {
1392
1393
1394 //cout << "processProductTerm(" << GetName() << ") case I " << endl ;
1395
1396 // Term factorizes
1397 return ret ;
1398 }
1399
1400 // CASE II: Dropped terms: if term is entirely unnormalized, it should be dropped
1401 // ------------------------------------------------------------------------------
1402 if (nset && termNSet.empty()) {
1403
1404 //cout << "processProductTerm(" << GetName() << ") case II " << endl ;
1405
1406 // Drop terms that are not asked to be normalized
1407 return ret ;
1408 }
1409
1410 if (iset && !termISet.empty()) {
1411 if (term->size()==1) {
1412
1413 // CASE IIIa: Normalized and partially integrated single PDF term
1414 //---------------------------------------------------------------
1415
1416 RooAbsPdf* pdf = static_cast<RooAbsPdf*>(term->first()) ;
1417
1418 RooAbsReal* partInt = std::unique_ptr<RooAbsReal>{pdf->createIntegral(termISet,termNSet,isetRangeName)}.release();
1419 partInt->setOperMode(operMode()) ;
1420 partInt->setStringAttribute("PROD_TERM_TYPE","IIIa") ;
1421
1422 isOwned=true ;
1423
1424 //cout << "processProductTerm(" << GetName() << ") case IIIa func = " << partInt->GetName() << endl ;
1425
1426 ret[0] = partInt ;
1427
1428 // Split mode results
1429 ret[1] = std::unique_ptr<RooAbsReal>{pdf->createIntegral(termISet,isetRangeName)}.release();
1430 ret[2] = std::unique_ptr<RooAbsReal>{pdf->createIntegral(termNSet,normRange())}.release();
1431
1432 return ret ;
1433
1434
1435 } else {
1436
1437 // CASE IIIb: Normalized and partially integrated composite PDF term
1438 //---------------------------------------------------------------
1439
1440 // Use auxiliary class RooGenProdProj to calculate this term
1441 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1442 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName) ;
1443 partInt->setStringAttribute("PROD_TERM_TYPE","IIIb") ;
1444 partInt->setOperMode(operMode()) ;
1445
1446 //cout << "processProductTerm(" << GetName() << ") case IIIb func = " << partInt->GetName() << endl ;
1447
1448 isOwned=true ;
1449 ret[0] = partInt ;
1450
1451 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),nullptr) ;
1452
1453 // WVE FIX THIS
1454 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1455
1456 ret[1] = std::unique_ptr<RooAbsReal>{tmp_prod->createIntegral(termISet,isetRangeName)}.release();
1457 ret[2] = std::unique_ptr<RooAbsReal>{tmp_prod->createIntegral(termNSet,normRange())}.release();
1458
1459 return ret ;
1460 }
1461 }
1462
1463 // CASE IVa: Normalized non-integrated composite PDF term
1464 // -------------------------------------------------------
1465 if (nset && !nset->empty() && term->size()>1) {
1466 // Composite term needs normalized integration
1467
1468 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1469 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName,normRange()) ;
1471
1472 partInt->setStringAttribute("PROD_TERM_TYPE","IVa") ;
1473 partInt->setOperMode(operMode()) ;
1474
1475 //cout << "processProductTerm(" << GetName() << ") case IVa func = " << partInt->GetName() << endl ;
1476
1477 isOwned=true ;
1478 ret[0] = partInt ;
1479
1480 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),nullptr) ;
1481
1482 // WVE FIX THIS
1483 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1484
1485 ret[1] = std::unique_ptr<RooAbsReal>{tmp_prod->createIntegral(termISet,isetRangeName)}.release();
1486 ret[2] = std::unique_ptr<RooAbsReal>{tmp_prod->createIntegral(termNSet,normRange())}.release();
1487
1488 return ret ;
1489 }
1490
1491 // CASE IVb: Normalized, non-integrated single PDF term
1492 // -----------------------------------------------------
1493 for (auto* pdf : static_range_cast<RooAbsPdf*>(*term)) {
1494
1495 if (forceWrap) {
1496
1497 // Construct representative name of normalization wrapper
1498 TString name(pdf->GetName()) ;
1499 name.Append("_NORM[") ;
1500 bool first(true) ;
1501 for (auto const* arg : termNSet) {
1502 if (!first) {
1503 name.Append(",") ;
1504 } else {
1505 first=false ;
1506 }
1507 name.Append(arg->GetName()) ;
1508 }
1509 if (normRange()) {
1510 name.Append("|") ;
1511 name.Append(normRange()) ;
1512 }
1513 name.Append("]") ;
1514
1515 RooAbsReal* partInt = new RooRealIntegral(name.Data(),name.Data(),*pdf,RooArgSet(),&termNSet) ;
1516 partInt->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1517 isOwned=true ;
1518
1519 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << partInt->GetName() << endl ;
1520
1521 ret[0] = partInt ;
1522
1523 ret[1] = std::unique_ptr<RooAbsReal>{pdf->createIntegral(RooArgSet())}.release();
1524 ret[2] = std::unique_ptr<RooAbsReal>{pdf->createIntegral(termNSet,normRange())}.release();
1525
1526 return ret ;
1527
1528
1529 } else {
1530 isOwned=false ;
1531
1532 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << pdf->GetName() << endl ;
1533
1534
1535 pdf->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1536 ret[0] = pdf ;
1537
1538 ret[1] = std::unique_ptr<RooAbsReal>{pdf->createIntegral(RooArgSet())}.release();
1539 ret[2] = !termNSet.empty() ? std::unique_ptr<RooAbsReal>{pdf->createIntegral(termNSet,normRange())}.release()
1540 : (static_cast<RooAbsReal*>(RooFit::RooConst(1).clone("1")));
1541 return ret ;
1542 }
1543 }
1544
1545 coutE(Eval) << "RooProdPdf::processProductTerm(" << GetName() << ") unidentified term!!!" << endl ;
1546 return ret ;
1547}
1548
1549
1550
1551
1552////////////////////////////////////////////////////////////////////////////////
1553/// Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1554
1555std::string RooProdPdf::makeRGPPName(const char* pfx, const RooArgSet& term, const RooArgSet& iset,
1556 const RooArgSet& nset, const char* isetRangeName) const
1557{
1558 // Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1559
1560 std::ostringstream os(pfx);
1561 os << "[";
1562
1563 // Encode component names
1564 bool first(true) ;
1565 for (auto const* pdf : static_range_cast<RooAbsPdf*>(term)) {
1566 if (!first) os << "_X_";
1567 first = false;
1568 os << pdf->GetName();
1569 }
1570 os << "]" << integralNameSuffix(iset,&nset,isetRangeName,true);
1571
1572 return os.str();
1573}
1574
1575
1576
1577////////////////////////////////////////////////////////////////////////////////
1578/// Force RooRealIntegral to offer all observables for internal integration
1579
1581{
1582 return true ;
1583}
1584
1585
1586
1587////////////////////////////////////////////////////////////////////////////////
1588/// Determine which part (if any) of given integral can be performed analytically.
1589/// If any analytical integration is possible, return integration scenario code.
1590///
1591/// RooProdPdf implements two strategies in implementing analytical integrals
1592///
1593/// First, PDF components whose entire set of dependents are requested to be integrated
1594/// can be dropped from the product, as they will integrate out to 1 by construction
1595///
1596/// Second, RooProdPdf queries each remaining component PDF for its analytical integration
1597/// capability of the requested set ('allVars'). It finds the largest common set of variables
1598/// that can be integrated by all remaining components. If such a set exists, it reconfirms that
1599/// each component is capable of analytically integrating the common set, and combines the components
1600/// individual integration codes into a single integration code valid for RooProdPdf.
1601
1603 const RooArgSet* normSet, const char* rangeName) const
1604{
1605 if (_forceNumInt) return 0 ;
1606
1607 // Declare that we can analytically integrate all requested observables
1608 analVars.add(allVars) ;
1609
1610 // Retrieve (or create) the required partial integral list
1611 Int_t code = getPartIntList(normSet,&allVars,rangeName);
1612
1613 return code+1 ;
1614}
1615
1616
1617
1618
1619////////////////////////////////////////////////////////////////////////////////
1620/// Return analytical integral defined by given scenario code
1621
1622double RooProdPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
1623{
1624 // No integration scenario
1625 if (code==0) {
1626 return getVal(normSet) ;
1627 }
1628
1629
1630 // WVE needs adaptation for rangename feature
1631
1632 // Partial integration scenarios
1633 CacheElem* cache = static_cast<CacheElem*>(_cacheMgr.getObjByIndex(code-1)) ;
1634
1635 // If cache has been sterilized, revive this slot
1636 if (cache==nullptr) {
1637 std::unique_ptr<RooArgSet> vars{getParameters(RooArgSet())} ;
1638 RooArgSet nset = _cacheMgr.selectFromSet1(*vars, code-1) ;
1639 RooArgSet iset = _cacheMgr.selectFromSet2(*vars, code-1) ;
1640
1641 Int_t code2 = getPartIntList(&nset, &iset, rangeName) ;
1642
1643 // preceding call to getPartIntList guarantees non-null return
1644 // coverity[NULL_RETURNS]
1645 cache = static_cast<CacheElem*>(_cacheMgr.getObj(&nset,&iset,&code2,rangeName)) ;
1646 }
1647
1648 double val = calculate(*cache,true) ;
1649// cout << "RPP::aIWN(" << GetName() << ") ,code = " << code << ", value = " << val << endl ;
1650
1651 return val ;
1652}
1653
1654
1655
1656////////////////////////////////////////////////////////////////////////////////
1657/// If this product contains exactly one extendable p.d.f return the extension abilities of
1658/// that p.d.f, otherwise return CanNotBeExtended
1659
1661{
1662 return (_extendedIndex>=0) ? (static_cast<RooAbsPdf*>(_pdfList.at(_extendedIndex)))->extendMode() : CanNotBeExtended ;
1663}
1664
1665
1666
1667////////////////////////////////////////////////////////////////////////////////
1668/// Return the expected number of events associated with the extendable input PDF
1669/// in the product. If there is no extendable term, abort.
1670
1671double RooProdPdf::expectedEvents(const RooArgSet* nset) const
1672{
1673 if (_extendedIndex<0) {
1674 coutF(Generation) << "Requesting expected number of events from a RooProdPdf that does not contain an extended p.d.f" << endl ;
1675 throw std::logic_error(std::string("RooProdPdf ") + GetName() + " could not be extended.");
1676 }
1677
1678 return static_cast<RooAbsPdf*>(_pdfList.at(_extendedIndex))->expectedEvents(nset) ;
1679}
1680
1681std::unique_ptr<RooAbsReal> RooProdPdf::createExpectedEventsFunc(const RooArgSet* nset) const
1682{
1683 if (_extendedIndex<0) {
1684 coutF(Generation) << "Requesting expected number of events from a RooProdPdf that does not contain an extended p.d.f" << endl ;
1685 throw std::logic_error(std::string("RooProdPdf ") + GetName() + " could not be extended.");
1686 }
1687
1688 return static_cast<RooAbsPdf*>(_pdfList.at(_extendedIndex))->createExpectedEventsFunc(nset);
1689}
1690
1691
1692////////////////////////////////////////////////////////////////////////////////
1693/// Return generator context optimized for generating events from product p.d.f.s
1694
1696 const RooArgSet* auxProto, bool verbose) const
1697{
1698 if (_useDefaultGen) return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ;
1699 return new RooProdGenContext(*this,vars,prototype,auxProto,verbose) ;
1700}
1701
1702
1703
1704////////////////////////////////////////////////////////////////////////////////
1705/// Query internal generation capabilities of component p.d.f.s and aggregate capabilities
1706/// into master configuration passed to the generator context
1707
1708Int_t RooProdPdf::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
1709{
1710 if (!_useDefaultGen) return 0 ;
1711
1712 // Find the subset directVars that only depend on a single PDF in the product
1713 RooArgSet directSafe ;
1714 for (auto const* arg : directVars) {
1715 if (isDirectGenSafe(*arg)) directSafe.add(*arg) ;
1716 }
1717
1718
1719 // Now find direct integrator for relevant components ;
1720 std::vector<Int_t> code;
1721 code.reserve(64);
1722 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1723 RooArgSet pdfDirect ;
1724 Int_t pdfCode = pdf->getGenerator(directSafe,pdfDirect,staticInitOK);
1725 code.push_back(pdfCode);
1726 if (pdfCode != 0) {
1727 generateVars.add(pdfDirect) ;
1728 }
1729 }
1730
1731
1732 if (!generateVars.empty()) {
1733 Int_t masterCode = _genCode.store(code) ;
1734 return masterCode+1 ;
1735 } else {
1736 return 0 ;
1737 }
1738}
1739
1740
1741
1742////////////////////////////////////////////////////////////////////////////////
1743/// Forward one-time initialization call to component generation initialization
1744/// methods.
1745
1747{
1748 if (!_useDefaultGen) return ;
1749
1750 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1751 Int_t i(0) ;
1752 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1753 if (codeList[i]!=0) {
1754 pdf->initGenerator(codeList[i]) ;
1755 }
1756 i++ ;
1757 }
1758}
1759
1760
1761
1762////////////////////////////////////////////////////////////////////////////////
1763/// Generate a single event with configuration specified by 'code'
1764/// Defer internal generation to components as encoded in the _genCode
1765/// registry for given generator code.
1766
1768{
1769 if (!_useDefaultGen) return ;
1770
1771 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1772 Int_t i(0) ;
1773 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1774 if (codeList[i]!=0) {
1775 pdf->generateEvent(codeList[i]) ;
1776 }
1777 i++ ;
1778 }
1779}
1780
1781
1782
1783////////////////////////////////////////////////////////////////////////////////
1784/// Return RooAbsArg components contained in the cache
1785
1787{
1788 RooArgList ret ;
1789 ret.add(_partList) ;
1790 ret.add(_numList) ;
1791 ret.add(_denList) ;
1792 if (_rearrangedNum) ret.add(*_rearrangedNum) ;
1793 if (_rearrangedDen) ret.add(*_rearrangedDen) ;
1794 return ret ;
1795
1796}
1797
1798
1799
1800////////////////////////////////////////////////////////////////////////////////
1801/// Hook function to print cache contents in tree printing of RooProdPdf
1802
1803void RooProdPdf::CacheElem::printCompactTreeHook(ostream& os, const char* indent, Int_t curElem, Int_t maxElem)
1804{
1805 if (curElem==0) {
1806 os << indent << "RooProdPdf begin partial integral cache" << endl ;
1807 }
1808
1809 auto indent2 = std::string(indent) + "[" + std::to_string(curElem) + "]";
1810 for(auto const& arg : _partList) {
1811 arg->printCompactTree(os,indent2.c_str()) ;
1812 }
1813
1814 if (curElem==maxElem) {
1815 os << indent << "RooProdPdf end partial integral cache" << endl ;
1816 }
1817}
1818
1819
1820
1821////////////////////////////////////////////////////////////////////////////////
1822/// Forward determination of safety of internal generator code to
1823/// component p.d.f that would generate the given observable
1824
1826{
1827 // Only override base class behaviour if default generator method is enabled
1828 if (!_useDefaultGen) return RooAbsPdf::isDirectGenSafe(arg) ;
1829
1830 // Argument may appear in only one PDF component
1831 RooAbsPdf* thePdf(nullptr) ;
1832 for (auto* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
1833
1834 if (pdf->dependsOn(arg)) {
1835 // Found PDF depending on arg
1836
1837 // If multiple PDFs depend on arg directGen is not safe
1838 if (thePdf) return false ;
1839
1840 thePdf = pdf ;
1841 }
1842 }
1843 // Forward call to relevant component PDF
1844 return thePdf?(thePdf->isDirectGenSafe(arg)):false ;
1845}
1846
1847
1848
1849////////////////////////////////////////////////////////////////////////////////
1850/// Look up user specified normalization set for given input PDF component
1851
1853{
1854 Int_t idx = _pdfList.index(&pdf) ;
1855 if (idx<0) return nullptr;
1856 return _pdfNSetList[idx].get() ;
1857}
1858
1859
1860
1861/// Add some full PDFs to the factors of this RooProdPdf.
1863{
1864 size_t numExtended = (_extendedIndex==-1) ? 0 : 1;
1865
1866 for(auto arg : pdfs) {
1867 RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg);
1868 if (!pdf) {
1869 coutW(InputArguments) << "RooProdPdf::addPdfs(" << GetName() << ") list arg "
1870 << arg->GetName() << " is not a PDF, ignored" << endl ;
1871 continue;
1872 }
1873 if(pdf->canBeExtended()) {
1874 if (_extendedIndex == -1) {
1876 } else {
1877 numExtended++;
1878 }
1879 }
1880 _pdfList.add(*pdf);
1881 _pdfNSetList.emplace_back(std::make_unique<RooArgSet>("nset"));
1882 }
1883
1884 // Protect against multiple extended terms
1885 if (numExtended>1) {
1886 coutW(InputArguments) << "RooProdPdf::addPdfs(" << GetName()
1887 << ") WARNING: multiple components with extended terms detected,"
1888 << " product will not be extendible." << endl ;
1889 _extendedIndex = -1 ;
1890 }
1891
1892 // Reset cache
1893 _cacheMgr.reset() ;
1894
1895}
1896
1897/// Remove some PDFs from the factors of this RooProdPdf.
1899{
1900 // Remember what the extended PDF is
1901 RooAbsArg const* extPdf = _extendedIndex >= 0 ? &_pdfList[_extendedIndex] : nullptr;
1902
1903 // Actually remove the PDFs and associated nsets
1904 for(size_t i=0;i < _pdfList.size(); i++) {
1905 if(pdfs.contains(_pdfList[i])) {
1907 _pdfNSetList.erase(_pdfNSetList.begin()+i);
1908 i--;
1909 }
1910 }
1911
1912 // Since we may have removed PDFs from the list, the index of the extended
1913 // PDF in the list needs to be updated. The new index might also be -1 if the
1914 // extended PDF got removed.
1915 if(extPdf) {
1916 _extendedIndex = _pdfList.index(*extPdf);
1917 }
1918
1919 // Reset cache
1920 _cacheMgr.reset() ;
1921}
1922
1923
1924namespace {
1925
1926std::vector<TNamed const*> sortedNamePtrs(RooAbsCollection const& col)
1927{
1928 std::vector<TNamed const*> ptrs;
1929 ptrs.reserve(col.size());
1930 for(RooAbsArg* arg : col) {
1931 ptrs.push_back(arg->namePtr());
1932 }
1933 std::sort(ptrs.begin(), ptrs.end());
1934 return ptrs;
1935}
1936
1937bool sortedNamePtrsOverlap(std::vector<TNamed const*> const& ptrsA, std::vector<TNamed const*> const& ptrsB)
1938{
1939 auto pA = ptrsA.begin();
1940 auto pB = ptrsB.begin();
1941 while (pA != ptrsA.end() && pB != ptrsB.end()) {
1942 if (*pA < *pB) {
1943 ++pA;
1944 } else if (*pB < *pA) {
1945 ++pB;
1946 } else {
1947 return true;
1948 }
1949 }
1950 return false;
1951}
1952
1953} // namespace
1954
1955
1956////////////////////////////////////////////////////////////////////////////////
1957/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
1958/// The observables set is required to distinguish unambiguously p.d.f in terms
1959/// of observables and parameters, which are not constraints, and p.d.fs in terms
1960/// of parameters only, which can serve as constraints p.d.f.s
1961
1962RooArgSet* RooProdPdf::getConstraints(const RooArgSet& observables, RooArgSet& constrainedParams,
1963 bool stripDisconnected, bool removeConstraintsFromPdf) const
1964{
1965 RooArgSet constraints ;
1966 RooArgSet pdfParams;
1967 RooArgSet conParams;
1968
1969 // For the optimized implementation of checking if two collections overlap by name.
1970 auto observablesNamePtrs = sortedNamePtrs(observables);
1971 auto constrainedParamsNamePtrs = sortedNamePtrs(constrainedParams);
1972
1973 // Loop over PDF components
1974 for (std::size_t iPdf = 0; iPdf < _pdfList.size(); ++iPdf) {
1975 auto * pdf = static_cast<RooAbsPdf*>(&_pdfList[iPdf]);
1976
1977 RooArgSet tmp;
1978 pdf->getParameters(nullptr, tmp);
1979
1980 // A constraint term is a p.d.f that doesn't contribute to the
1981 // expectedEvents() and does not depend on any of the listed observables
1982 // but does depends on any of the parameters that should be constrained
1983 bool isConstraint = false;
1984
1985 if(static_cast<int>(iPdf) != _extendedIndex) {
1986 auto tmpNamePtrs = sortedNamePtrs(tmp);
1987 // Before, there were calls to `pdf->dependsOn()` here, but they were very
1988 // expensive for large computation graphs! Given that we have to traverse
1989 // the computation graph with a call to `pdf->getParameters()` anyway, we
1990 // can just check if the set of all variables operlaps with the observables
1991 // or constraind parameters.
1992 //
1993 // We are using an optimized implementation of overlap checking. Because
1994 // the overlap is checked by name, we can check overlap of the
1995 // corresponding name pointers. The optimization can't be in
1996 // RooAbsCollection itself, because it is crucial that the memory for the
1997 // non-tmp name pointers is not reallocated for each pdf.
1998 isConstraint = !sortedNamePtrsOverlap(tmpNamePtrs, observablesNamePtrs) &&
1999 sortedNamePtrsOverlap(tmpNamePtrs, constrainedParamsNamePtrs);
2000 }
2001 if (isConstraint) {
2002 constraints.add(*pdf) ;
2003 conParams.add(tmp,true) ;
2004 } else {
2005 // We only want to add parameter, not observables. Since a call like
2006 // `pdf->getParameters(&observables)` would be expensive, we take the set
2007 // of all variables and remove the ovservables, which is much cheaper. In
2008 // a call to `pdf->getParameters(&observables)`, the observables are
2009 // matched by name, so we have to pass the `matchByNameOnly` here.
2010 tmp.remove(observables, /*silent=*/false, /*matchByNameOnly=*/true);
2011 pdfParams.add(tmp,true) ;
2012 }
2013 }
2014
2015 // Remove the constraints now from the PDF if the caller requested it
2016 if(removeConstraintsFromPdf) {
2017 const_cast<RooProdPdf*>(this)->removePdfs(constraints);
2018 }
2019
2020 // Strip any constraints that are completely decoupled from the other product terms
2021 RooArgSet* finalConstraints = new RooArgSet("constraints") ;
2022 for(auto * pdf : static_range_cast<RooAbsPdf*>(constraints)) {
2023 if (pdf->dependsOnValue(pdfParams) || !stripDisconnected) {
2024 finalConstraints->add(*pdf) ;
2025 } else {
2026 coutI(Minimization) << "RooProdPdf::getConstraints(" << GetName() << ") omitting term " << pdf->GetName()
2027 << " as constraint term as it does not share any parameters with the other pdfs in product. "
2028 << "To force inclusion in likelihood, add an explicit Constrain() argument for the target parameter" << endl ;
2029 }
2030 }
2031
2032 // Now remove from constrainedParams all parameters that occur exclusively in constraint term and not in regular pdf term
2033
2034 RooArgSet cexl;
2035 conParams.selectCommon(constrainedParams, cexl);
2036 cexl.remove(pdfParams,true,true) ;
2037 constrainedParams.remove(cexl,true,true) ;
2038
2039 return finalConstraints ;
2040}
2041
2042
2043
2044
2045////////////////////////////////////////////////////////////////////////////////
2046/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams.
2047/// The observables set is required to distinguish unambiguously p.d.f in terms
2048/// of observables and parameters, which are not constraints, and p.d.fs in terms
2049/// of parameters only, which can serve as constraints p.d.f.s
2050
2052{
2053 RooArgSet* connectedPars = new RooArgSet("connectedPars") ;
2054 for (std::size_t iPdf = 0; iPdf < _pdfList.size(); ++iPdf) {
2055 auto * pdf = static_cast<RooAbsPdf*>(&_pdfList[iPdf]);
2056 // Check if term is relevant, either because it provides a propablity
2057 // density in the observables or because it is used for the expected
2058 // events.
2059 if (static_cast<int>(iPdf) == _extendedIndex || pdf->dependsOn(observables)) {
2060 RooArgSet tmp;
2061 pdf->getParameters(&observables, tmp);
2062 connectedPars->add(tmp) ;
2063 }
2064 }
2065 return connectedPars ;
2066}
2067
2068
2069
2070
2071////////////////////////////////////////////////////////////////////////////////
2072
2073void RooProdPdf::getParametersHook(const RooArgSet* nset, RooArgSet* params, bool stripDisconnected) const
2074{
2075 if (!stripDisconnected) return ;
2076 if (!nset || nset->empty()) return ;
2077
2078 // Get/create appropriate term list for this normalization set
2079 Int_t code = getPartIntList(nset, nullptr);
2080 RooArgList & plist = static_cast<CacheElem*>(_cacheMgr.getObj(nset, &code))->_partList;
2081
2082 // Strip any terms from params that do not depend on any term
2083 RooArgSet tostrip ;
2084 for (auto param : *params) {
2085 bool anyDep(false) ;
2086 for (auto term : plist) {
2087 if (term->dependsOnValue(*param)) {
2088 anyDep=true ;
2089 }
2090 }
2091 if (!anyDep) {
2092 tostrip.add(*param) ;
2093 }
2094 }
2095
2096 if (!tostrip.empty()) {
2097 params->remove(tostrip,true,true);
2098 }
2099
2100}
2101
2102
2103
2104////////////////////////////////////////////////////////////////////////////////
2105/// Interface function used by test statistics to freeze choice of range
2106/// for interpretation of conditional product terms
2107
2108void RooProdPdf::selectNormalizationRange(const char* rangeName, bool force)
2109{
2110 if (!force && _refRangeName) {
2111 return ;
2112 }
2113
2114 fixRefRange(rangeName) ;
2115}
2116
2117
2118
2119
2120////////////////////////////////////////////////////////////////////////////////
2121
2122void RooProdPdf::fixRefRange(const char* rangeName)
2123{
2124 _refRangeName = const_cast<TNamed*>(RooNameReg::ptr(rangeName));
2125}
2126
2127
2128
2129////////////////////////////////////////////////////////////////////////////////
2130/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2131
2132std::list<double>* RooProdPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
2133{
2134 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2135 if (std::list<double>* hint = pdf->plotSamplingHint(obs,xlo,xhi)) {
2136 return hint ;
2137 }
2138 }
2139
2140 return nullptr;
2141}
2142
2143
2144
2145////////////////////////////////////////////////////////////////////////////////
2146/// If all components that depend on obs are binned that so is the product
2147
2149{
2150 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2151 if (pdf->dependsOn(obs) && !pdf->isBinnedDistribution(obs)) {
2152 return false ;
2153 }
2154 }
2155
2156 return true ;
2157}
2158
2159
2160
2161
2162
2163
2164////////////////////////////////////////////////////////////////////////////////
2165/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2166
2167std::list<double>* RooProdPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
2168{
2169 for (auto const* pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
2170 if (std::list<double>* hint = pdf->binBoundaries(obs,xlo,xhi)) {
2171 return hint ;
2172 }
2173 }
2174
2175 return nullptr;
2176}
2177
2178
2179////////////////////////////////////////////////////////////////////////////////
2180/// Label OK'ed components of a RooProdPdf with cache-and-track, _and_ label all RooProdPdf
2181/// descendants with extra information about (conditional) normalization, needed to be able
2182/// to Cache-And-Track them outside the RooProdPdf context.
2183
2185{
2186 for (const auto parg : _pdfList) {
2187
2188 if (parg->canNodeBeCached()==Always) {
2189 trackNodes.add(*parg) ;
2190// cout << "tracking node RooProdPdf component " << parg << " " << parg->ClassName() << "::" << parg->GetName() << endl ;
2191
2192 // Additional processing to fix normalization sets in case product defines conditional observables
2193 if (RooArgSet* pdf_nset = findPdfNSet(static_cast<RooAbsPdf&>(*parg))) {
2194 // Check if conditional normalization is specified
2196 if (string("nset")==pdf_nset->GetName() && !pdf_nset->empty()) {
2197 parg->setStringAttribute("CATNormSet",getColonSeparatedNameString(*pdf_nset).c_str()) ;
2198 }
2199 if (string("cset")==pdf_nset->GetName()) {
2200 parg->setStringAttribute("CATCondSet",getColonSeparatedNameString(*pdf_nset).c_str()) ;
2201 }
2202 } else {
2203 coutW(Optimization) << "RooProdPdf::setCacheAndTrackHints(" << GetName() << ") WARNING product pdf does not specify a normalization set for component " << parg->GetName() << endl ;
2204 }
2205 }
2206 }
2207}
2208
2209
2210
2211////////////////////////////////////////////////////////////////////////////////
2212/// Customized printing of arguments of a RooProdPdf to more intuitively reflect the contents of the
2213/// product operator construction
2214
2215void RooProdPdf::printMetaArgs(ostream& os) const
2216{
2217 for (std::size_t i=0 ; i<_pdfList.size() ; i++) {
2218 if (i>0) os << " * " ;
2219 RooArgSet* ncset = _pdfNSetList[i].get() ;
2220 os << _pdfList.at(i)->GetName() ;
2221 if (!ncset->empty()) {
2222 if (string("nset")==ncset->GetName()) {
2223 os << *ncset ;
2224 } else {
2225 os << "|" ;
2226 bool first(true) ;
2227 for (auto const* arg : *ncset) {
2228 if (!first) {
2229 os << "," ;
2230 } else {
2231 first = false ;
2232 }
2233 os << arg->GetName() ;
2234 }
2235 }
2236 }
2237 }
2238 os << " " ;
2239}
2240
2241
2242
2243////////////////////////////////////////////////////////////////////////////////
2244/// Implement support for node removal
2245
2246bool RooProdPdf::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
2247{
2248 if (nameChange && _pdfList.find("REMOVAL_DUMMY")) {
2249
2250 cxcoutD(LinkStateMgmt) << "RooProdPdf::redirectServersHook(" << GetName() << "): removing REMOVAL_DUMMY" << endl ;
2251
2252 // Remove node from _pdfList proxy and remove corresponding entry from normset list
2253 RooAbsArg* pdfDel = _pdfList.find("REMOVAL_DUMMY") ;
2254
2255 _pdfNSetList.erase(_pdfNSetList.begin() + _pdfList.index("REMOVAL_DUMMY")) ;
2256 _pdfList.remove(*pdfDel) ;
2257
2258 // Clear caches
2259 _cacheMgr.reset() ;
2260 }
2261
2262 // If the replaced server is an observable that is used in any of the
2263 // normalization sets for conditional fits, replace the element in the
2264 // normalization set too.
2265 for(std::unique_ptr<RooArgSet> const& normSet : _pdfNSetList) {
2266 for(RooAbsArg * arg : *normSet) {
2267 if(RooAbsArg * newArg = arg->findNewServer(newServerList, nameChange)) {
2268 // Since normSet is owning, the original arg is now deleted.
2269 normSet->replace(arg, std::unique_ptr<RooAbsArg>{newArg->cloneTree()});
2270 }
2271 }
2272 }
2273
2274 return RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
2275}
2276
2277void RooProdPdf::CacheElem::writeToStream(std::ostream& os) const {
2278 using namespace RooHelpers;
2279 os << "_partList\n";
2280 os << getColonSeparatedNameString(_partList) << "\n";
2281 os << "_numList\n";
2282 os << getColonSeparatedNameString(_numList) << "\n";
2283 os << "_denList\n";
2284 os << getColonSeparatedNameString(_denList) << "\n";
2285 os << "_ownedList\n";
2286 os << getColonSeparatedNameString(_ownedList) << "\n";
2287 os << "_normList\n";
2288 for(auto const& set : _normList) {
2289 os << getColonSeparatedNameString(*set) << "\n";
2290 }
2291 os << "_isRearranged" << "\n";
2292 os << _isRearranged << "\n";
2293 os << "_rearrangedNum" << "\n";
2294 if(_rearrangedNum) {
2295 os << getColonSeparatedNameString(*_rearrangedNum) << "\n";
2296 } else {
2297 os << "nullptr" << "\n";
2298 }
2299 os << "_rearrangedDen" << "\n";
2300 if(_rearrangedDen) {
2301 os << getColonSeparatedNameString(*_rearrangedDen) << "\n";
2302 } else {
2303 os << "nullptr" << "\n";
2304 }
2305}
2306
2307std::unique_ptr<RooArgSet> RooProdPdf::fillNormSetForServer(RooArgSet const &normSet, RooAbsArg const &server) const
2308{
2309 if (normSet.empty())
2310 return nullptr;
2311 auto *pdfNset = findPdfNSet(static_cast<RooAbsPdf const &>(server));
2312 if (pdfNset && !pdfNset->empty()) {
2313 std::unique_ptr<RooArgSet> out;
2314 if (0 == strcmp("cset", pdfNset->GetName())) {
2315 // If the name of the normalization set is "cset", it doesn't contain the
2316 // normalization set but the conditional observables that should *not* be
2317 // normalized over.
2318 out = std::make_unique<RooArgSet>(normSet);
2319 RooArgSet common;
2320 out->selectCommon(*pdfNset, common);
2321 out->remove(common);
2322 } else {
2323 out = std::make_unique<RooArgSet>(*pdfNset);
2324 }
2325 // prefix also the arguments in the normSets if they have not already been
2326 if (auto prefix = getStringAttribute("__prefix__")) {
2327 for (RooAbsArg *arg : *out) {
2328 if (!arg->getStringAttribute("__prefix__")) {
2329 arg->SetName((std::string(prefix) + arg->GetName()).c_str());
2330 arg->setStringAttribute("__prefix__", prefix);
2331 }
2332 }
2333 }
2334 return out;
2335 } else {
2336 return nullptr;
2337 }
2338}
2339
2340/// A RooProdPdf with a fixed normalization set can be replaced by this class.
2341/// Its purpose is to provide the right client-server interface for the
2342/// evaluation of RooProdPdf cache elements that were created for a given
2343/// normalization set.
2345public:
2346 RooFixedProdPdf(std::unique_ptr<RooProdPdf> &&prodPdf, RooArgSet const &normSet)
2347 : RooAbsPdf(prodPdf->GetName(), prodPdf->GetTitle()), _normSet{normSet},
2348 _servers("!servers", "List of servers", this), _prodPdf{std::move(prodPdf)}
2349 {
2350 initialize();
2351 }
2352 RooFixedProdPdf(const RooFixedProdPdf &other, const char *name = nullptr)
2353 : RooAbsPdf(other, name), _normSet{other._normSet},
2354 _servers("!servers", "List of servers", this), _prodPdf{static_cast<RooProdPdf *>(other._prodPdf->Clone())}
2355 {
2356 initialize();
2357 }
2358 TObject *clone(const char *newname) const override { return new RooFixedProdPdf(*this, newname); }
2359
2360 bool selfNormalized() const override { return true; }
2361
2362 inline bool canComputeBatchWithCuda() const override { return true; }
2363
2364 void doEval(RooFit::EvalContext &ctx) const override
2365 {
2366 _prodPdf->doEvalImpl(this, *_cache, ctx);
2367 }
2368
2370 {
2371 if (_cache->_isRearranged) {
2372 ctx.addResult(this, ctx.buildCall("RooFit::Detail::MathFuncs::ratio", *_cache->_rearrangedNum, *_cache->_rearrangedDen));
2373 } else {
2374 ctx.addResult(this, ctx.buildCall("RooFit::Detail::MathFuncs::product", _cache->_partList, _cache->_partList.size()));
2375 }
2376 }
2377
2378 ExtendMode extendMode() const override { return _prodPdf->extendMode(); }
2379 double expectedEvents(const RooArgSet * /*nset*/) const override { return _prodPdf->expectedEvents(&_normSet); }
2380 std::unique_ptr<RooAbsReal> createExpectedEventsFunc(const RooArgSet * /*nset*/) const override
2381 {
2382 return _prodPdf->createExpectedEventsFunc(&_normSet);
2383 }
2384
2385 // Analytical Integration handling
2386 bool forceAnalyticalInt(const RooAbsArg &dep) const override { return _prodPdf->forceAnalyticalInt(dep); }
2387 Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet,
2388 const char *rangeName = nullptr) const override
2389 {
2390 return _prodPdf->getAnalyticalIntegralWN(allVars, analVars, normSet, rangeName);
2391 }
2392 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &numVars, const char *rangeName = nullptr) const override
2393 {
2394 return _prodPdf->getAnalyticalIntegral(allVars, numVars, rangeName);
2395 }
2396 double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName) const override
2397 {
2398 return _prodPdf->analyticalIntegralWN(code, normSet, rangeName);
2399 }
2400 double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override
2401 {
2402 return _prodPdf->analyticalIntegral(code, rangeName);
2403 }
2404
2405private:
2407 {
2408 _cache = _prodPdf->createCacheElem(&_normSet, nullptr);
2409 auto &cache = *_cache;
2410
2411 // The actual servers for a given normalization set depend on whether the
2412 // cache is rearranged or not. See RooProdPdf::calculateBatch to see
2413 // which args in the cache are used directly.
2414 if (cache._isRearranged) {
2415 _servers.add(*cache._rearrangedNum);
2416 _servers.add(*cache._rearrangedDen);
2417 } else {
2418 for (std::size_t i = 0; i < cache._partList.size(); ++i) {
2419 _servers.add(cache._partList[i]);
2420 }
2421 }
2422 }
2423
2424 double evaluate() const override { return _prodPdf->calculate(*_cache); }
2425
2427 std::unique_ptr<RooProdPdf::CacheElem> _cache;
2429 std::unique_ptr<RooProdPdf> _prodPdf;
2430};
2431
2432std::unique_ptr<RooAbsArg>
2434{
2435 if (ctx.likelihoodMode()) {
2436 auto binnedInfo = RooHelpers::getBinnedL(*this);
2437 if (binnedInfo.binnedPdf && binnedInfo.binnedPdf != this) {
2438 return binnedInfo.binnedPdf->compileForNormSet(normSet, ctx);
2439 }
2440 }
2441
2442 std::unique_ptr<RooProdPdf> prodPdfClone{static_cast<RooProdPdf *>(this->Clone())};
2443 ctx.markAsCompiled(*prodPdfClone);
2444
2445 for (const auto server : prodPdfClone->servers()) {
2446 auto nsetForServer = fillNormSetForServer(normSet, *server);
2447 RooArgSet const &nset = nsetForServer ? *nsetForServer : normSet;
2448
2449 RooArgSet depList;
2450 server->getObservables(&nset, depList);
2451
2452 ctx.compileServer(*server, *prodPdfClone, depList);
2453 }
2454
2455 auto fixedProdPdf = std::make_unique<RooFixedProdPdf>(std::move(prodPdfClone), normSet);
2456 ctx.markAsCompiled(*fixedProdPdf);
2457
2458 return fixedProdPdf;
2459}
bool initialize()
#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:377
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:2489
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...
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooExpensiveObjectCache & expensiveObjectCache() const
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:501
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:631
RooFit::OwningPtr< RooArgSet > getComponents() const
Create a RooArgSet with all components (branch nodes) of the expression tree headed by this object.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:89
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:482
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.
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
void clear()
Clear contents. If the collection is owning, it will also delete the contents.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
TString _normRange
Normalization range.
Definition RooAbsPdf.h:343
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:321
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:219
@ CanNotBeExtended
Definition RooAbsPdf.h:213
const char * normRange() const
Definition RooAbsPdf.h:251
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.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:546
TString integralNameSuffix(const RooArgSet &iset, const RooArgSet *nset=nullptr, const char *rangeName=nullptr, bool omitEmpty=false) const
Construct string with unique suffix name to give to integral object that encodes integrated observabl...
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:191
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:186
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 containing the objects that are both in the cached set 1 with a given index and an i...
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 2 with a given index and an i...
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.
Named container for two doubles, two integers two object points and three string pointers that can be...
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:29
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 occurrence 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...
Container class to hold unbinned data.
Definition RooDataSet.h:57
A class to maintain the context for squashing of RooFit models into code.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
void markAsCompiled(RooAbsArg &arg) const
void compileServer(RooAbsArg &server, RooAbsArg &arg, RooArgSet const &normSet)
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
A RooProdPdf with a fixed normalization set can be replaced by this class.
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
std::unique_ptr< RooProdPdf::CacheElem > _cache
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 ...
std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *) const override
Returns an object that represents the expected number of events for a given normalization set,...
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.
void doEval(RooFit::EvalContext &ctx) const override
Base function for computing multiple values of a RooAbsReal.
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...
General form of projected integral of product of PDFs, utility class for RooProdPdf.
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
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:39
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
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
std::unique_ptr< RooAbsReal > specializeRatio(RooFormulaVar &input, const char *targetRangeName) const
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
std::unique_ptr< RooAbsReal > specializeIntegral(RooAbsReal &orig, const char *targetRangeName) const
bool forceAnalyticalInt(const RooAbsArg &dep) const override
Force RooRealIntegral to offer all observables for internal integration.
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.
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...
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
std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const override
Returns an object that represents the expected number of events for a given normalization set,...
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 doEvalImpl(RooAbsArg const *caller, const RooProdPdf::CacheElem &cache, RooFit::EvalContext &) const
Evaluate product of PDFs in batch mode.
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
Represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
static TClass * Class()
RooArgList components()
Definition RooProduct.h:48
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooArgSet *nset=nullptr, const RooNumIntConfig *cfg=nullptr, const char *rangeName=nullptr) const override
Create an object that represents the integral of the function over one or more observables listed in ...
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:417
const char * Data() const
Definition TString.h:376
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
RooConstVar & RooConst(double val)
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
BinnedLOutput getBinnedL(RooAbsPdf const &pdf)
std::string getColonSeparatedNameString(RooArgSet const &argSet, char delim=':')
TLine l
Definition textangle.C:4