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