Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAddModel.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 RooAddModel
19///
20/// RooAddModel is an efficient implementation of a sum of PDFs of the form
21/// \f[
22/// c_1 \cdot \mathrm{PDF}_1 + c_2 \cdot \mathrm{PDF}_2 + ... + c_n \cdot \mathrm{PDF}_n
23/// \f]
24/// or
25/// \f[
26/// 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
27/// \f]
28/// The first form is for extended likelihood fits, where the
29/// expected number of events is \f$ \sum_i c_i \f$. The coefficients \f$ c_i \f$
30/// can either be explicitly provided, or, if all components support
31/// extended likelihood fits, they can be calculated from the contribution
32/// of each PDF to the total number of expected events.
33///
34/// In the second form, the sum of the coefficients is enforced to be one,
35/// and the coefficient of the last PDF is calculated from that condition.
36///
37/// RooAddModel relies on each component PDF to be normalized, and will perform
38/// no normalization other than calculating the proper last coefficient \f$ c_n \f$, if requested.
39/// An (enforced) condition for this assumption is that each \f$ \mathrm{PDF}_i \f$ is independent
40/// of each coefficient \f$ i \f$.
41///
42///
43
44#include "RooAddModel.h"
45
46#include "RooAddHelpers.h"
47#include "RooMsgService.h"
48#include "RooDataSet.h"
49#include "RooRealProxy.h"
50#include "RooPlot.h"
51#include "RooRealVar.h"
52#include "RooAddGenContext.h"
53#include "RooNameReg.h"
54#include "RooBatchCompute.h"
55
56using std::endl, std::ostream;
57
58
59
60////////////////////////////////////////////////////////////////////////////////
61
63 : _refCoefNorm("!refCoefNorm", "Reference coefficient normalization set", this, false, false),
64 _projCacheMgr(this, 10),
65 _intCacheMgr(this, 10),
66 _coefErrCount(_errorCount)
67{
68}
69
70
71
72////////////////////////////////////////////////////////////////////////////////
73/// Generic constructor from list of PDFs and list of coefficients.
74/// Each pdf list element (i) is paired with coefficient list element (i).
75/// The number of coefficients must be either equal to the number of PDFs,
76/// in which case extended MLL fitting is enabled, or be one less.
77///
78/// All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal.
79
80RooAddModel::RooAddModel(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, bool ownPdfList) :
82 _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,false,false),
83 _projCacheMgr(this,10),
84 _intCacheMgr(this,10),
85 _codeReg(10),
86 _pdfList("!pdfs","List of PDFs",this),
87 _coefList("!coefficients","List of coefficients",this)
88{
89 const std::string ownName(GetName() ? GetName() : "");
90 if (inPdfList.size() > inCoefList.size() + 1 || inPdfList.size() < inCoefList.size()) {
91 std::stringstream msgSs;
92 msgSs << "RooAddModel::RooAddModel(" << ownName
93 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
94 const std::string msgStr = msgSs.str();
95 coutE(InputArguments) << msgStr << "\n";
96 throw std::runtime_error(msgStr);
97 }
98
99 // Constructor with N PDFs and N or N-1 coefs
100 std::size_t i = 0;
101 for (auto const &coef : inCoefList) {
102 auto pdf = inPdfList.at(i);
103 if (!pdf) {
104 std::stringstream msgSs;
105 msgSs << "RooAddModel::RooAddModel(" << ownName
106 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
107 const std::string msgStr = msgSs.str();
108 coutE(InputArguments) << msgStr << "\n";
109 throw std::runtime_error(msgStr);
110 }
111 if (!coef) {
112 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName
113 << ") encountered and undefined coefficient, ignored\n";
114 continue;
115 }
116 if (!dynamic_cast<RooAbsReal *>(coef)) {
117 auto coefName = coef->GetName();
118 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName << ") coefficient "
119 << (coefName != nullptr ? coefName : "") << " is not of type RooAbsReal, ignored\n";
120 continue;
121 }
122 if (!dynamic_cast<RooAbsPdf *>(pdf)) {
123 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName << ") pdf "
124 << (pdf->GetName() ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored\n";
125 continue;
126 }
127 _pdfList.add(*pdf);
128 _coefList.add(*coef);
129 i++;
130 }
131
132 if (i < inPdfList.size()) {
133 auto pdf = inPdfList.at(i);
134 if (!dynamic_cast<RooAbsPdf *>(pdf)) {
135 std::stringstream msgSs;
136 msgSs << "RooAddModel::RooAddModel(" << ownName << ") last pdf " << (pdf->GetName() ? pdf->GetName() : "")
137 << " is not of type RooAbsPdf, fatal error";
138 const std::string msgStr = msgSs.str();
139 coutE(InputArguments) << msgStr << "\n";
140 throw std::runtime_error(msgStr);
141 }
142 _pdfList.add(*pdf);
143 } else {
144 _haveLastCoef = true;
145 }
146
148
149 if (ownPdfList) {
151 }
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// Copy constructor
156
159 _refCoefNorm("!refCoefNorm", this, other._refCoefNorm),
160 _refCoefRangeName((TNamed *)other._refCoefRangeName),
161 _projCacheMgr(other._projCacheMgr, this),
162 _intCacheMgr(other._intCacheMgr, this),
163 _codeReg(other._codeReg),
164 _pdfList("!pdfs", this, other._pdfList),
165 _coefList("!coefficients", this, other._coefList),
166 _haveLastCoef(other._haveLastCoef),
167 _allExtendable(other._allExtendable),
168 _coefErrCount(_errorCount)
169{
170}
171
172
173
174////////////////////////////////////////////////////////////////////////////////
175/// By default the interpretation of the fraction coefficients is
176/// performed in the contextual choice of observables. This makes the
177/// shape of the p.d.f explicitly dependent on the choice of
178/// observables. This method instructs RooAddModel to freeze the
179/// interpretation of the coefficients to be done in the given set of
180/// observables. If frozen, fractions are automatically transformed
181/// from the reference normalization set to the contextual normalization
182/// set by ratios of integrals
183
185{
186 if (refCoefNorm.empty()) {
187 return ;
188 }
189
192
194}
195
196
197
198////////////////////////////////////////////////////////////////////////////////
199/// By default the interpretation of the fraction coefficients is
200/// performed in the default range. This make the shape of a RooAddModel
201/// explicitly dependent on the range of the observables. To allow
202/// a range independent definition of the fraction this function
203/// instructs RooAddModel to freeze its interpretation in the given
204/// named range. If the current normalization range is different
205/// from the reference range, the appropriate fraction coefficients
206/// are automatically calculated from the reference fractions using
207/// ratios of integrals.
208
213
214
215
216////////////////////////////////////////////////////////////////////////////////
217/// Instantiate a clone of this resolution model representing a convolution with given
218/// basis function. The owners object name is incorporated in the clones name
219/// to avoid multiple convolution objects with the same name in complex PDF structures.
220///
221/// RooAddModel will clone all the component models to create a composite convolution object
222
224{
225 // Check that primary variable of basis functions is our convolution variable
226 if (inBasis->getParameter(0) != x.absArg()) {
227 coutE(InputArguments) << "RooAddModel::convolution(" << GetName()
228 << ") convolution parameter of basis function and PDF don't match" << std::endl ;
229 ccoutE(InputArguments) << "basis->findServer(0) = " << inBasis->findServer(0) << " " << inBasis->findServer(0)->GetName() << std::endl ;
230 ccoutE(InputArguments) << "x.absArg() = " << x.absArg() << " " << x.absArg()->GetName() << std::endl ;
231 inBasis->Print("v") ;
232 return nullptr ;
233 }
234
236 newName.Append("_conv_") ;
237 newName.Append(inBasis->GetName()) ;
238 newName.Append("_[") ;
239 newName.Append(owner->GetName()) ;
240 newName.Append("]") ;
241
243 newTitle.Append(" convoluted with basis function ") ;
244 newTitle.Append(inBasis->GetName()) ;
245
248 // Create component convolution
249 RooResolutionModel* conv = model->convolution(inBasis,owner) ;
250 modelList.add(*conv) ;
251 }
252
254 for (auto coef : _coefList) {
255 theCoefList.add(*coef) ;
256 }
257
259 for (std::set<std::string>::const_iterator attrIt = _boolAttrib.begin();
260 attrIt != _boolAttrib.end(); ++attrIt) {
261 convSum->setAttribute((*attrIt).c_str()) ;
262 }
263 for (std::map<std::string,std::string>::const_iterator attrIt = _stringAttrib.begin();
264 attrIt != _stringAttrib.end(); ++attrIt) {
265 convSum->setStringAttribute((attrIt->first).c_str(), (attrIt->second).c_str()) ;
266 }
267 convSum->changeBasis(inBasis) ;
268 return convSum ;
269}
270
271
272
273////////////////////////////////////////////////////////////////////////////////
274/// Return code for basis function representing by 'name' string.
275/// The basis code of the first component model will be returned,
276/// if the basis is supported by all components. Otherwise 0
277/// is returned
278
280{
281 bool first(true);
282 bool code(false);
284 Int_t subCode = model->basisCode(name) ;
285 if (first) {
286 code = subCode ;
287 first = false ;
288 } else if (subCode==0) {
289 code = false ;
290 }
291 }
292
293 return code ;
294}
295
296
297
298////////////////////////////////////////////////////////////////////////////////
299/// Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated over iset
300/// in range 'rangeName'. If cache element does not exist, create and fill it on the fly. The cache contains
301/// suplemental normalization terms (in case not all added p.d.f.s have the same observables), projection
302/// integrals to calculated transformed fraction coefficients when a frozen reference frame is provided
303/// and projection integrals for similar transformations when a frozen reference range is provided.
304
305AddCacheElem* RooAddModel::getProjCache(const RooArgSet* nset, const RooArgSet* iset) const
306{
307 // Check if cache already exists
308 auto cache = static_cast<AddCacheElem*>(_projCacheMgr.getObj(nset,iset,nullptr,normRange()));
309 if (cache) {
310 return cache ;
311 }
312
313 //Create new cache
314 cache = new AddCacheElem{*this, _pdfList, _coefList, nset, iset, _refCoefNorm,
317
319
320 return cache;
321}
322
323
324////////////////////////////////////////////////////////////////////////////////
325/// Update the coefficient values in the given cache element: calculate new remainder
326/// fraction, normalize fractions obtained from extended ML terms to unity, and
327/// multiply the various range and dimensional corrections needed in the
328/// current use context.
329
330void RooAddModel::updateCoefficients(AddCacheElem &cache, const RooArgSet *nset) const
331{
332 _coefCache.resize(_pdfList.size());
333 for (std::size_t i = 0; i < _coefList.size(); ++i) {
334 _coefCache[i] = static_cast<RooAbsReal const &>(_coefList[i]).getVal(nset);
335 }
336 if (_allExtendable) {
337 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
338 auto &pdf = static_cast<RooAbsPdf &>(_pdfList[i]);
339 _coefCache[i] = pdf.expectedEvents(!_refCoefNorm.empty() ? &_refCoefNorm : nset);
340 }
341 }
342
343 RooAddHelpers::updateCoefficients(*this, _pdfList.size(), _coefCache, _haveLastCoef || _allExtendable, cache,
345}
346
347
348////////////////////////////////////////////////////////////////////////////////
349/// Calculate the current value
350
352{
353 const RooArgSet* nset = _normSet ;
354 AddCacheElem* cache = getProjCache(nset) ;
355
356 updateCoefficients(*cache,nset) ;
357
358
359 // Do running sum of coef/pdf pairs, calculate lastCoef.
360 double snormVal ;
361 double value(0) ;
362 Int_t i(0) ;
363 for (auto *pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
364
365 if (_coefCache[i]!=0.) {
366 snormVal = nset ? cache->suppNormVal(i) : 1.0 ;
367 double pdfVal = pdf->getVal(nset) ;
368 // double pdfNorm = pdf->getNorm(nset) ;
369 if (pdf->isSelectedComp()) {
371 cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
372 << pdf->GetName() << "] " << pdfVal << " * " << _coefCache[i] << " / " << snormVal << std::endl ;
373 }
374 }
375 i++ ;
376 }
377
378 return value ;
379}
380
382{
383 // Like many other functions in this class, the implementation was copy-pasted from the RooAddPdf
384 RooBatchCompute::Config config = ctx.config(this);
385
386 _coefCache.resize(_pdfList.size());
387 for (std::size_t i = 0; i < _coefList.size(); ++i) {
388 auto coefVals = ctx.at(&_coefList[i]);
389 // We don't support per-event coefficients in this function. If the CPU
390 // mode is used, we can just fall back to the RooAbsReal implementation.
391 // With CUDA, we can't do that because the inputs might be on the device.
392 // That's why we throw an exception then.
393 if (coefVals.size() > 1) {
394 if (config.useCuda()) {
395 throw std::runtime_error("The RooAddPdf doesn't support per-event coefficients in CUDA mode yet!");
396 }
398 return;
399 }
400 _coefCache[i] = coefVals[0];
401 }
402
403 std::vector<std::span<const double>> pdfs;
404 std::vector<double> coefs;
405 AddCacheElem *cache = getProjCache(nullptr);
406 updateCoefficients(*cache, nullptr);
407
408 for (unsigned int pdfNo = 0; pdfNo < _pdfList.size(); ++pdfNo) {
409 auto pdf = static_cast<RooAbsPdf *>(&_pdfList[pdfNo]);
410 if (pdf->isSelectedComp()) {
411 pdfs.push_back(ctx.at(pdf));
412 coefs.push_back(_coefCache[pdfNo] / cache->suppNormVal(pdfNo));
413 }
414 }
416}
417
418
419////////////////////////////////////////////////////////////////////////////////
420/// Reset error counter to given value, limiting the number
421/// of future error messages for this pdf to 'resetValue'
422
428
429
430
431////////////////////////////////////////////////////////////////////////////////
432/// Check if PDF is valid for given normalization set.
433/// Coefficient and PDF must be non-overlapping, but pdf-coefficient
434/// pairs may overlap each other
435
437{
438 bool ret(false) ;
439
440 for (unsigned int i = 0; i < _coefList.size(); ++i) {
441 auto pdf = &_pdfList[i];
442 auto coef = &_coefList[i];
443
444 if (pdf->observableOverlaps(nset,*coef)) {
445 coutE(InputArguments) << "RooAddModel::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName()
446 << " and PDF " << pdf->GetName() << " have one or more dependents in common" << std::endl ;
447 ret = true ;
448 }
449 }
450
451 return ret ;
452}
453
454
455
456////////////////////////////////////////////////////////////////////////////////
457
459 const RooArgSet* normSet, const char* rangeName) const
460{
461 if (_forceNumInt) return 0 ;
462
463 // Declare that we can analytically integrate all requested observables
464 analVars.add(allVars) ;
465
466 // Retrieve (or create) the required component integral list
467 Int_t code ;
469 getCompIntList(normSet,&allVars,cilist,code,rangeName) ;
470
471 return code+1 ;
472
473}
474
475
476
477////////////////////////////////////////////////////////////////////////////////
478/// Check if this configuration was created before
479
481{
482 Int_t sterileIdx(-1) ;
483
485 if (cache) {
486 code = _intCacheMgr.lastIndex() ;
487 compIntList = &cache->_intList ;
488
489 return ;
490 }
491
492 // Create containers for partial integral components to be generated
493 cache = new IntCacheElem ;
494
495 // Fill Cache
497
498 cache->_intList.addOwned(std::unique_ptr<RooAbsReal>{model->createIntegral(*iset,nset,nullptr,isetRangeName)});
499 }
500
501 // Store the partial integral list and return the assigned code ;
503
504 // Fill references to be returned
505 compIntList = &cache->_intList ;
506}
507
508
509
510////////////////////////////////////////////////////////////////////////////////
511/// Return analytical integral defined by given scenario code
512
514{
515 // No integration scenario
516 if (code==0) {
517 return getVal(normSet) ;
518 }
519
520 // Partial integration scenarios
521 IntCacheElem* cache = static_cast<IntCacheElem*>(_intCacheMgr.getObjByIndex(code-1)) ;
522
524
525 // If cache has been sterilized, revive this slot
526 if (cache==nullptr) {
527 std::unique_ptr<RooArgSet> vars{getParameters(RooArgSet())} ;
528 RooArgSet nset = _intCacheMgr.selectFromSet1(*vars, code-1) ;
529 RooArgSet iset = _intCacheMgr.selectFromSet2(*vars, code-1) ;
530
531 int code2 = -1 ;
533 } else {
534
535 compIntList = &cache->_intList ;
536
537 }
538
539 // Calculate the current value
540 const RooArgSet* nset = _normSet ;
542
544
545 // Do running sum of coef/pdf pairs, calculate lastCoef.
546 double snormVal ;
547 double value(0) ;
548 Int_t i(0) ;
550 if (_coefCache[i]!=0.) {
551 snormVal = nset ? pcache->suppNormVal(i) : 1.0 ;
552 double intVal = pdfInt->getVal(nset) ;
554 cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
555 << pdfInt->GetName() << "] " << intVal << " * " << _coefCache[i] << " / " << snormVal << std::endl ;
556 }
557 i++ ;
558 }
559
560 return value ;
561
562}
563
564
565
566////////////////////////////////////////////////////////////////////////////////
567/// Return the number of expected events, which is either the sum of all coefficients
568/// or the sum of the components extended terms
569
570double RooAddModel::expectedEvents(const RooArgSet* nset) const
571{
572 double expectedTotal(0.0);
573
574 if (_allExtendable) {
575
576 // Sum of the extended terms
577 for (auto *pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
578 expectedTotal += pdf->expectedEvents(nset) ;
579 }
580
581 } else {
582
583 // Sum the coefficients
584 for (auto *coef : static_range_cast<RooAbsReal*>(_coefList)) {
585 expectedTotal += coef->getVal() ;
586 }
587 }
588
589 return expectedTotal;
590}
591
592
593
594////////////////////////////////////////////////////////////////////////////////
595/// Interface function used by test statistics to freeze choice of observables
596/// for interpretation of fraction coefficients
597
599{
600 if (!force && !_refCoefNorm.empty()) {
601 return ;
602 }
603
604 if (!depSet) {
606 return ;
607 }
608
612}
613
614
615
616////////////////////////////////////////////////////////////////////////////////
617/// Interface function used by test statistics to freeze choice of range
618/// for interpretation of fraction coefficients
619
621{
622 if (!force && _refCoefRangeName) {
623 return ;
624 }
625
627}
628
629
630
631////////////////////////////////////////////////////////////////////////////////
632/// Return specialized context to efficiently generate toy events from RooAddModels.
633
635 const RooArgSet* auxProto, bool verbose) const
636{
637 return RooAddGenContext::create(*this,vars,prototype,auxProto,verbose).release();
638}
639
640
641
642////////////////////////////////////////////////////////////////////////////////
643/// Direct generation is safe if all components say so
644
646{
647 for (auto *pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
648
649 if (!pdf->isDirectGenSafe(arg)) {
650 return false ;
651 }
652 }
653 return true ;
654}
655
656
657
658////////////////////////////////////////////////////////////////////////////////
659/// Return pseud-code that indicates if all components can do internal generation (1) or not (0)
660
661Int_t RooAddModel::getGenerator(const RooArgSet& directVars, RooArgSet &/*generateVars*/, bool /*staticInitOK*/) const
662{
663 for (auto *pdf : static_range_cast<RooAbsPdf*>(_pdfList)) {
664
665 RooArgSet tmp ;
666 if (pdf->getGenerator(directVars,tmp)==0) {
667 return 0 ;
668 }
669 }
670 return 1 ;
671}
672
673
674
675
676////////////////////////////////////////////////////////////////////////////////
677/// This function should never be called as RooAddModel implements a custom generator context
678
680{
681 assert(0) ;
682}
683
684
685////////////////////////////////////////////////////////////////////////////////
686/// List all RooAbsArg derived contents in this cache element
687
693
694
695////////////////////////////////////////////////////////////////////////////////
696/// Customized printing of arguments of a RooAddModel to more intuitively reflect the contents of the
697/// product operator construction
698
699void RooAddModel::printMetaArgs(ostream& os) const
700{
701 bool first(true) ;
702
703 os << "(" ;
704 for (unsigned int i=0; i < _coefList.size(); ++i) {
705 auto coef = &_coefList[i];
706 auto pdf = &_pdfList[i];
707 if (!first) {
708 os << " + " ;
709 } else {
710 first = false ;
711 }
712 os << coef->GetName() << " * " << pdf->GetName() ;
713 }
714 if (_pdfList.size() > _coefList.size()) {
715 os << " + [%] * " << _pdfList[_pdfList.size()-1].GetName() ;
716 }
717 os << ") " ;
718}
719
#define ccoutE(a)
#define cxcoutD(a)
#define coutE(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:148
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
std::map< std::string, std::string > _stringAttrib
Definition RooAbsArg.h:589
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
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.
std::set< std::string > _boolAttrib
Definition RooAbsArg.h:588
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
const char * GetName() const override
Returns name of object.
Storage_t::size_type size() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
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...
RooArgSet const * _normSet
! Normalization set with for above integral
Definition RooAbsPdf.h:314
Int_t _errorCount
Number of errors remaining to print.
Definition RooAbsPdf.h:328
const char * normRange() const
Definition RooAbsPdf.h:246
static Int_t _verboseEval
Definition RooAbsPdf.h:308
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:542
friend class AddCacheElem
Definition RooAbsReal.h:407
virtual void doEval(RooFit::EvalContext &) const
Base function for computing multiple values of a RooAbsReal.
static std::unique_ptr< RooAbsGenContext > create(const Pdf_t &pdf, const RooArgSet &vars, const RooDataSet *prototype, const RooArgSet *auxProto, bool verbose)
Returns a RooAddGenContext if possible, or, if the RooAddGenContext doesn't support this particular R...
RooArgList containedArgs(Action) override
List all RooAbsArg derived contents in this cache element.
RooArgList _intList
List of component integrals.
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:27
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 RooAddModels.
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooAddModel to more intuitively reflect the contents of the pro...
RooObjCacheManager _projCacheMgr
! Manager of cache with coefficient projections and transformations
void getCompIntList(const RooArgSet *nset, const RooArgSet *iset, pRooArgList &compIntList, Int_t &code, const char *isetRangeName) const
Check if this configuration was created before.
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...
RooSetProxy _refCoefNorm
! Reference observable set for coefficient interpretation
Definition RooAddModel.h:96
bool _allExtendable
Flag indicating if all PDF components are extendable.
Int_t _coefErrCount
! Coefficient error counter
RooArgSet _ownedComps
! Owned components
RooListProxy _coefList
List of coefficients.
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 basisCode(const char *name) const override
Return code for basis function representing by 'name' string.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral defined by given scenario code.
RooResolutionModel * convolution(RooFormulaVar *basis, RooAbsArg *owner) const override
Instantiate a clone of this resolution model representing a convolution with given basis function.
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
bool checkObservables(const RooArgSet *nset) const override
Check if PDF is valid for given normalization set.
void generateEvent(Int_t code) override
This function should never be called as RooAddModel implements a custom generator context.
bool _haveLastCoef
Flag indicating if last PDFs coefficient was supplied in the constructor.
double expectedEvents(const RooArgSet *nset) const override
Return expected number of events for extended likelihood calculation, which is the sum of all coeffic...
RooListProxy _pdfList
List of component PDFs.
RooObjCacheManager _intCacheMgr
! Manager of cache with integrals
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...
void fixCoefNormalization(const RooArgSet &refCoefNorm)
By default the interpretation of the fraction coefficients is performed in the contextual choice of o...
void fixCoefRange(const char *rangeName)
By default the interpretation of the fraction coefficients is performed in the default range.
TNamed * _refCoefRangeName
! Reference range name for coefficient interpretation
Definition RooAddModel.h:97
std::vector< double > _coefCache
! Transient cache with transformed values of coefficients
Definition RooAddModel.h:99
double evaluate() const override
Calculate the current value.
void updateCoefficients(AddCacheElem &cache, const RooArgSet *nset) const
Update the coefficient values in the given cache element: calculate new remainder fraction,...
bool isDirectGenSafe(const RooAbsArg &arg) const override
Direct generation is safe if all components say so.
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Return pseud-code that indicates if all components can do internal generation (1) or not (0)
AddCacheElem * getProjCache(const RooArgSet *nset, const RooArgSet *iset=nullptr) const
Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated...
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * absArg() const
Return pointer to contained argument.
Definition RooArgProxy.h:46
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Minimal configuration struct to steer the evaluation of a single node with the RooBatchCompute librar...
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
RooArgSet selectFromSet1(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 1 with a given index and an i...
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 2 with a given index and an i...
void reset()
Clear the cache.
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Container class to hold unbinned data.
Definition RooDataSet.h:32
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:39
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
virtual RooResolutionModel * convolution(RooFormulaVar *basis, RooAbsArg *owner) const
Instantiate a clone of this resolution model representing a convolution with given basis function.
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
Basic string class.
Definition TString.h:138
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})