Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAddPdf.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/** \class RooAddPdf
19 \ingroup Roofitcore
20
21Efficient implementation of a sum of PDFs of the form
22
23\f[
24 \sum_{i=1}^{n} c_i \cdot \mathrm{PDF}_i
25\f]
26
27or
28\f[
29 c_1\cdot\mathrm{PDF}_1 + c_2\cdot\mathrm{PDF}_2 \; + \; ... \; + \; \left( 1-\sum_{i=1}^{n-1}c_i \right) \cdot \mathrm{PDF}_n
30\f]
31
32The first form is for extended likelihood fits, where the
33expected number of events is \f$ \sum_i c_i \f$. The coefficients \f$ c_i \f$
34can either be explicitly provided, or, if all components support
35extended likelihood fits, they can be calculated from the contribution
36of each PDF to the total expected number of events.
37
38In the second form, the sum of the coefficients is required to be 1 or less,
39and the coefficient of the last PDF is calculated automatically from the condition
40that the sum of all coefficients has to be 1.
41
42### Recursive coefficients
43It is also possible to parameterise the coefficients recursively
44
45\f[
46 \sum_{i=1}^n c_i \prod_{j=1}^{i-1} \left[ (1-c_j) \right] \cdot \mathrm{PDF}_i \\
47 = c_1 \cdot \mathrm{PDF}_1 + (1-c_1)\, c_2 \cdot \mathrm{PDF}_2 + \ldots + (1-c_1)\ldots(1-c_{n-1}) \cdot 1 \cdot \mathrm{PDF}_n \\
48\f]
49
50In this form the sum of the coefficients is always less than 1.0
51for all possible values of the individual coefficients between 0 and 1.
52\note Don't pass the \f$ n^\mathrm{th} \f$ coefficient. It is always 1, since the normalisation condition removes one degree of freedom.
53
54RooAddPdf relies on each component PDF to be normalized and will perform
55no normalization other than calculating the proper last coefficient \f$ c_n \f$, if requested.
56An (enforced) condition for this assumption is that each \f$ \mathrm{PDF}_i \f$ is independent of each \f$ c_i \f$.
57
58## Difference between RooAddPdf / RooRealSumFunc / RooRealSumPdf
59- RooAddPdf is a PDF of PDFs, *i.e.* its components need to be normalised and non-negative.
60- RooRealSumPdf is a PDF of functions, *i.e.*, its components can be negative, but their sum cannot be. The normalisation
61 is computed automatically, unless the PDF is extended (see above).
62- RooRealSumFunc is a sum of functions. It is neither normalised, nor need it be positive.
63
64*/
65
66#include <RooAddPdf.h>
67
68#include <RooAddGenContext.h>
69#include <RooAddition.h>
70#include <RooBatchCompute.h>
71#include <RooDataSet.h>
72#include <RooGenericPdf.h>
73#include <RooGlobalFunc.h>
74#include <RooProdPdf.h>
75#include <RooProduct.h>
76#include <RooRatio.h>
77#include <RooRealConstant.h>
78#include <RooRealProxy.h>
79#include <RooRealSumFunc.h>
80#include <RooRealSumPdf.h>
81#include <RooRealVar.h>
83
84#include "RooAddHelpers.h"
85#include "RooFitImplHelpers.h"
86
87#include <ROOT/StringUtils.hxx>
88
89#include <algorithm>
90#include <memory>
91#include <set>
92#include <sstream>
93
94
95
96////////////////////////////////////////////////////////////////////////////////
97/// Dummy constructor
98
99RooAddPdf::RooAddPdf(const char *name, const char *title) :
100 RooAbsPdf(name,title),
101 _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,false,false),
102 _projCacheMgr(this,10),
103 _pdfList("!pdfs","List of PDFs",this),
104 _coefList("!coefficients","List of coefficients",this),
105 _coefErrCount{_errorCount}
106{
108}
109
110
112
113 // Two pdfs with the same name are only allowed in the input list if they are
114 // actually the same object.
115 using PdfInfo = std::pair<std::string,RooAbsArg*>;
116 std::set<PdfInfo> seen;
117 for(auto const& pdf : _pdfList) {
118 PdfInfo elem{pdf->GetName(), pdf};
119 auto comp = [&](PdfInfo const& p){ return p.first == elem.first && p.second != elem.second; };
120 auto found = std::find_if(seen.begin(), seen.end(), comp);
121 if(found != seen.end()) {
122 std::stringstream errorMsg;
123 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
124 << ") pdf list contains pdfs with duplicate name \"" << pdf->GetName() << "\".";
125 coutE(InputArguments) << errorMsg.str() << std::endl;
126 throw std::invalid_argument(errorMsg.str().c_str());
127 }
128 seen.insert(elem);
129 }
130}
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Constructor with two PDFs and one coefficient
135
136RooAddPdf::RooAddPdf(const char *name, const char *title,
137 RooAbsPdf& pdf1, RooAbsPdf& pdf2, RooAbsReal& coef1) :
138 RooAddPdf(name, title)
139{
140 _pdfList.add(pdf1) ;
141 _pdfList.add(pdf2) ;
143
145}
146
147
148////////////////////////////////////////////////////////////////////////////////
149/// Generic constructor from list of PDFs and list of coefficients.
150/// Each pdf list element (i) is paired with coefficient list element (i).
151/// The number of coefficients must be either equal to the number of PDFs,
152/// in which case extended MLL fitting is enabled, or be one less.
153///
154/// All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal
155///
156/// If the recursiveFraction flag is true, the coefficients are interpreted as recursive
157/// coefficients as explained in the class description.
158
159RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList &inPdfList, const RooArgList &inCoefList,
161 : RooAddPdf(name, title)
162{
164
165 if (inPdfList.size()>inCoefList.size()+1 || inPdfList.size()<inCoefList.size()) {
166 std::stringstream errorMsg;
167 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
168 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1.";
169 coutE(InputArguments) << errorMsg.str() << std::endl;
170 throw std::invalid_argument(errorMsg.str().c_str());
171 }
172
173 if (recursiveFractions && inPdfList.size()!=inCoefList.size()+1) {
174 std::stringstream errorMsg;
175 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
176 << "): Recursive fractions option can only be used if Npdf=Ncoef+1.";
177 coutE(InputArguments) << errorMsg.str() << std::endl;
178 throw std::invalid_argument(errorMsg.str());
179 }
180
181 // Constructor with N PDFs and N or N-1 coefs
183
184 auto addRecursiveCoef = [this,&partinCoefList](RooAbsPdf& pdf, RooAbsReal& coef) -> RooAbsReal & {
185 partinCoefList.add(coef) ;
186 if(partinCoefList.size() == 1) {
187 // The first fraction is the first plain fraction
188 return coef;
189 }
190 // The i-th recursive fraction = (1-f1)*(1-f2)*...(fi) and is calculated from the list (f1,...,fi) by RooRecursiveFraction)
191 std::stringstream rfracName;
192 rfracName << GetName() << "_recursive_fraction_" << pdf.GetName() << "_" << partinCoefList.size();
193 auto rfrac = std::make_unique<RooRecursiveFraction>(rfracName.str().c_str(),"Recursive Fraction",partinCoefList) ;
194 auto & rfracRef = *rfrac;
195 addOwnedComponents(std::move(rfrac)) ;
196 return rfracRef;
197 };
198
199 for (auto i = 0u; i < inCoefList.size(); ++i) {
200 auto coef = dynamic_cast<RooAbsReal*>(inCoefList.at(i));
201 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(i));
202 if (inPdfList.at(i) == nullptr) {
203 std::stringstream errorMsg;
204 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
205 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
206 coutE(InputArguments) << errorMsg.str() << std::endl;
207 throw std::invalid_argument(errorMsg.str());
208 }
209 if (!coef) {
210 std::stringstream errorMsg;
211 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") coefficient " << (coef ? coef->GetName() : "") << " is not of type RooAbsReal, ignored";
212 coutE(InputArguments) << errorMsg.str() << std::endl;
213 throw std::invalid_argument(errorMsg.str());
214 }
215 if (!pdf) {
216 std::stringstream errorMsg;
217 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored";
218 coutE(InputArguments) << errorMsg.str() << std::endl;
219 throw std::invalid_argument(errorMsg.str());
220 }
221 _pdfList.add(*pdf) ;
222
223 // Process recursive fraction mode separately
224 _coefList.add(recursiveFractions ? addRecursiveCoef(*pdf, *coef) : *coef);
225 }
226
227 if (inPdfList.size() == inCoefList.size() + 1) {
228 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(inCoefList.size()));
229
230 if (!pdf) {
231 coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") last argument " << inPdfList.at(inCoefList.size())->GetName() << " is not of type RooAbsPdf." << std::endl;
232 throw std::invalid_argument("Last argument for RooAddPdf is not a PDF.");
233 }
234 _pdfList.add(*pdf) ;
235
236 // Process recursive fractions mode. Above, we verified that we don't have a last coefficient
237 if (recursiveFractions) {
239 // In recursive mode we always have Ncoef=Npdf, since we added it just above
241 }
242
243 } else {
245 }
246
248}
249
250
251////////////////////////////////////////////////////////////////////////////////
252/// Generic constructor from list of extended PDFs. There are no coefficients as the expected
253/// number of events from each components determine the relative weight of the PDFs.
254///
255/// All PDFs must inherit from RooAbsPdf.
256
257RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList &inPdfList)
258 : RooAddPdf(name, title)
259{
260 _allExtendable = true;
261
262 // Constructor with N PDFs
263 for (const auto pdfArg : inPdfList) {
264 auto pdf = dynamic_cast<const RooAbsPdf*>(pdfArg);
265
266 if (!pdf) {
267 std::stringstream errorMsg;
268 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "")
269 << " is not of type RooAbsPdf, RooAddPdf constructor call is invalid!";
270 coutE(InputArguments) << errorMsg.str() << std::endl;
271 throw std::invalid_argument(errorMsg.str().c_str());
272 }
273 if (!pdf->canBeExtended()) {
274 std::stringstream errorMsg;
275 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName()
276 << " is not extendable, RooAddPdf constructor call is invalid!";
277 coutE(InputArguments) << errorMsg.str() << std::endl;
278 throw std::invalid_argument(errorMsg.str().c_str());
279 }
280 _pdfList.add(*pdf) ;
281 }
282
284}
285
286
287////////////////////////////////////////////////////////////////////////////////
288/// Copy constructor
289
291 : RooAbsPdf(other, name),
292 _refCoefNorm("!refCoefNorm", this, other._refCoefNorm),
293 _refCoefRangeName((TNamed *)other._refCoefRangeName),
294 _projCacheMgr(other._projCacheMgr, this),
295 _codeReg(other._codeReg),
296 _pdfList("!pdfs", this, other._pdfList),
297 _coefList("!coefficients", this, other._coefList),
298 _haveLastCoef(other._haveLastCoef),
299 _allExtendable(other._allExtendable),
300 _recursive(other._recursive),
301 _coefErrCount(_errorCount)
302{
303
306}
307
308
309////////////////////////////////////////////////////////////////////////////////
310/// By default the interpretation of the fraction coefficients is
311/// performed in the contextual choice of observables. This makes the
312/// shape of the p.d.f explicitly dependent on the choice of
313/// observables. This method instructs RooAddPdf to freeze the
314/// interpretation of the coefficients to be done in the given set of
315/// observables. If frozen, fractions are automatically transformed
316/// from the reference normalization set to the contextual normalization
317/// set by ratios of integrals.
318
320{
321 if (refCoefNorm.empty()) {
322 return ;
323 }
324
325 // Also set an attribute with this information, which is the easiest way to
326 // preserve this in the JSON IO.
328
331
333}
334
340
341// For the JSON IO, we are not storing the _refCoefNorm directly. Instead, it
342// is stored by names in a string attribute. This function should be called
343// internally before _refCoefNorm is used to materialize it from the attribute
344// if necessary.
346{
347 // _refCoefNorm was already materialized
348 if (!_refCoefNorm.empty())
349 return;
350
351 std::vector<std::string> names;
352 if (auto attrib = getStringAttribute("ref_coef_norm")) {
353 names = ROOT::Split(attrib, ",", /*skipEmpty=*/true);
354 } else {
355 return;
356 }
357
359
362 for (std::string const &name : names) {
363 if (RooAbsArg *arg = serverSet.find(name.c_str())) {
364 refCoefNorm.add(*arg);
365 } else {
366 throw std::runtime_error("Internal logic error in RooAddPdf::materializeRefCoefNormFromAttribute()");
367 }
368 }
369
370 const_cast<RooAddPdf *>(this)->fixCoefNormalization(refCoefNorm);
371}
372
373
374////////////////////////////////////////////////////////////////////////////////
375/// By default, fraction coefficients are assumed to refer to the default
376/// fit range. This makes the shape of a RooAddPdf
377/// explicitly dependent on the range of the observables. Calling this function
378/// allows for a range-independent definition of the fractions, because it
379/// ties all coefficients to the given
380/// named range. If the normalisation range is different
381/// from this reference range, the appropriate fraction coefficients
382/// are automatically calculated from the reference fractions by
383/// integrating over the ranges, and comparing these integrals.
384
386{
387 auto* newNamePtr = const_cast<TNamed*>(RooNameReg::ptr(rangeName));
390 }
392}
393
394
395
396////////////////////////////////////////////////////////////////////////////////
397/// Retrieve cache element for the computation of the PDF normalisation.
398/// \param[in] nset Current normalisation set (integration over these variables yields 1).
399/// \param[in] iset Integration set. Variables to be integrated over (if integrations are performed).
400///
401/// If a cache element does not exist, create and fill it on the fly. The cache also contains
402/// - Supplemental normalization terms (in case not all added p.d.f.s have the same observables)
403/// - Projection integrals to calculate transformed fraction coefficients when a frozen reference frame is provided
404/// - Projection integrals for similar transformations when a frozen reference range is provided.
405
406AddCacheElem* RooAddPdf::getProjCache(const RooArgSet* nset, const RooArgSet* iset) const
407{
408 // Check if cache already exists
409 auto cache = static_cast<AddCacheElem*>(_projCacheMgr.getObj(nset,iset,nullptr,normRange()));
410 if (cache) {
411 return cache ;
412 }
413
414 // Make sure _refCoefNorm is defined
416
417 //Create new cache
418 cache = new AddCacheElem{*this, _pdfList, _coefList, nset, iset, _refCoefNorm,
421
423
424 return cache;
425}
426
427
428////////////////////////////////////////////////////////////////////////////////
429/// Update the coefficient values in the given cache element: calculate new remainder
430/// fraction, normalize fractions obtained from extended ML terms to unity, and
431/// multiply the various range and dimensional corrections needed in the
432/// current use context.
433///
434/// param[in] cache The cache element for the given normalization set that
435/// stores the supplementary normalization values and
436/// projection-related objects.
437/// param[in] nset The set of variables to normalize over.
438/// param[in] syncCoefValues If the initial values of the coefficients still
439/// need to be copied from the `_coefList` elements to
440/// the `_coefCache`. True by default.
441
442void RooAddPdf::updateCoefficients(AddCacheElem &cache, const RooArgSet *nset, bool syncCoefValues) const
443{
444 _coefCache.resize(_pdfList.size());
445 if (syncCoefValues) {
446 for (std::size_t i = 0; i < _coefList.size(); ++i) {
447 _coefCache[i] = static_cast<RooAbsReal const &>(_coefList[i]).getVal(nset);
448 }
449 }
450 if (_allExtendable) {
451 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
452 auto &pdf = static_cast<RooAbsPdf &>(_pdfList[i]);
453 _coefCache[i] = pdf.expectedEvents(!_refCoefNorm.empty() ? &_refCoefNorm : nset);
454 }
455 }
456
457 RooAddHelpers::updateCoefficients(*this, _pdfList.size(), _coefCache, _haveLastCoef || _allExtendable, cache,
459}
460
461////////////////////////////////////////////////////////////////////////////////
462/// Look up projection cache and per-PDF norm sets. If a PDF doesn't have a special
463/// norm set, use the `defaultNorm`. If `defaultNorm == nullptr`, use the member
464/// _normSet.
465std::pair<const RooArgSet*, AddCacheElem*> RooAddPdf::getNormAndCache(const RooArgSet* nset) const {
466
467 // Treat empty normalization set and nullptr the same way.
468 if(nset && nset->empty()) nset = nullptr;
469
470 if (nset == nullptr) {
471 // Make sure _refCoefNorm is defined
473
474 if (!_refCoefNorm.empty()) {
475 nset = &_refCoefNorm ;
476 }
477 }
478
479 // A RooAddPdf needs to have a normalization set defined, otherwise its
480 // coefficient will not be uniquely defined. Its shape depends on the
481 // normalization provided. Un-normalized calls to RooAddPdf can happen in
482 // Roofit, when printing the pdf's or when computing integrals. In these case,
483 // if the pdf has a normalization set previously defined (i.e. stored as a
484 // datamember in _copyOfLastNormSet) it should use it by default when the pdf
485 // is evaluated without passing a normalizations set (in pdf->getVal(nullptr) )
486 // In the case of no pre-defined normalization set exists, a warning will be
487 // produced, since the obtained value will be arbitrary. Note that to avoid
488 // unnecessary warning messages, when calling RooAbsPdf::printValue or
489 // RooAbsPdf::graphVizTree, the printing of the warning messages for the
490 // RooFit::Eval topic is explicitly disabled.
491 {
492 // If nset is still nullptr, get the pointer to a copy of the last-used
493 // normalization set. It nset is not nullptr, check whether the copy of
494 // the last-used normalization set needs an update.
495 if(nset == nullptr) {
496 nset = _copyOfLastNormSet.get();
498 _copyOfLastNormSet = std::make_unique<const RooArgSet>(*nset);
500 }
501
502 // If nset is STILL nullptr, print a warning.
503 if (nset == nullptr) {
504 coutW(Eval) << "Evaluating RooAddPdf " << GetName() << " without a defined normalization set. This can lead to ambiguous "
505 "coefficients definition and incorrect results."
506 << " Use RooAddPdf::fixCoefNormalization(nset) to provide a normalization set for "
507 "defining uniquely RooAddPdf coefficients!"
508 << std::endl;
509 }
510 }
511
512
513 AddCacheElem* cache = getProjCache(nset) ;
514
515 return {nset, cache};
516}
517
518
519////////////////////////////////////////////////////////////////////////////////
520/// Calculate and return the current value
521
523{
525 const RooArgSet* nset = normAndCache.first;
526 AddCacheElem* cache = normAndCache.second;
527 updateCoefficients(*cache, nset);
528
529 // Process change in last data set used
530 bool nsetChanged(false) ;
531 if (!isActiveNormSet(nset) || _norm==nullptr) {
533 }
534
535 // Do running sum of coef/pdf pairs, calculate lastCoef.
536 if (isValueDirty() || nsetChanged) {
537 _value = 0.0;
538
539 for (unsigned int i=0; i < _pdfList.size(); ++i) {
540 auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
541 double snormVal = 1.;
542 snormVal = cache->suppNormVal(i);
543
544 double pdfVal = pdf.getVal(nset);
545 if (pdf.isSelectedComp()) {
547 }
548 }
550 }
551
552 return _value;
553}
554
555////////////////////////////////////////////////////////////////////////////////
556/// Compute addition of PDFs in batches.
558{
559 std::span<double> output = ctx.output();
560
561 RooBatchCompute::Config config = ctx.config(this);
562
563 _coefCache.resize(_pdfList.size());
564 for(std::size_t i = 0; i < _coefList.size(); ++i) {
565 auto coefVals = ctx.at(&_coefList[i]);
566 // We don't support per-event coefficients in this function. If the CPU
567 // mode is used, we can just fall back to the RooAbsReal implementation.
568 // With CUDA, we can't do that because the inputs might be on the device.
569 // That's why we throw an exception then.
570 if(coefVals.size() > 1) {
571 if (config.useCuda()) {
572 throw std::runtime_error("The RooAddPdf doesn't support per-event coefficients in CUDA mode yet!");
573 }
575 return;
576 }
577 _coefCache[i] = coefVals[0];
578 }
579
580 std::vector<std::span<const double>> pdfs;
581 std::vector<double> coefs;
582 AddCacheElem* cache = getProjCache(nullptr);
583 RooAddHelpers::updateCoefficients(*this, _pdfList.size(), _coefCache, _haveLastCoef || _allExtendable, *cache,
585
586 for (unsigned int pdfNo = 0; pdfNo < _pdfList.size(); ++pdfNo)
587 {
588 auto pdf = static_cast<RooAbsPdf*>(&_pdfList[pdfNo]);
589 if (pdf->isSelectedComp())
590 {
591 pdfs.push_back(ctx.at(pdf));
592 coefs.push_back(_coefCache[pdfNo] / cache->suppNormVal(pdfNo) );
593 }
594 }
596}
597
598
599////////////////////////////////////////////////////////////////////////////////
600/// Reset error counter to given value, limiting the number
601/// of future error messages for this pdf to 'resetValue'
602
608
609
610
611////////////////////////////////////////////////////////////////////////////////
612/// Check if PDF is valid for given normalization set.
613/// Coefficient and PDF must be non-overlapping, but pdf-coefficient
614/// pairs may overlap each other
615
617{
619}
620
621
622////////////////////////////////////////////////////////////////////////////////
623/// Determine which part (if any) of given integral can be performed analytically.
624/// If any analytical integration is possible, return integration scenario code
625///
626/// RooAddPdf queries each component PDF for its analytical integration capability of the requested
627/// set ('allVars'). It finds the largest common set of variables that can be integrated
628/// by all components. If such a set exists, it reconfirms that each component is capable of
629/// analytically integrating the common set, and combines the components individual integration
630/// codes into a single integration code valid for RooAddPdf.
631
633 const RooArgSet* normSet, const char* rangeName) const
634{
635 // Make sure _refCoefNorm is defined
637
638 RooArgSet allAnalVars(*std::unique_ptr<RooArgSet>{getObservables(allVars)}) ;
639
640 Int_t n(0) ;
641
642 // First iteration, determine what each component can integrate analytically
643 for (const auto pdfArg : _pdfList) {
644 auto pdf = static_cast<const RooAbsPdf *>(pdfArg);
646 pdf->getAnalyticalIntegralWN(allVars,subAnalVars,normSet,rangeName) ;
647
648 // Observables that cannot be integrated analytically by this component are dropped from the common list
649 for (const auto arg : allVars) {
650 if (!subAnalVars.find(arg->GetName()) && pdf->dependsOn(*arg)) {
651 allAnalVars.remove(*arg,true,true) ;
652 }
653 }
654 n++ ;
655 }
656
657 // If no observables can be integrated analytically, return code 0 here
658 if (allAnalVars.empty()) {
659 return 0 ;
660 }
661
662
663 // Now retrieve codes for integration over common set of analytically integrable observables for each component
664 n=0 ;
665 std::vector<Int_t> subCode(_pdfList.size());
666 bool allOK(true) ;
667 for (const auto arg : _pdfList) {
668 auto pdf = static_cast<const RooAbsPdf *>(arg);
670 auto allAnalVars2 = std::unique_ptr<RooArgSet>{pdf->getObservables(allAnalVars)} ;
671 subCode[n] = pdf->getAnalyticalIntegralWN(*allAnalVars2,subAnalVars,normSet,rangeName) ;
672 if (subCode[n]==0 && !allAnalVars2->empty()) {
673 coutE(InputArguments) << "RooAddPdf::getAnalyticalIntegral(" << GetName() << ") WARNING: component PDF " << pdf->GetName()
674 << " advertises inconsistent set of integrals (e.g. (X,Y) but not X or Y individually."
675 << " Distributed analytical integration disabled. Please fix PDF" << std::endl ;
676 allOK = false ;
677 }
678 n++ ;
679 }
680 if (!allOK) {
681 return 0 ;
682 }
683
684 // Mare all analytically integrated observables as such
685 analVars.add(allAnalVars) ;
686
687 // Store set of variables analytically integrated
690
691 return masterCode ;
692}
693
694
695
696////////////////////////////////////////////////////////////////////////////////
697/// Return analytical integral defined by given scenario code
698
699double RooAddPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
700{
701 // WVE needs adaptation to handle new rangeName feature
702 if (code==0) {
703 return getVal(normSet) ;
704 }
705
706 // Retrieve analytical integration subCodes and set of observabels integrated over
707 RooArgSet* intSet = nullptr;
708 const std::vector<Int_t>& subCode = _codeReg.retrieve(code-1,intSet) ;
709 if (subCode.empty()) {
710 std::stringstream errorMsg;
711 errorMsg << "RooAddPdf::analyticalIntegral(" << GetName() << "): ERROR unrecognized integration code, " << code;
712 coutE(InputArguments) << errorMsg.str() << std::endl;
713 throw std::invalid_argument(errorMsg.str().c_str());
714 }
715
716 cxcoutD(Caching) << "RooAddPdf::aiWN(" << GetName() << ") calling getProjCache with nset = " << (normSet?*normSet:RooArgSet()) << std::endl ;
717
718 if ((normSet==nullptr || normSet->empty()) && !_refCoefNorm.empty()) {
719// std::cout << "WVE integration of RooAddPdf without normalization, but have reference set, using ref set for normalization" << std::endl ;
721 }
722
725
726 // Calculate the current value of this object
727 double value(0) ;
728
729 // Do running sum of coef/pdf pairs, calculate lastCoef.
730 double snormVal ;
731
732 //cout << "ROP::aIWN updateCoefCache with rangeName = " << (rangeName?rangeName:"<null>") << std::endl ;
733 for (std::size_t i = 0; i < _pdfList.size(); ++i ) {
734 auto pdf = static_cast<const RooAbsPdf*>(_pdfList.at(i));
735
736 if (_coefCache[i]) {
737 snormVal = cache->suppNormVal(i);
738
739 // WVE swap this?
740 double val = pdf->analyticalIntegralWN(subCode[i],normSet,rangeName) ;
741 if (pdf->isSelectedComp()) {
742 value += val*_coefCache[i]/snormVal ;
743 }
744 }
745 }
746
747 return value ;
748}
749
750
751
752////////////////////////////////////////////////////////////////////////////////
753/// Return the number of expected events, which is either the sum of all coefficients
754/// or the sum of the components extended terms, multiplied with the fraction that
755/// is in the current range w.r.t the reference range
756
757double RooAddPdf::expectedEvents(const RooArgSet* nset) const
758{
759 double expectedTotal{0.0};
760
761 cxcoutD(Caching) << "RooAddPdf::expectedEvents(" << GetName() << ") calling getProjCache with nset = " << (nset?*nset:RooArgSet()) << std::endl ;
762 AddCacheElem& cache = *getProjCache(nset) ;
763 updateCoefficients(cache, nset);
764
765 if (cache.doProjection()) {
766
767 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
768 double ncomp = _allExtendable ? static_cast<RooAbsPdf&>(_pdfList[i]).expectedEvents(nset)
769 : static_cast<RooAbsReal&>(_coefList[i]).getVal(nset);
770 expectedTotal += cache.rangeProjScaleFactor(i) * ncomp ;
771
772 }
773
774 } else {
775
776 if (_allExtendable) {
777 for(auto const& arg : _pdfList) {
778 expectedTotal += static_cast<RooAbsPdf*>(arg)->expectedEvents(nset) ;
779 }
780 } else {
781 for(auto const& arg : _coefList) {
782 expectedTotal += static_cast<RooAbsReal*>(arg)->getVal(nset) ;
783 }
784 }
785
786 }
787 return expectedTotal ;
788}
789
790
791std::unique_ptr<RooAbsReal> RooAddPdf::createExpectedEventsFunc(const RooArgSet *nset) const
792{
793 std::unique_ptr<RooAbsReal> out;
794
795 auto name = std::string(GetName()) + "_expectedEvents";
796 if (_allExtendable) {
798 for (auto *pdf : static_range_cast<RooAbsPdf *>(_pdfList)) {
799 sumSet.addOwned(pdf->createExpectedEventsFunc(nset));
800 }
801 out = std::make_unique<RooAddition>(name.c_str(), name.c_str(), sumSet);
802 out->addOwnedComponents(std::move(sumSet));
803 } else {
804 out = std::make_unique<RooAddition>(name.c_str(), name.c_str(), _coefList);
805 }
806
808
809 // Make sure _refCoefNorm is defined
811
812 if (!_allExtendable) {
813 // If the _refCoefNorm is empty or it's equal to normSet anyway, this is not
814 // a conditional pdf and we don't need to do any transformation. See also
815 // RooAddPdf::compileForNormSet() for more explanations on a similar logic.
816 if (!_refCoefNorm.empty() && !nset->equals(_refCoefNorm)) {
817 prodList.addOwned(std::unique_ptr<RooAbsReal>{createIntegral(*nset, _refCoefNorm)});
818 }
819
820 // Optionally multiply with fractional normalization. I this case, we
821 // replace the original factor stored in "out".
822 if (!_normRange.IsNull()) {
823 std::unique_ptr<RooAbsReal> owner;
824 RooArgList terms;
825 // The integrals own each other in a chain. We do this because it's
826 // not possible to add two objects with the same name via
827 // addOwnedComponents(), and it happens in some user models that some
828 // component pdfs are the same. Hence, the integrals might share names
829 // too and we can't add them all in one go as owned objects of the
830 // final integral sum.
831 for (auto *pdf : static_range_cast<RooAbsPdf *>(_pdfList)) {
832 auto integrl = std::unique_ptr<RooAbsReal>{pdf->createIntegral(*nset, *nset)};
833 auto formulaName = std::string(pdf->GetName()) + "_formulaVar";
834 auto next = std::make_unique<RooFormulaVar>(formulaName.c_str(), "1./x[0]", RooArgList{*integrl});
835 next->addOwnedComponents(std::move(integrl));
836 terms.add(*next);
837 if (owner)
838 next->addOwnedComponents(std::move(owner));
839 owner = std::move(next);
840 }
841 auto fracIntegName = std::string(GetName()) + "_integSum";
842 auto fracInteg =
843 std::make_unique<RooRealSumFunc>(fracIntegName.c_str(), fracIntegName.c_str(), _coefList, terms);
844 fracInteg->addOwnedComponents(std::move(owner));
845
846 out = std::move(fracInteg);
847 }
848 }
849
850 std::string finalName = std::string(out->GetName()) + "_finalized";
851 if (prodList.empty()) {
852 // If there are no additional factors, just return the single factor we have
853 return out;
854 } else {
855 prodList.addOwned(std::move(out));
856 }
857 auto finalOut = std::make_unique<RooProduct>(finalName.c_str(), finalName.c_str(), prodList);
858 finalOut->addOwnedComponents(std::move(prodList));
859 return finalOut;
860}
861
862
863////////////////////////////////////////////////////////////////////////////////
864/// Interface function used by test statistics to freeze choice of observables
865/// for interpretation of fraction coefficients
866
868{
869 // Make sure _refCoefNorm is defined
871
872 if (!force && !_refCoefNorm.empty()) {
873 return ;
874 }
875
876 if (!depSet) {
878 return ;
879 }
880
881 fixCoefNormalization(*std::unique_ptr<RooArgSet>{getObservables(depSet)}) ;
882}
883
884
885
886////////////////////////////////////////////////////////////////////////////////
887/// Interface function used by test statistics to freeze choice of range
888/// for interpretation of fraction coefficients
889
891{
892 if (!force && _refCoefRangeName) {
893 return ;
894 }
895
897}
898
899
900
901////////////////////////////////////////////////////////////////////////////////
902/// Return specialized context to efficiently generate toy events from RooAddPdfs
903/// return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ; // WVE DEBUG
904
906 const RooArgSet* auxProto, bool verbose) const
907{
908 return RooAddGenContext::create(*this,vars,prototype,auxProto,verbose).release();
909}
910
911
912
913////////////////////////////////////////////////////////////////////////////////
914/// Loop over components for plot sampling hints and merge them if there are multiple
915
916std::list<double>* RooAddPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
917{
918 return RooRealSumPdf::plotSamplingHint(_pdfList, obs, xlo, xhi);
919}
920
921
922////////////////////////////////////////////////////////////////////////////////
923/// Loop over components for plot sampling hints and merge them if there are multiple
924
925std::list<double>* RooAddPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
926{
927 return RooRealSumPdf::binBoundaries(_pdfList, obs, xlo, xhi);
928}
929
930
931////////////////////////////////////////////////////////////////////////////////
932/// If all components that depend on obs are binned, so is their sum.
937
938
939////////////////////////////////////////////////////////////////////////////////
940/// Label OK'ed components of a RooAddPdf with cache-and-track
941
946
947
948
949////////////////////////////////////////////////////////////////////////////////
950/// Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the
951/// product operator construction
952
953void RooAddPdf::printMetaArgs(std::ostream& os) const
954{
956}
957
958
960 bool nameChange, bool isRecursiveStep)
961{
962 // If a server is redirected, the cached normalization set might not point
963 // to the right observables anymore. We need to reset it.
964 _copyOfLastNormSet.reset();
966}
967
968
969std::unique_ptr<RooAbsArg>
971{
972 // Make sure _refCoefNorm is defined
974
975 // Instead of cloning this RooAddPdf directly, we have to massage it a bit.
976 // In the case of extended pdfs, the coefficients should be set to functions
977 // representing the expected number of events, so we don't have to fall back
978 // to legacy code paths that don't support evaluation with the
979 // RooFit::EvalContext, like RooAbsPdf::expectedEvents().
981 if (_allExtendable) {
982 for (auto *pdf : static_range_cast<RooAbsPdf *>(_pdfList)) {
983 coefListNew.addOwned(pdf->createExpectedEventsFunc(!_refCoefNorm.empty() ? &_refCoefNorm : &normSet));
984 }
985 } else {
986 coefListNew.add(coefList());
987 }
988 auto newArg = std::make_unique<RooAddPdf>(GetName(), GetTitle(), pdfList(), coefListNew, _recursive);
989 // Copy some other info that the RooAddPdf copy constructor would otherwise take care of.
990 newArg->setNormRange(normRange());
991 newArg->_codeReg = _codeReg;
992 if (!_refCoefNorm.empty()) {
993 newArg->fixCoefNormalization(_refCoefNorm);
994 }
995 if (_refCoefRangeName) {
996 newArg->fixCoefRange(getCoefRange());
997 }
998
1000
1001 // If we set the normalization ranges of the component pdfs to the
1002 // normalization range of the RooAddPdf, the RooAddPdf doesn't need to
1003 // compute as many correction factor integrals internally.
1004 // Remember the previous normalizations ranges to reset later.
1005 class ResetNormRangesRAII {
1006 public:
1007 ResetNormRangesRAII(RooAbsCollection const &pdfs, std::string const &normRange)
1008 {
1009 _componentPdfs.reserve(pdfs.size());
1010 _oldNormRanges.reserve(pdfs.size());
1011
1012 bool isMultiRange = normRange.find(',') != std::string::npos;
1013
1015
1016 // The RooProdPdf is not able to deal with a multi-range _normRange
1017 // for now, so skip the changing the range in this case.
1018 bool changeRange = !(isMultiRange && dynamic_cast<RooProdPdf *>(componentPdf));
1019
1020 const char *old = componentPdf->normRange();
1021 const char *newVal = changeRange ? normRange.c_str() : old;
1022 componentPdf->setNormRange(newVal);
1023 _componentPdfs.emplace_back(componentPdf);
1024 _oldNormRanges.emplace_back(old ? old : "");
1025 }
1026 }
1028 {
1029 for (std::size_t i = 0; i < _componentPdfs.size(); ++i) {
1030 _componentPdfs[i]->setNormRange(_oldNormRanges[i].c_str());
1031 }
1032 }
1033
1034 private:
1035 std::vector<RooAbsPdf *> _componentPdfs;
1036 std::vector<std::string> _oldNormRanges;
1038
1039 // In case conditional observables, e.g. p(x|y), the _refCoefNorm is set to
1040 // all observables (x, y) and the normSet doesn't contain the conditional
1041 // observables (so it only contains x in this example).
1042
1043 // If the _refCoefNorm is empty or it's equal to normSet anyway, this is not
1044 // a conditional pdf and we don't need to do any transformation.
1045 if (_refCoefNorm.empty() || normSet.equals(_refCoefNorm)) {
1047 return newArg;
1048 }
1049
1050 // In the conditional case, things become more complicated. The original
1051 // getValV() method is covering this case with very complicated logic,
1052 // caching multiple new RooFit objects to scale the individual coefficients
1053 // of the RooAddPdf.
1054 //
1055 // However, it's not complicated what we need to do mathematically:
1056 //
1057 // Since:
1058 // 1. p(x, y) = p(x | y) * p(y)
1059 // 2. p(y) = Integral of p(x, y) over x
1060 //
1061 // We conclude:
1062 // p(x, y)
1063 // p(x | y) = --------------------------
1064 // Integral of p(x, y) over x
1065 //
1066 // What follows is the implementation of this formula in RooFit. By doing
1067 // this here in compileForNormSet(), we don't invoke the old RooAddPdf
1068 // projection caches (note that no conditional pdfs are on the right hand
1069 // side of the equation).
1070 std::string finalName = std::string(GetName()) + "_conditional";
1071 std::unique_ptr<RooAbsReal> denom{newArg->createIntegral(normSet, _refCoefNorm)};
1072 auto finalArg = std::make_unique<RooGenericPdf>(finalName.c_str(), "@0/@1", RooArgList{*newArg, *denom});
1074 ctx.markAsCompiled(*denom);
1077 finalArg->addOwnedComponents(std::move(newArg));
1078 finalArg->addOwnedComponents(std::move(denom));
1079 return finalArg;
1080}
#define cxcoutD(a)
#define coutW(a)
#define coutE(a)
#define TRACE_CREATE
Definition RooTrace.h:23
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
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
const_iterator begin() const
const_iterator end() const
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:76
void clearValueAndShapeDirty() const
Definition RooAbsArg.h:535
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
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'.
bool isValueDirty() const
Definition RooAbsArg.h:356
Abstract container object that can hold multiple RooAbsArg objects.
bool equals(const RooAbsCollection &otherColl) const
Check if this and other collection have identically-named contents.
RooFit::UniqueId< RooAbsCollection > const & uniqueId() const
Returns a unique ID that is different for every instantiated RooAbsCollection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
virtual bool syncNormalization(const RooArgSet *dset, bool adjustProxies=true) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
bool isActiveNormSet(RooArgSet const *normSet) const
Checks if normSet is the currently active normalization set of this PDF, meaning is exactly the same ...
Definition RooAbsPdf.h:295
TString _normRange
Normalization range.
Definition RooAbsPdf.h:336
RooAbsReal * _norm
! Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:313
const char * normRange() const
Definition RooAbsPdf.h:246
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Hook function intercepting redirectServer calls.
static Int_t _verboseEval
Definition RooAbsPdf.h:308
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:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
friend class AddCacheElem
Definition RooAbsReal.h:416
double _value
Cache for current value of object.
Definition RooAbsReal.h:545
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 ...
virtual void doEval(RooFit::EvalContext &) const
Base function for computing multiple values of a RooAbsReal.
static std::unique_ptr< RooAbsGenContext > create(const Pdf_t &pdf, const RooArgSet &vars, const RooDataSet *prototype, const RooArgSet *auxProto, bool verbose)
Returns a RooAddGenContext if possible, or, if the RooAddGenContext doesn't support this particular R...
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooListProxy _coefList
List of coefficients.
Definition RooAddPdf.h:130
bool _allExtendable
Flag indicating if all PDF components are extendable.
Definition RooAddPdf.h:134
void doEval(RooFit::EvalContext &) const override
Compute addition of PDFs in batches.
RooAICRegistry _codeReg
! Registry of component analytical integration codes
Definition RooAddPdf.h:127
RooFit::UniqueId< RooArgSet >::Value_t _idOfLastUsedNormSet
!
Definition RooAddPdf.h:143
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral defined by given scenario code.
const char * getCoefRange() const
Definition RooAddPdf.h:83
std::unique_ptr< const RooArgSet > _copyOfLastNormSet
!
Definition RooAddPdf.h:144
void updateCoefficients(AddCacheElem &cache, const RooArgSet *nset, bool syncCoefValues=true) const
Update the coefficient values in the given cache element: calculate new remainder fraction,...
Int_t _coefErrCount
! Coefficient error counter
Definition RooAddPdf.h:137
bool _haveLastCoef
Flag indicating if last PDFs coefficient was supplied in the constructor.
Definition RooAddPdf.h:133
void selectNormalization(const RooArgSet *depSet=nullptr, bool force=false) override
Interface function used by test statistics to freeze choice of observables for interpretation of frac...
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the produ...
void finalizeConstruction()
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components of a RooAddPdf with cache-and-track.
const RooArgList & coefList() const
Definition RooAddPdf.h:74
bool _recursive
Flag indicating is fractions are treated recursively.
Definition RooAddPdf.h:135
RooObjCacheManager _projCacheMgr
! Manager of cache with coefficient projections and transformations
Definition RooAddPdf.h:108
void materializeRefCoefNormFromAttribute() const
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Return specialized context to efficiently generate toy events from RooAddPdfs return RooAbsPdf::genCo...
bool checkObservables(const RooArgSet *nset) const override
Check if PDF is valid for given normalization set.
void fixCoefNormalization(const RooArgSet &refCoefNorm)
By default the interpretation of the fraction coefficients is performed in the contextual choice of o...
std::pair< const RooArgSet *, AddCacheElem * > getNormAndCache(const RooArgSet *nset) const
Look up projection cache and per-PDF norm sets.
RooSetProxy _refCoefNorm
Reference observable set for coefficient interpretation.
Definition RooAddPdf.h:102
void selectNormalizationRange(const char *rangeName=nullptr, bool force=false) override
Interface function used by test statistics to freeze choice of range for interpretation of fraction c...
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.
void resetErrorCounters(Int_t resetValue=10) override
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
double expectedEvents(const RooArgSet *nset) const override
Return expected number of events for extended likelihood calculation, which is the sum of all coeffic...
double getValV(const RooArgSet *set=nullptr) const override
Calculate and return the current value.
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
void fixCoefRange(const char *rangeName)
By default, fraction coefficients are assumed to refer to the default fit range.
AddCacheElem * getProjCache(const RooArgSet *nset, const RooArgSet *iset=nullptr) const
Retrieve cache element for the computation of the PDF normalisation.
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Loop over components for plot sampling hints and merge them if there are multiple.
RooListProxy _pdfList
List of component PDFs.
Definition RooAddPdf.h:129
TNamed * _refCoefRangeName
Reference range name for coefficient interpretation.
Definition RooAddPdf.h:103
const RooArgList & pdfList() const
Definition RooAddPdf.h:70
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,...
bool isBinnedDistribution(const RooArgSet &obs) const override
If all components that depend on obs are binned, so is their sum.
bool redirectServersHook(const RooAbsCollection &, bool, bool, bool) override
Hook function intercepting redirectServer calls.
std::vector< double > _coefCache
! Transient cache with transformed values of coefficients
Definition RooAddPdf.h:105
const RooArgSet & getCoefNormalization() const
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Loop over components for plot sampling hints and merge them if there are multiple.
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:24
Minimal configuration struct to steer the evaluation of a single node with the RooBatchCompute librar...
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
void reset()
Clear the cache.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
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...
Container class to hold unbinned data.
Definition RooDataSet.h:32
void markAsCompiled(RooAbsArg &arg) const
void compileServers(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
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.
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:36
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components of a RooRealSumPdf with cache-and-track.
bool checkObservables(const RooArgSet *nset) const override
Check if FUNC is valid for given normalization set.
std::list< double > * plotSamplingHint(RooAbsRealLValue &, double, double) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Retrieve bin boundaries if this distribution is binned in obs.
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the p...
bool isBinnedDistribution(const RooArgSet &obs) const override
Check if all components that depend on obs are binned.
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:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
const char * Data() const
Definition TString.h:384
Bool_t IsNull() const
Definition TString.h:422
RooConstVar & RooConst(double val)
const Int_t n
Definition legend1.C:16
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
void getSortedComputationGraph(RooAbsArg const &func, RooArgSet &out)
std::string getColonSeparatedNameString(RooArgSet const &argSet, char delim=':')
unsigned long Value_t
Definition UniqueId.h:41
static void output()