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