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
59
60
61////////////////////////////////////////////////////////////////////////////////
62
64 : _refCoefNorm("!refCoefNorm", "Reference coefficient normalization set", this, false, false),
65 _projCacheMgr(this, 10),
66 _intCacheMgr(this, 10),
67 _coefErrCount(_errorCount)
68{
69}
70
71
72
73////////////////////////////////////////////////////////////////////////////////
74/// Generic constructor from list of PDFs and list of coefficients.
75/// Each pdf list element (i) is paired with coefficient list element (i).
76/// The number of coefficients must be either equal to the number of PDFs,
77/// in which case extended MLL fitting is enabled, or be one less.
78///
79/// All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal.
80
81RooAddModel::RooAddModel(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, bool ownPdfList) :
82 RooResolutionModel(name,title,(static_cast<RooResolutionModel*>(inPdfList.at(0)))->convVar()),
83 _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,false,false),
84 _projCacheMgr(this,10),
85 _intCacheMgr(this,10),
86 _codeReg(10),
87 _pdfList("!pdfs","List of PDFs",this),
88 _coefList("!coefficients","List of coefficients",this)
89{
90 const std::string ownName(GetName() ? GetName() : "");
91 if (inPdfList.size() > inCoefList.size() + 1 || inPdfList.size() < inCoefList.size()) {
92 std::stringstream msgSs;
93 msgSs << "RooAddModel::RooAddModel(" << ownName
94 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
95 const std::string msgStr = msgSs.str();
96 coutE(InputArguments) << msgStr << "\n";
97 throw std::runtime_error(msgStr);
98 }
99
100 // Constructor with N PDFs and N or N-1 coefs
101 std::size_t i = 0;
102 for (auto const &coef : inCoefList) {
103 auto pdf = inPdfList.at(i);
104 if (!pdf) {
105 std::stringstream msgSs;
106 msgSs << "RooAddModel::RooAddModel(" << ownName
107 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1";
108 const std::string msgStr = msgSs.str();
109 coutE(InputArguments) << msgStr << "\n";
110 throw std::runtime_error(msgStr);
111 }
112 if (!coef) {
113 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName
114 << ") encountered and undefined coefficient, ignored\n";
115 continue;
116 }
117 if (!dynamic_cast<RooAbsReal *>(coef)) {
118 auto coefName = coef->GetName();
119 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName << ") coefficient "
120 << (coefName != nullptr ? coefName : "") << " is not of type RooAbsReal, ignored\n";
121 continue;
122 }
123 if (!dynamic_cast<RooAbsPdf *>(pdf)) {
124 coutE(InputArguments) << "RooAddModel::RooAddModel(" << ownName << ") pdf "
125 << (pdf->GetName() ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored\n";
126 continue;
127 }
128 _pdfList.add(*pdf);
129 _coefList.add(*coef);
130 i++;
131 }
132
133 if (i < inPdfList.size()) {
134 auto pdf = inPdfList.at(i);
135 if (!dynamic_cast<RooAbsPdf *>(pdf)) {
136 std::stringstream msgSs;
137 msgSs << "RooAddModel::RooAddModel(" << ownName << ") last pdf " << (pdf->GetName() ? pdf->GetName() : "")
138 << " is not of type RooAbsPdf, fatal error";
139 const std::string msgStr = msgSs.str();
140 coutE(InputArguments) << msgStr << "\n";
141 throw std::runtime_error(msgStr);
142 }
143 _pdfList.add(*pdf);
144 } else {
145 _haveLastCoef = true;
146 }
147
149
150 if (ownPdfList) {
152 }
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// Copy constructor
157
158RooAddModel::RooAddModel(const RooAddModel &other, const char *name)
159 : RooResolutionModel(other, name),
160 _refCoefNorm("!refCoefNorm", this, other._refCoefNorm),
161 _refCoefRangeName((TNamed *)other._refCoefRangeName),
162 _projCacheMgr(other._projCacheMgr, this),
163 _intCacheMgr(other._intCacheMgr, this),
164 _codeReg(other._codeReg),
165 _pdfList("!pdfs", this, other._pdfList),
166 _coefList("!coefficients", this, other._coefList),
167 _haveLastCoef(other._haveLastCoef),
168 _allExtendable(other._allExtendable),
169 _coefErrCount(_errorCount)
170{
171}
172
173
174
175////////////////////////////////////////////////////////////////////////////////
176/// By default the interpretation of the fraction coefficients is
177/// performed in the contextual choice of observables. This makes the
178/// shape of the p.d.f explicitly dependent on the choice of
179/// observables. This method instructs RooAddModel to freeze the
180/// interpretation of the coefficients to be done in the given set of
181/// observables. If frozen, fractions are automatically transformed
182/// from the reference normalization set to the contextual normalization
183/// set by ratios of integrals
184
186{
187 if (refCoefNorm.empty()) {
188 return ;
189 }
190
192 _refCoefNorm.add(refCoefNorm) ;
193
195}
196
197
198
199////////////////////////////////////////////////////////////////////////////////
200/// By default the interpretation of the fraction coefficients is
201/// performed in the default range. This make the shape of a RooAddModel
202/// explicitly dependent on the range of the observables. To allow
203/// a range independent definition of the fraction this function
204/// instructs RooAddModel to freeze its interpretation in the given
205/// named range. If the current normalization range is different
206/// from the reference range, the appropriate fraction coefficients
207/// are automatically calculated from the reference fractions using
208/// ratios of integrals.
209
210void RooAddModel::fixCoefRange(const char* rangeName)
211{
212 _refCoefRangeName = const_cast<TNamed*>(RooNameReg::ptr(rangeName));
213}
214
215
216
217////////////////////////////////////////////////////////////////////////////////
218/// Instantiate a clone of this resolution model representing a convolution with given
219/// basis function. The owners object name is incorporated in the clones name
220/// to avoid multiple convolution objects with the same name in complex PDF structures.
221///
222/// RooAddModel will clone all the component models to create a composite convolution object
223
225{
226 // Check that primary variable of basis functions is our convolution variable
227 if (inBasis->getParameter(0) != x.absArg()) {
228 coutE(InputArguments) << "RooAddModel::convolution(" << GetName()
229 << ") convolution parameter of basis function and PDF don't match" << endl ;
230 ccoutE(InputArguments) << "basis->findServer(0) = " << inBasis->findServer(0) << " " << inBasis->findServer(0)->GetName() << endl ;
231 ccoutE(InputArguments) << "x.absArg() = " << x.absArg() << " " << x.absArg()->GetName() << endl ;
232 inBasis->Print("v") ;
233 return nullptr ;
234 }
235
236 TString newName(GetName()) ;
237 newName.Append("_conv_") ;
238 newName.Append(inBasis->GetName()) ;
239 newName.Append("_[") ;
240 newName.Append(owner->GetName()) ;
241 newName.Append("]") ;
242
243 TString newTitle(GetTitle()) ;
244 newTitle.Append(" convoluted with basis function ") ;
245 newTitle.Append(inBasis->GetName()) ;
246
247 RooArgList modelList ;
248 for (auto obj : _pdfList) {
249 auto model = static_cast<RooResolutionModel*>(obj);
250 // Create component convolution
251 RooResolutionModel* conv = model->convolution(inBasis,owner) ;
252 modelList.add(*conv) ;
253 }
254
255 RooArgList theCoefList ;
256 for (auto coef : _coefList) {
257 theCoefList.add(*coef) ;
258 }
259
260 RooAddModel* convSum = new RooAddModel(newName,newTitle,modelList,theCoefList,true) ;
261 for (std::set<std::string>::const_iterator attrIt = _boolAttrib.begin();
262 attrIt != _boolAttrib.end(); ++attrIt) {
263 convSum->setAttribute((*attrIt).c_str()) ;
264 }
265 for (std::map<std::string,std::string>::const_iterator attrIt = _stringAttrib.begin();
266 attrIt != _stringAttrib.end(); ++attrIt) {
267 convSum->setStringAttribute((attrIt->first).c_str(), (attrIt->second).c_str()) ;
268 }
269 convSum->changeBasis(inBasis) ;
270 return convSum ;
271}
272
273
274
275////////////////////////////////////////////////////////////////////////////////
276/// Return code for basis function representing by 'name' string.
277/// The basis code of the first component model will be returned,
278/// if the basis is supported by all components. Otherwise 0
279/// is returned
280
282{
283 bool first(true);
284 bool code(false);
285 for (auto obj : _pdfList) {
286 auto model = static_cast<RooResolutionModel*>(obj);
287 Int_t subCode = model->basisCode(name) ;
288 if (first) {
289 code = subCode ;
290 first = false ;
291 } else if (subCode==0) {
292 code = false ;
293 }
294 }
295
296 return code ;
297}
298
299
300
301////////////////////////////////////////////////////////////////////////////////
302/// Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated over iset
303/// in range 'rangeName'. If cache element does not exist, create and fill it on the fly. The cache contains
304/// suplemental normalization terms (in case not all added p.d.f.s have the same observables), projection
305/// integrals to calculated transformed fraction coefficients when a frozen reference frame is provided
306/// and projection integrals for similar transformations when a frozen reference range is provided.
307
309{
310 // Check if cache already exists
311 auto cache = static_cast<AddCacheElem*>(_projCacheMgr.getObj(nset,iset,nullptr,normRange()));
312 if (cache) {
313 return cache ;
314 }
315
316 //Create new cache
317 cache = new AddCacheElem{*this, _pdfList, _coefList, nset, iset, _refCoefNorm,
320
321 _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(normRange())) ;
322
323 return cache;
324}
325
326
327////////////////////////////////////////////////////////////////////////////////
328/// Update the coefficient values in the given cache element: calculate new remainder
329/// fraction, normalize fractions obtained from extended ML terms to unity, and
330/// multiply the various range and dimensional corrections needed in the
331/// current use context.
332
334{
335 _coefCache.resize(_pdfList.size());
336 for(std::size_t i = 0; i < _coefList.size(); ++i) {
337 _coefCache[i] = static_cast<RooAbsReal const&>(_coefList[i]).getVal(nset);
338 }
341}
342
343
344////////////////////////////////////////////////////////////////////////////////
345/// Calculate the current value
346
348{
349 const RooArgSet* nset = _normSet ;
350 AddCacheElem* cache = getProjCache(nset) ;
351
352 updateCoefficients(*cache,nset) ;
353
354
355 // Do running sum of coef/pdf pairs, calculate lastCoef.
356 double snormVal ;
357 double value(0) ;
358 Int_t i(0) ;
359 for (auto obj : _pdfList) {
360 auto pdf = static_cast<RooAbsPdf*>(obj);
361
362 if (_coefCache[i]!=0.) {
363 snormVal = nset ? cache->suppNormVal(i) : 1.0 ;
364 double pdfVal = pdf->getVal(nset) ;
365 // double pdfNorm = pdf->getNorm(nset) ;
366 if (pdf->isSelectedComp()) {
367 value += pdfVal*_coefCache[i]/snormVal ;
368 cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
369 << pdf->GetName() << "] " << pdfVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
370 }
371 }
372 i++ ;
373 }
374
375 return value ;
376}
377
379{
380 // Like many other functions in this class, the implementation was copy-pasted from the RooAddPdf
381 RooBatchCompute::Config config = ctx.config(this);
382
383 _coefCache.resize(_pdfList.size());
384 for (std::size_t i = 0; i < _coefList.size(); ++i) {
385 auto coefVals = ctx.at(&_coefList[i]);
386 // We don't support per-event coefficients in this function. If the CPU
387 // mode is used, we can just fall back to the RooAbsReal implementation.
388 // With CUDA, we can't do that because the inputs might be on the device.
389 // That's why we throw an exception then.
390 if (coefVals.size() > 1) {
391 if (config.useCuda()) {
392 throw std::runtime_error("The RooAddPdf doesn't support per-event coefficients in CUDA mode yet!");
393 }
395 return;
396 }
397 _coefCache[i] = coefVals[0];
398 }
399
400 std::vector<std::span<const double>> pdfs;
401 std::vector<double> coefs;
402 AddCacheElem *cache = getProjCache(nullptr);
403 updateCoefficients(*cache, nullptr);
404
405 for (unsigned int pdfNo = 0; pdfNo < _pdfList.size(); ++pdfNo) {
406 auto pdf = static_cast<RooAbsPdf *>(&_pdfList[pdfNo]);
407 if (pdf->isSelectedComp()) {
408 pdfs.push_back(ctx.at(pdf));
409 coefs.push_back(_coefCache[pdfNo] / cache->suppNormVal(pdfNo));
410 }
411 }
412 RooBatchCompute::compute(config, RooBatchCompute::AddPdf, ctx.output(), pdfs, coefs);
413}
414
415
416////////////////////////////////////////////////////////////////////////////////
417/// Reset error counter to given value, limiting the number
418/// of future error messages for this pdf to 'resetValue'
419
421{
423 _coefErrCount = resetValue ;
424}
425
426
427
428////////////////////////////////////////////////////////////////////////////////
429/// Check if PDF is valid for given normalization set.
430/// Coefficient and PDF must be non-overlapping, but pdf-coefficient
431/// pairs may overlap each other
432
434{
435 bool ret(false) ;
436
437 for (unsigned int i = 0; i < _coefList.size(); ++i) {
438 auto pdf = &_pdfList[i];
439 auto coef = &_coefList[i];
440
441 if (pdf->observableOverlaps(nset,*coef)) {
442 coutE(InputArguments) << "RooAddModel::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName()
443 << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
444 ret = true ;
445 }
446 }
447
448 return ret ;
449}
450
451
452
453////////////////////////////////////////////////////////////////////////////////
454
456 const RooArgSet* normSet, const char* rangeName) const
457{
458 if (_forceNumInt) return 0 ;
459
460 // Declare that we can analytically integrate all requested observables
461 analVars.add(allVars) ;
462
463 // Retrieve (or create) the required component integral list
464 Int_t code ;
465 RooArgList *cilist ;
466 getCompIntList(normSet,&allVars,cilist,code,rangeName) ;
467
468 return code+1 ;
469
470}
471
472
473
474////////////////////////////////////////////////////////////////////////////////
475/// Check if this configuration was created before
476
477void RooAddModel::getCompIntList(const RooArgSet* nset, const RooArgSet* iset, pRooArgList& compIntList, Int_t& code, const char* isetRangeName) const
478{
479 Int_t sterileIdx(-1) ;
480
481 IntCacheElem* cache = static_cast<IntCacheElem*>(_intCacheMgr.getObj(nset,iset,&sterileIdx,RooNameReg::ptr(isetRangeName))) ;
482 if (cache) {
483 code = _intCacheMgr.lastIndex() ;
484 compIntList = &cache->_intList ;
485
486 return ;
487 }
488
489 // Create containers for partial integral components to be generated
490 cache = new IntCacheElem ;
491
492 // Fill Cache
493 for (auto obj : _pdfList) {
494 auto model = static_cast<RooResolutionModel*>(obj);
495
496 cache->_intList.addOwned(std::unique_ptr<RooAbsReal>{model->createIntegral(*iset,nset,nullptr,isetRangeName)});
497 }
498
499 // Store the partial integral list and return the assigned code ;
500 code = _intCacheMgr.setObj(nset,iset,(RooAbsCacheElement*)cache,RooNameReg::ptr(isetRangeName)) ;
501
502 // Fill references to be returned
503 compIntList = &cache->_intList ;
504}
505
506
507
508////////////////////////////////////////////////////////////////////////////////
509/// Return analytical integral defined by given scenario code
510
511double RooAddModel::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
512{
513 // No integration scenario
514 if (code==0) {
515 return getVal(normSet) ;
516 }
517
518 // Partial integration scenarios
519 IntCacheElem* cache = static_cast<IntCacheElem*>(_intCacheMgr.getObjByIndex(code-1)) ;
520
521 RooArgList* compIntList ;
522
523 // If cache has been sterilized, revive this slot
524 if (cache==nullptr) {
525 std::unique_ptr<RooArgSet> vars{getParameters(RooArgSet())} ;
526 RooArgSet nset = _intCacheMgr.selectFromSet1(*vars, code-1) ;
527 RooArgSet iset = _intCacheMgr.selectFromSet2(*vars, code-1) ;
528
529 int code2 = -1 ;
530 getCompIntList(&nset,&iset,compIntList,code2,rangeName) ;
531 } else {
532
533 compIntList = &cache->_intList ;
534
535 }
536
537 // Calculate the current value
538 const RooArgSet* nset = _normSet ;
539 AddCacheElem* pcache = getProjCache(nset) ;
540
541 updateCoefficients(*pcache,nset) ;
542
543 // Do running sum of coef/pdf pairs, calculate lastCoef.
544 double snormVal ;
545 double value(0) ;
546 Int_t i(0) ;
547 for (const auto obj : *compIntList) {
548 auto pdfInt = static_cast<const RooAbsReal*>(obj);
549 if (_coefCache[i]!=0.) {
550 snormVal = nset ? pcache->suppNormVal(i) : 1.0 ;
551 double intVal = pdfInt->getVal(nset) ;
552 value += intVal*_coefCache[i]/snormVal ;
553 cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
554 << pdfInt->GetName() << "] " << intVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
555 }
556 i++ ;
557 }
558
559 return value ;
560
561}
562
563
564
565////////////////////////////////////////////////////////////////////////////////
566/// Return the number of expected events, which is either the sum of all coefficients
567/// or the sum of the components extended terms
568
569double RooAddModel::expectedEvents(const RooArgSet* nset) const
570{
571 double expectedTotal(0.0);
572
573 if (_allExtendable) {
574
575 // Sum of the extended terms
576 for (auto obj : _pdfList) {
577 auto pdf = static_cast<RooAbsPdf*>(obj);
578 expectedTotal += pdf->expectedEvents(nset) ;
579 }
580
581 } else {
582
583 // Sum the coefficients
584 for (const auto obj : _coefList) {
585 auto coef = static_cast<RooAbsReal*>(obj);
586 expectedTotal += coef->getVal() ;
587 }
588 }
589
590 return expectedTotal;
591}
592
593
594
595////////////////////////////////////////////////////////////////////////////////
596/// Interface function used by test statistics to freeze choice of observables
597/// for interpretation of fraction coefficients
598
599void RooAddModel::selectNormalization(const RooArgSet* depSet, bool force)
600{
601 if (!force && !_refCoefNorm.empty()) {
602 return ;
603 }
604
605 if (!depSet) {
607 return ;
608 }
609
610 RooArgSet myDepSet;
611 getObservables(depSet, myDepSet);
612 fixCoefNormalization(myDepSet);
613}
614
615
616
617////////////////////////////////////////////////////////////////////////////////
618/// Interface function used by test statistics to freeze choice of range
619/// for interpretation of fraction coefficients
620
621void RooAddModel::selectNormalizationRange(const char* rangeName, bool force)
622{
623 if (!force && _refCoefRangeName) {
624 return ;
625 }
626
627 fixCoefRange(rangeName) ;
628}
629
630
631
632////////////////////////////////////////////////////////////////////////////////
633/// Return specialized context to efficiently generate toy events from RooAddModels.
634
636 const RooArgSet* auxProto, bool verbose) const
637{
638 return RooAddGenContext::create(*this,vars,prototype,auxProto,verbose).release();
639}
640
641
642
643////////////////////////////////////////////////////////////////////////////////
644/// Direct generation is safe if all components say so
645
647{
648 for (auto obj : _pdfList) {
649 auto pdf = static_cast<RooAbsPdf*>(obj);
650
651 if (!pdf->isDirectGenSafe(arg)) {
652 return false ;
653 }
654 }
655 return true ;
656}
657
658
659
660////////////////////////////////////////////////////////////////////////////////
661/// Return pseud-code that indicates if all components can do internal generation (1) or not (0)
662
663Int_t RooAddModel::getGenerator(const RooArgSet& directVars, RooArgSet &/*generateVars*/, bool /*staticInitOK*/) const
664{
665 for (auto obj : _pdfList) {
666 auto pdf = static_cast<RooAbsPdf*>(obj);
667
668 RooArgSet tmp ;
669 if (pdf->getGenerator(directVars,tmp)==0) {
670 return 0 ;
671 }
672 }
673 return 1 ;
674}
675
676
677
678
679////////////////////////////////////////////////////////////////////////////////
680/// This function should never be called as RooAddModel implements a custom generator context
681
683{
684 assert(0) ;
685}
686
687
688////////////////////////////////////////////////////////////////////////////////
689/// List all RooAbsArg derived contents in this cache element
690
692{
693 RooArgList allNodes(_intList) ;
694 return allNodes ;
695}
696
697
698////////////////////////////////////////////////////////////////////////////////
699/// Customized printing of arguments of a RooAddModel to more intuitively reflect the contents of the
700/// product operator construction
701
702void RooAddModel::printMetaArgs(ostream& os) const
703{
704 bool first(true) ;
705
706 os << "(" ;
707 for (unsigned int i=0; i < _coefList.size(); ++i) {
708 auto coef = &_coefList[i];
709 auto pdf = &_pdfList[i];
710 if (!first) {
711 os << " + " ;
712 } else {
713 first = false ;
714 }
715 os << coef->GetName() << " * " << pdf->GetName() ;
716 }
717 if (_pdfList.size() > _coefList.size()) {
718 os << " + [%] * " << _pdfList[_pdfList.size()-1].GetName() ;
719 }
720 os << ") " ;
721}
722
#define ccoutE(a)
#define cxcoutD(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:382
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 suppNormVal(std::size_t idx) const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:263
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
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:600
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
std::map< std::string, std::string > _stringAttrib
Definition RooAbsArg.h:601
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition RooAbsArg.h:153
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
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
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:40
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 integral (owned by _normMgr)
Definition RooAbsPdf.h:320
Int_t _errorCount
Number of errors remaining to print.
Definition RooAbsPdf.h:334
const char * normRange() const
Definition RooAbsPdf.h:250
static Int_t _verboseEval
Definition RooAbsPdf.h:314
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:539
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...
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.
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 * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
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:33
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...
RooAbsArg * getParameter(const char *name) const
Return pointer to parameter with given name.
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 void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
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:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:572
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})