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 <RooProduct.h>
75#include <RooRatio.h>
76#include <RooRealConstant.h>
77#include <RooRealProxy.h>
78#include <RooRealSumFunc.h>
79#include <RooRealSumPdf.h>
80#include <RooRealVar.h>
82
83#include "RooAddHelpers.h"
84#include "RooFitImplHelpers.h"
85
86#include <ROOT/StringUtils.hxx>
87
88#include <algorithm>
89#include <memory>
90#include <set>
91#include <sstream>
92
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) ;
142 _coefList.add(coef1) ;
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, bool recursiveFractions) :
160 RooAddPdf(name,title)
161{
162 _recursive = recursiveFractions;
163
164 if (inPdfList.size()>inCoefList.size()+1 || inPdfList.size()<inCoefList.size()) {
165 std::stringstream errorMsg;
166 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
167 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1.";
168 coutE(InputArguments) << errorMsg.str() << std::endl;
169 throw std::invalid_argument(errorMsg.str().c_str());
170 }
171
172 if (recursiveFractions && inPdfList.size()!=inCoefList.size()+1) {
173 std::stringstream errorMsg;
174 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
175 << "): Recursive fractions option can only be used if Npdf=Ncoef+1.";
176 coutE(InputArguments) << errorMsg.str() << std::endl;
177 throw std::invalid_argument(errorMsg.str());
178 }
179
180 // Constructor with N PDFs and N or N-1 coefs
181 RooArgList partinCoefList ;
182
183 auto addRecursiveCoef = [this,&partinCoefList](RooAbsPdf& pdf, RooAbsReal& coef) -> RooAbsReal & {
184 partinCoefList.add(coef) ;
185 if(partinCoefList.size() == 1) {
186 // The first fraction is the first plain fraction
187 return coef;
188 }
189 // The i-th recursive fraction = (1-f1)*(1-f2)*...(fi) and is calculated from the list (f1,...,fi) by RooRecursiveFraction)
190 std::stringstream rfracName;
191 rfracName << GetName() << "_recursive_fraction_" << pdf.GetName() << "_" << partinCoefList.size();
192 auto rfrac = std::make_unique<RooRecursiveFraction>(rfracName.str().c_str(),"Recursive Fraction",partinCoefList) ;
193 auto & rfracRef = *rfrac;
194 addOwnedComponents(std::move(rfrac)) ;
195 return rfracRef;
196 };
197
198 for (auto i = 0u; i < inCoefList.size(); ++i) {
199 auto coef = dynamic_cast<RooAbsReal*>(inCoefList.at(i));
200 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(i));
201 if (inPdfList.at(i) == nullptr) {
202 std::stringstream errorMsg;
203 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
204 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
205 coutE(InputArguments) << errorMsg.str() << std::endl;
206 throw std::invalid_argument(errorMsg.str());
207 }
208 if (!coef) {
209 std::stringstream errorMsg;
210 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") coefficient " << (coef ? coef->GetName() : "") << " is not of type RooAbsReal, ignored";
211 coutE(InputArguments) << errorMsg.str() << std::endl;
212 throw std::invalid_argument(errorMsg.str());
213 }
214 if (!pdf) {
215 std::stringstream errorMsg;
216 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored";
217 coutE(InputArguments) << errorMsg.str() << std::endl;
218 throw std::invalid_argument(errorMsg.str());
219 }
220 _pdfList.add(*pdf) ;
221
222 // Process recursive fraction mode separately
223 _coefList.add(recursiveFractions ? addRecursiveCoef(*pdf, *coef) : *coef);
224 }
225
226 if (inPdfList.size() == inCoefList.size() + 1) {
227 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(inCoefList.size()));
228
229 if (!pdf) {
230 coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") last argument " << inPdfList.at(inCoefList.size())->GetName() << " is not of type RooAbsPdf." << std::endl;
231 throw std::invalid_argument("Last argument for RooAddPdf is not a PDF.");
232 }
233 _pdfList.add(*pdf) ;
234
235 // Process recursive fractions mode. Above, we verified that we don't have a last coefficient
236 if (recursiveFractions) {
237 _coefList.add(addRecursiveCoef(*pdf, RooFit::RooConst(1)));
238 // In recursive mode we always have Ncoef=Npdf, since we added it just above
239 _haveLastCoef=true ;
240 }
241
242 } else {
243 _haveLastCoef=true ;
244 }
245
247}
248
249
250////////////////////////////////////////////////////////////////////////////////
251/// Generic constructor from list of extended PDFs. There are no coefficients as the expected
252/// number of events from each components determine the relative weight of the PDFs.
253///
254/// All PDFs must inherit from RooAbsPdf.
255
256RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList& inPdfList) :
257 RooAddPdf(name,title)
258{
259 _allExtendable = true;
260
261 // Constructor with N PDFs
262 for (const auto pdfArg : inPdfList) {
263 auto pdf = dynamic_cast<const RooAbsPdf*>(pdfArg);
264
265 if (!pdf) {
266 std::stringstream errorMsg;
267 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "")
268 << " is not of type RooAbsPdf, RooAddPdf constructor call is invalid!";
269 coutE(InputArguments) << errorMsg.str() << std::endl;
270 throw std::invalid_argument(errorMsg.str().c_str());
271 }
272 if (!pdf->canBeExtended()) {
273 std::stringstream errorMsg;
274 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName()
275 << " is not extendable, RooAddPdf constructor call is invalid!";
276 coutE(InputArguments) << errorMsg.str() << std::endl;
277 throw std::invalid_argument(errorMsg.str().c_str());
278 }
279 _pdfList.add(*pdf) ;
280 }
281
283}
284
285
286////////////////////////////////////////////////////////////////////////////////
287/// Copy constructor
288
289RooAddPdf::RooAddPdf(const RooAddPdf& other, const char* name) :
290 RooAbsPdf(other,name),
291 _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
292 _refCoefRangeName((TNamed*)other._refCoefRangeName),
293 _projCacheMgr(other._projCacheMgr,this),
294 _codeReg(other._codeReg),
295 _pdfList("!pdfs",this,other._pdfList),
296 _coefList("!coefficients",this,other._coefList),
297 _haveLastCoef(other._haveLastCoef),
298 _allExtendable(other._allExtendable),
299 _recursive(other._recursive)
300{
304}
305
306
307////////////////////////////////////////////////////////////////////////////////
308/// By default the interpretation of the fraction coefficients is
309/// performed in the contextual choice of observables. This makes the
310/// shape of the p.d.f explicitly dependent on the choice of
311/// observables. This method instructs RooAddPdf to freeze the
312/// interpretation of the coefficients to be done in the given set of
313/// observables. If frozen, fractions are automatically transformed
314/// from the reference normalization set to the contextual normalization
315/// set by ratios of integrals.
316
318{
319 if (refCoefNorm.empty()) {
320 return ;
321 }
322
323 // Also set an attribute with this information, which is the easiest way to
324 // preserve this in the JSON IO.
325 setStringAttribute("ref_coef_norm", RooHelpers::getColonSeparatedNameString(refCoefNorm, ',').c_str());
326
328 _refCoefNorm.add(refCoefNorm) ;
329
331}
332
334{
336 return _refCoefNorm;
337}
338
339// For the JSON IO, we are not storing the _refCoefNorm directly. Instead, it
340// is stored by names in a string attribute. This function should be called
341// internally before _refCoefNorm is used to materialize it from the attribute
342// if necessary.
344{
345 // _refCoefNorm was already materialized
346 if (!_refCoefNorm.empty())
347 return;
348
349 std::vector<std::string> names;
350 if (auto attrib = getStringAttribute("ref_coef_norm")) {
351 names = ROOT::Split(attrib, ",", /*skipEmpty=*/true);
352 } else {
353 return;
354 }
355
356 RooArgSet refCoefNorm;
357
358 RooArgSet serverSet;
360 for (std::string const &name : names) {
361 if (RooAbsArg *arg = serverSet.find(name.c_str())) {
362 refCoefNorm.add(*arg);
363 } else {
364 throw std::runtime_error("Internal logic error in RooAddPdf::materializeRefCoefNormFromAttribute()");
365 }
366 }
367
368 const_cast<RooAddPdf *>(this)->fixCoefNormalization(refCoefNorm);
369}
370
371
372////////////////////////////////////////////////////////////////////////////////
373/// By default, fraction coefficients are assumed to refer to the default
374/// fit range. This makes the shape of a RooAddPdf
375/// explicitly dependent on the range of the observables. Calling this function
376/// allows for a range-independent definition of the fractions, because it
377/// ties all coefficients to the given
378/// named range. If the normalisation range is different
379/// from this reference range, the appropriate fraction coefficients
380/// are automatically calculated from the reference fractions by
381/// integrating over the ranges, and comparing these integrals.
382
383void RooAddPdf::fixCoefRange(const char* rangeName)
384{
385 auto* newNamePtr = const_cast<TNamed*>(RooNameReg::ptr(rangeName));
386 if(newNamePtr != _refCoefRangeName) {
388 }
389 _refCoefRangeName = newNamePtr;
390}
391
392
393
394////////////////////////////////////////////////////////////////////////////////
395/// Retrieve cache element for the computation of the PDF normalisation.
396/// \param[in] nset Current normalisation set (integration over these variables yields 1).
397/// \param[in] iset Integration set. Variables to be integrated over (if integrations are performed).
398/// \param[in] rangeName Reference range for the integrals.
399///
400/// If a cache element does not exist, create and fill it on the fly. The cache also contains
401/// - Supplemental normalization terms (in case not all added p.d.f.s have the same observables)
402/// - Projection integrals to calculate transformed fraction coefficients when a frozen reference frame is provided
403/// - Projection integrals for similar transformations when a frozen reference range is provided.
404
406{
407 // Check if cache already exists
408 auto cache = static_cast<AddCacheElem*>(_projCacheMgr.getObj(nset,iset,nullptr,normRange()));
409 if (cache) {
410 return cache ;
411 }
412
413 // Make sure _refCoefNorm is defined
415
416 //Create new cache
417 cache = new AddCacheElem{*this, _pdfList, _coefList, nset, iset, _refCoefNorm,
420
421 _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(normRange())) ;
422
423 return cache;
424}
425
426
427////////////////////////////////////////////////////////////////////////////////
428/// Update the coefficient values in the given cache element: calculate new remainder
429/// fraction, normalize fractions obtained from extended ML terms to unity, and
430/// multiply the various range and dimensional corrections needed in the
431/// current use context.
432///
433/// param[in] cache The cache element for the given normalization set that
434/// stores the supplementary normalization values and
435/// projection-related objects.
436/// param[in] nset The set of variables to normalize over.
437/// param[in] syncCoefValues If the initial values of the coefficients still
438/// need to be copied from the `_coefList` elements to
439/// the `_coefCache`. True by default.
440
441void RooAddPdf::updateCoefficients(AddCacheElem& cache, const RooArgSet* nset, bool syncCoefValues) const
442{
443 _coefCache.resize(_pdfList.size());
444 if(syncCoefValues) {
445 for(std::size_t i = 0; i < _coefList.size(); ++i) {
446 _coefCache[i] = static_cast<RooAbsReal const&>(_coefList[i]).getVal(nset);
447 }
448 }
451}
452
453////////////////////////////////////////////////////////////////////////////////
454/// Look up projection cache and per-PDF norm sets. If a PDF doesn't have a special
455/// norm set, use the `defaultNorm`. If `defaultNorm == nullptr`, use the member
456/// _normSet.
457std::pair<const RooArgSet*, AddCacheElem*> RooAddPdf::getNormAndCache(const RooArgSet* nset) const {
458
459 // Treat empty normalization set and nullptr the same way.
460 if(nset && nset->empty()) nset = nullptr;
461
462 if (nset == nullptr) {
463 // Make sure _refCoefNorm is defined
465
466 if (!_refCoefNorm.empty()) {
467 nset = &_refCoefNorm ;
468 }
469 }
470
471 // A RooAddPdf needs to have a normalization set defined, otherwise its
472 // coefficient will not be uniquely defined. Its shape depends on the
473 // normalization provided. Un-normalized calls to RooAddPdf can happen in
474 // Roofit, when printing the pdf's or when computing integrals. In these case,
475 // if the pdf has a normalization set previously defined (i.e. stored as a
476 // datamember in _copyOfLastNormSet) it should use it by default when the pdf
477 // is evaluated without passing a normalizations set (in pdf->getVal(nullptr) )
478 // In the case of no pre-defined normalization set exists, a warning will be
479 // produced, since the obtained value will be arbitrary. Note that to avoid
480 // unnecessary warning messages, when calling RooAbsPdf::printValue or
481 // RooAbsPdf::graphVizTree, the printing of the warning messages for the
482 // RooFit::Eval topic is explicitly disabled.
483 {
484 // If nset is still nullptr, get the pointer to a copy of the last-used
485 // normalization set. It nset is not nullptr, check whether the copy of
486 // the last-used normalization set needs an update.
487 if(nset == nullptr) {
488 nset = _copyOfLastNormSet.get();
490 _copyOfLastNormSet = std::make_unique<const RooArgSet>(*nset);
492 }
493
494 // If nset is STILL nullptr, print a warning.
495 if (nset == nullptr) {
496 coutW(Eval) << "Evaluating RooAddPdf " << GetName() << " without a defined normalization set. This can lead to ambiguous "
497 "coefficients definition and incorrect results."
498 << " Use RooAddPdf::fixCoefNormalization(nset) to provide a normalization set for "
499 "defining uniquely RooAddPdf coefficients!"
500 << std::endl;
501 }
502 }
503
504
505 AddCacheElem* cache = getProjCache(nset) ;
506
507 return {nset, cache};
508}
509
510
511////////////////////////////////////////////////////////////////////////////////
512/// Calculate and return the current value
513
514double RooAddPdf::getValV(const RooArgSet* normSet) const
515{
516 auto normAndCache = getNormAndCache(normSet);
517 const RooArgSet* nset = normAndCache.first;
518 AddCacheElem* cache = normAndCache.second;
519 updateCoefficients(*cache, nset);
520
521 // Process change in last data set used
522 bool nsetChanged(false) ;
523 if (!isActiveNormSet(nset) || _norm==nullptr) {
524 nsetChanged = syncNormalization(nset) ;
525 }
526
527 // Do running sum of coef/pdf pairs, calculate lastCoef.
528 if (isValueDirty() || nsetChanged) {
529 _value = 0.0;
530
531 for (unsigned int i=0; i < _pdfList.size(); ++i) {
532 auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
533 double snormVal = 1.;
534 snormVal = cache->suppNormVal(i);
535
536 double pdfVal = pdf.getVal(nset);
537 if (pdf.isSelectedComp()) {
538 _value += pdfVal*_coefCache[i]/snormVal;
539 }
540 }
542 }
543
544 return _value;
545}
546
548{
550}
551
552////////////////////////////////////////////////////////////////////////////////
553/// Compute addition of PDFs in batches.
554void RooAddPdf::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
555{
556 RooBatchCompute::Config config = dataMap.config(this);
557
558 _coefCache.resize(_pdfList.size());
559 for(std::size_t i = 0; i < _coefList.size(); ++i) {
560 auto coefVals = dataMap.at(&_coefList[i]);
561 // We don't support per-event coefficients in this function. If the CPU
562 // mode is used, we can just fall back to the RooAbsReal implementation.
563 // With CUDA, we can't do that because the inputs might be on the device.
564 // That's why we throw an exception then.
565 if(coefVals.size() > 1) {
566 if (config.useCuda()) {
567 throw std::runtime_error("The RooAddPdf doesn't support per-event coefficients in CUDA mode yet!");
568 }
569 RooAbsReal::computeBatch(output, nEvents, dataMap);
570 return;
571 }
572 _coefCache[i] = coefVals[0];
573 }
574
577 AddCacheElem* cache = getProjCache(nullptr);
578 // We don't sync the coefficient values from the _coefList to the _coefCache
579 // because we have already done it using the dataMap.
580 updateCoefficients(*cache, nullptr, /*syncCoefValues=*/false);
581
582 for (unsigned int pdfNo = 0; pdfNo < _pdfList.size(); ++pdfNo)
583 {
584 auto pdf = static_cast<RooAbsPdf*>(&_pdfList[pdfNo]);
585 if (pdf->isSelectedComp())
586 {
587 pdfs.push_back(dataMap.at(pdf));
588 coefs.push_back(_coefCache[pdfNo] / cache->suppNormVal(pdfNo) );
589 }
590 }
591 RooBatchCompute::compute(config, RooBatchCompute::AddPdf, output, nEvents, pdfs, coefs);
592}
593
594
595////////////////////////////////////////////////////////////////////////////////
596/// Reset error counter to given value, limiting the number
597/// of future error messages for this pdf to 'resetValue'
598
600{
602 _coefErrCount = resetValue ;
603}
604
605
606
607////////////////////////////////////////////////////////////////////////////////
608/// Check if PDF is valid for given normalization set.
609/// Coefficient and PDF must be non-overlapping, but pdf-coefficient
610/// pairs may overlap each other
611
613{
615}
616
617
618////////////////////////////////////////////////////////////////////////////////
619/// Determine which part (if any) of given integral can be performed analytically.
620/// If any analytical integration is possible, return integration scenario code
621///
622/// RooAddPdf queries each component PDF for its analytical integration capability of the requested
623/// set ('allVars'). It finds the largest common set of variables that can be integrated
624/// by all components. If such a set exists, it reconfirms that each component is capable of
625/// analytically integrating the common set, and combines the components individual integration
626/// codes into a single integration code valid for RooAddPdf.
627
629 const RooArgSet* normSet, const char* rangeName) const
630{
631 // Make sure _refCoefNorm is defined
633
634 RooArgSet allAnalVars(*std::unique_ptr<RooArgSet>{getObservables(allVars)}) ;
635
636 Int_t n(0) ;
637
638 // First iteration, determine what each component can integrate analytically
639 for (const auto pdfArg : _pdfList) {
640 auto pdf = static_cast<const RooAbsPdf *>(pdfArg);
641 RooArgSet subAnalVars ;
642 pdf->getAnalyticalIntegralWN(allVars,subAnalVars,normSet,rangeName) ;
643
644 // Observables that cannot be integrated analytically by this component are dropped from the common list
645 for (const auto arg : allVars) {
646 if (!subAnalVars.find(arg->GetName()) && pdf->dependsOn(*arg)) {
647 allAnalVars.remove(*arg,true,true) ;
648 }
649 }
650 n++ ;
651 }
652
653 // If no observables can be integrated analytically, return code 0 here
654 if (allAnalVars.empty()) {
655 return 0 ;
656 }
657
658
659 // Now retrieve codes for integration over common set of analytically integrable observables for each component
660 n=0 ;
661 std::vector<Int_t> subCode(_pdfList.size());
662 bool allOK(true) ;
663 for (const auto arg : _pdfList) {
664 auto pdf = static_cast<const RooAbsPdf *>(arg);
665 RooArgSet subAnalVars ;
666 auto allAnalVars2 = std::unique_ptr<RooArgSet>{pdf->getObservables(allAnalVars)} ;
667 subCode[n] = pdf->getAnalyticalIntegralWN(*allAnalVars2,subAnalVars,normSet,rangeName) ;
668 if (subCode[n]==0 && !allAnalVars2->empty()) {
669 coutE(InputArguments) << "RooAddPdf::getAnalyticalIntegral(" << GetName() << ") WARNING: component PDF " << pdf->GetName()
670 << " advertises inconsistent set of integrals (e.g. (X,Y) but not X or Y individually."
671 << " Distributed analytical integration disabled. Please fix PDF" << std::endl ;
672 allOK = false ;
673 }
674 n++ ;
675 }
676 if (!allOK) {
677 return 0 ;
678 }
679
680 // Mare all analytically integrated observables as such
681 analVars.add(allAnalVars) ;
682
683 // Store set of variables analytically integrated
684 RooArgSet* intSet = new RooArgSet(allAnalVars) ;
685 Int_t masterCode = _codeReg.store(subCode,intSet)+1 ;
686
687 return masterCode ;
688}
689
690
691
692////////////////////////////////////////////////////////////////////////////////
693/// Return analytical integral defined by given scenario code
694
695double RooAddPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
696{
697 // WVE needs adaptation to handle new rangeName feature
698 if (code==0) {
699 return getVal(normSet) ;
700 }
701
702 // Retrieve analytical integration subCodes and set of observabels integrated over
703 RooArgSet* intSet ;
704 const std::vector<Int_t>& subCode = _codeReg.retrieve(code-1,intSet) ;
705 if (subCode.empty()) {
706 std::stringstream errorMsg;
707 errorMsg << "RooAddPdf::analyticalIntegral(" << GetName() << "): ERROR unrecognized integration code, " << code;
708 coutE(InputArguments) << errorMsg.str() << std::endl;
709 throw std::invalid_argument(errorMsg.str().c_str());
710 }
711
712 cxcoutD(Caching) << "RooAddPdf::aiWN(" << GetName() << ") calling getProjCache with nset = " << (normSet?*normSet:RooArgSet()) << std::endl ;
713
714 if ((normSet==nullptr || normSet->empty()) && !_refCoefNorm.empty()) {
715// cout << "WVE integration of RooAddPdf without normalization, but have reference set, using ref set for normalization" << std::endl ;
716 normSet = &_refCoefNorm ;
717 }
718
719 AddCacheElem* cache = getProjCache(normSet,intSet);
720 updateCoefficients(*cache,normSet);
721
722 // Calculate the current value of this object
723 double value(0) ;
724
725 // Do running sum of coef/pdf pairs, calculate lastCoef.
726 double snormVal ;
727
728 //cout << "ROP::aIWN updateCoefCache with rangeName = " << (rangeName?rangeName:"<null>") << std::endl ;
729 for (std::size_t i = 0; i < _pdfList.size(); ++i ) {
730 auto pdf = static_cast<const RooAbsPdf*>(_pdfList.at(i));
731
732 if (_coefCache[i]) {
733 snormVal = cache->suppNormVal(i);
734
735 // WVE swap this?
736 double val = pdf->analyticalIntegralWN(subCode[i],normSet,rangeName) ;
737 if (pdf->isSelectedComp()) {
738 value += val*_coefCache[i]/snormVal ;
739 }
740 }
741 }
742
743 return value ;
744}
745
746
747
748////////////////////////////////////////////////////////////////////////////////
749/// Return the number of expected events, which is either the sum of all coefficients
750/// or the sum of the components extended terms, multiplied with the fraction that
751/// is in the current range w.r.t the reference range
752
753double RooAddPdf::expectedEvents(const RooArgSet* nset) const
754{
755 double expectedTotal{0.0};
756
757 cxcoutD(Caching) << "RooAddPdf::expectedEvents(" << GetName() << ") calling getProjCache with nset = " << (nset?*nset:RooArgSet()) << std::endl ;
758 AddCacheElem& cache = *getProjCache(nset) ;
759 updateCoefficients(cache, nset);
760
761 if (cache.doProjection()) {
762
763 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
764 double ncomp = _allExtendable ? static_cast<RooAbsPdf&>(_pdfList[i]).expectedEvents(nset)
765 : static_cast<RooAbsReal&>(_coefList[i]).getVal(nset);
766 expectedTotal += cache.rangeProjScaleFactor(i) * ncomp ;
767
768 }
769
770 } else {
771
772 if (_allExtendable) {
773 for(auto const& arg : _pdfList) {
774 expectedTotal += static_cast<RooAbsPdf*>(arg)->expectedEvents(nset) ;
775 }
776 } else {
777 for(auto const& arg : _coefList) {
778 expectedTotal += static_cast<RooAbsReal*>(arg)->getVal(nset) ;
779 }
780 }
781
782 }
783 return expectedTotal ;
784}
785
786
787std::unique_ptr<RooAbsReal> RooAddPdf::createExpectedEventsFunc(const RooArgSet *nset) const
788{
789 std::unique_ptr<RooAbsReal> out;
790
791 auto name = std::string(GetName()) + "_expectedEvents";
792 if (_allExtendable) {
793 RooArgSet sumSet;
794 for (auto *pdf : static_range_cast<RooAbsPdf *>(_pdfList)) {
795 sumSet.addOwned(pdf->createExpectedEventsFunc(nset));
796 }
797 out = std::make_unique<RooAddition>(name.c_str(), name.c_str(), sumSet);
798 out->addOwnedComponents(std::move(sumSet));
799 } else {
800 out = std::make_unique<RooAddition>(name.c_str(), name.c_str(), _coefList);
801 }
802
803 RooArgList prodList;
804
805 // Make sure _refCoefNorm is defined
807
808 if (!_allExtendable) {
809 // If the _refCoefNorm is empty or it's equal to normSet anyway, this is not
810 // a conditional pdf and we don't need to do any transformation. See also
811 // RooAddPdf::compleForNormSet() for more explanations on a similar logic.
812 if (!_refCoefNorm.empty() && !nset->equals(_refCoefNorm)) {
813 prodList.addOwned(std::unique_ptr<RooAbsReal>{createIntegral(*nset, _refCoefNorm)});
814 }
815
816 // Optionally multiply with fractional normalization. I this case, we
817 // replace the original factor stored in "out".
818 if (!_normRange.IsNull()) {
819 std::unique_ptr<RooAbsReal> owner;
820 RooArgList terms;
821 // The integrals own each other in a chain. We do this because it's
822 // not possible to add two objects with the same name via
823 // addOwnedComponents(), and it happens in some user models that some
824 // component pdfs are the same. Hence, the integrals might share names
825 // too and we can't add them all in one go as owned objects of the
826 // final integral sum.
827 for (auto *pdf : static_range_cast<RooAbsPdf *>(_pdfList)) {
828 auto next = std::unique_ptr<RooAbsReal>{pdf->createIntegral(*nset, *nset, _normRange)};
829 terms.add(*next);
830 if (owner)
831 next->addOwnedComponents(std::move(owner));
832 owner = std::move(next);
833 }
834 auto fracIntegName = std::string(GetName()) + "_integSum";
835 auto fracInteg =
836 std::make_unique<RooRealSumFunc>(fracIntegName.c_str(), fracIntegName.c_str(), _coefList, terms);
837 fracInteg->addOwnedComponents(std::move(owner));
838
839 out = std::move(fracInteg);
840 }
841 }
842
843 std::string finalName = std::string(out->GetName()) + "_finalized";
844 if (prodList.empty()) {
845 // If there are no additional factors, just return the single factor we have
846 return out;
847 } else {
848 prodList.addOwned(std::move(out));
849 }
850 auto finalOut = std::make_unique<RooProduct>(finalName.c_str(), finalName.c_str(), prodList);
851 finalOut->addOwnedComponents(std::move(prodList));
852 return finalOut;
853}
854
855
856////////////////////////////////////////////////////////////////////////////////
857/// Interface function used by test statistics to freeze choice of observables
858/// for interpretation of fraction coefficients
859
860void RooAddPdf::selectNormalization(const RooArgSet* depSet, bool force)
861{
862 // Make sure _refCoefNorm is defined
864
865 if (!force && !_refCoefNorm.empty()) {
866 return ;
867 }
868
869 if (!depSet) {
871 return ;
872 }
873
874 fixCoefNormalization(*std::unique_ptr<RooArgSet>{getObservables(depSet)}) ;
875}
876
877
878
879////////////////////////////////////////////////////////////////////////////////
880/// Interface function used by test statistics to freeze choice of range
881/// for interpretation of fraction coefficients
882
883void RooAddPdf::selectNormalizationRange(const char* rangeName, bool force)
884{
885 if (!force && _refCoefRangeName) {
886 return ;
887 }
888
889 fixCoefRange(rangeName) ;
890}
891
892
893
894////////////////////////////////////////////////////////////////////////////////
895/// Return specialized context to efficiently generate toy events from RooAddPdfs
896/// return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ; // WVE DEBUG
897
899 const RooArgSet* auxProto, bool verbose) const
900{
901 return RooAddGenContext::create(*this,vars,prototype,auxProto,verbose).release();
902}
903
904
905
906////////////////////////////////////////////////////////////////////////////////
907/// Loop over components for plot sampling hints and merge them if there are multiple
908
909std::list<double>* RooAddPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
910{
911 return RooRealSumPdf::plotSamplingHint(_pdfList, obs, xlo, xhi);
912}
913
914
915////////////////////////////////////////////////////////////////////////////////
916/// Loop over components for plot sampling hints and merge them if there are multiple
917
918std::list<double>* RooAddPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
919{
920 return RooRealSumPdf::binBoundaries(_pdfList, obs, xlo, xhi);
921}
922
923
924////////////////////////////////////////////////////////////////////////////////
925/// If all components that depend on obs are binned, so is their sum.
927{
929}
930
931
932////////////////////////////////////////////////////////////////////////////////
933/// Label OK'ed components of a RooAddPdf with cache-and-track
934
936{
938}
939
940
941
942////////////////////////////////////////////////////////////////////////////////
943/// Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the
944/// product operator construction
945
946void RooAddPdf::printMetaArgs(std::ostream& os) const
947{
949}
950
951
952bool RooAddPdf::redirectServersHook(const RooAbsCollection & newServerList, bool mustReplaceAll,
953 bool nameChange, bool isRecursiveStep)
954{
955 // If a server is redirected, the cached normalization set might not point
956 // to the right observables anymore. We need to reset it.
957 _copyOfLastNormSet.reset();
958 return RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursiveStep);
959}
960
961
962std::unique_ptr<RooAbsArg>
964{
965 // Make sure _refCoefNorm is defined
967
968 auto newArg = std::unique_ptr<RooAbsReal>{static_cast<RooAbsReal *>(Clone())};
969 ctx.markAsCompiled(*newArg);
970
971 // In case conditional observables, e.g. p(x|y), the _refCoefNorm is set to
972 // all observables (x, y) and the normSet doesn't contain the conditional
973 // observables (so it only contains x in this example).
974
975 // If the _refCoefNorm is empty or it's equal to normSet anyway, this is not
976 // a conditional pdf and we don't need to do any transformation.
977 if(_refCoefNorm.empty() || normSet.equals(_refCoefNorm)) {
978 ctx.compileServers(*newArg, normSet);
979 return newArg;
980 }
981
982 // In the conditional case, things become more complicated. The original
983 // getValV() method is covering this case with very complicated logic,
984 // caching multiple new RooFit objects to scale the individual coefficients
985 // of the RooAddPdf.
986 //
987 // However, it's not complicated what we need to do mathematically:
988 //
989 // Since:
990 // 1. p(x, y) = p(x | y) * p(y)
991 // 2. p(y) = Integral of p(x, y) over x
992 //
993 // We conclude:
994 // p(x, y)
995 // p(x | y) = --------------------------
996 // Integral of p(x, y) over x
997 //
998 // What follows is the implementation of this formula in RooFit. By doing
999 // this here in compileForNormSet(), we don't invoke the old RooAddPdf
1000 // projection caches (note that no conditional pdfs are on the right hand
1001 // side of the equation).
1002 std::string finalName = std::string(GetName()) + "_conditional";
1003 std::unique_ptr<RooAbsReal> denom{newArg->createIntegral(normSet, _refCoefNorm)};
1004 auto finalArg = std::make_unique<RooGenericPdf>(finalName.c_str(), "@0/@1", RooArgList{*newArg, *denom});
1005 ctx.compileServers(*denom, _refCoefNorm);
1006 ctx.markAsCompiled(*denom);
1007 ctx.markAsCompiled(*finalArg);
1008 ctx.compileServers(*newArg, _refCoefNorm);
1009 finalArg->addOwnedComponents(std::move(newArg));
1010 finalArg->addOwnedComponents(std::move(denom));
1011 return finalArg;
1012}
#define cxcoutD(a)
#define coutW(a)
#define coutE(a)
#define TRACE_CREATE
Definition RooTrace.h:23
#define ClassImp(name)
Definition Rtypes.h:377
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
double rangeProjScaleFactor(std::size_t idx) const
bool doProjection() const
double suppNormVal(std::size_t idx) 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:79
void clearValueAndShapeDirty() const
Definition RooAbsArg.h:599
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
void 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:421
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:91
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 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.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * first() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
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
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:299
TString _normRange
Normalization range.
Definition RooAbsPdf.h:343
RooAbsReal * _norm
Definition RooAbsPdf.h:319
Int_t _errorCount
Number of errors remaining to print.
Definition RooAbsPdf.h:335
const char * normRange() const
Definition RooAbsPdf.h:250
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
static Int_t _verboseEval
Definition RooAbsPdf.h:314
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
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
virtual void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
double _value
Cache for current value of object.
Definition RooAbsReal.h:543
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...
static void updateCoefficients(RooAbsPdf const &addPdf, std::vector< double > &coefCache, RooArgList const &pdfList, bool haveLastCoef, AddCacheElem &cache, const RooArgSet *nset, RooArgSet const &refCoefNormSet, bool allExtendable, int &coefErrCount)
Update the RooAddPdf coefficients for a given normalization set and projection configuration.
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooListProxy _coefList
List of coefficients.
Definition RooAddPdf.h:132
bool _allExtendable
Flag indicating if all PDF components are extendable.
Definition RooAddPdf.h:136
RooAICRegistry _codeReg
! Registry of component analytical integration codes
Definition RooAddPdf.h:129
RooFit::UniqueId< RooArgSet >::Value_t _idOfLastUsedNormSet
!
Definition RooAddPdf.h:145
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral defined by given scenario code.
std::unique_ptr< const RooArgSet > _copyOfLastNormSet
!
Definition RooAddPdf.h:146
void updateCoefficients(AddCacheElem &cache, const RooArgSet *nset, bool syncCoefValues=true) const
Update the coefficient values in the given cache element: calculate new remainder fraction,...
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 ...
Int_t _coefErrCount
! Coefficient error counter
Definition RooAddPdf.h:139
bool _haveLastCoef
Flag indicating if last PDFs coefficient was supplied in the ctor.
Definition RooAddPdf.h:135
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.
bool _recursive
Flag indicating is fractions are treated recursively.
Definition RooAddPdf.h:137
RooObjCacheManager _projCacheMgr
Definition RooAddPdf.h:110
void materializeRefCoefNormFromAttribute() const
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute addition of PDFs in batches.
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:104
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
Manager of cache with coefficient projections and transformations.
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:131
TNamed * _refCoefRangeName
Reference range name for coefficient interpretation.
Definition RooAddPdf.h:105
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
The cache manager.
std::vector< double > _coefCache
! Transient cache with transformed values of coefficients
Definition RooAddPdf.h:107
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:55
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...
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
A class to maintain the context for squashing of RooFit models into code.
void markAsCompiled(RooAbsArg &arg) const
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition DataMap.cxx:40
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:22
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:37
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
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...
static void translateImpl(RooFit::Detail::CodeSquashContext &ctx, RooAbsArg const *klass, RooArgList const &funcList, RooArgList const &coefList)
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:47
Bool_t IsNull() const
Definition TString.h:418
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.
std::vector< std::span< const double > > VarVector
std::vector< double > ArgVector
void compute(Config cfg, Computer comp, RestrictArr output, size_t size, const VarVector &vars, ArgVector &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()