Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsAnaConvPdf.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 RooAbsAnaConvPdf
19/// \ingroup Roofitcore
20///
21/// RooAbsAnaConvPdf is the base class for PDFs that represent a
22/// physics model that can be analytically convolved with a resolution model.
23///
24/// To achieve factorization between the physics model and the resolution
25/// model, each physics model must be able to be written in the form
26/// \f[
27/// \mathrm{Phys}(x, \bar{a}, \bar{b}) = \sum_k \mathrm{coef}_k(\bar{a}) * \mathrm{basis}_k(x,\bar{b})
28/// \f]
29///
30/// where \f$ \mathrm{basis}_k \f$ are a limited number of functions in terms of the variable
31/// to be convoluted, and \f$ \mathrm{coef}_k \f$ are coefficients independent of the convolution
32/// variable.
33///
34/// Classes derived from RooResolutionModel implement
35/// \f[
36/// R_k(x,\bar{b},\bar{c}) = \int \mathrm{basis}_k(x', \bar{b}) \cdot \mathrm{resModel}(x-x',\bar{c}) \; \mathrm{d}x',
37/// \f]
38///
39/// which RooAbsAnaConvPdf uses to construct the pdf for [ Phys (x) R ] :
40/// \f[
41/// \mathrm{PDF}(x,\bar{a},\bar{b},\bar{c}) = \sum_k \mathrm{coef}_k(\bar{a}) * R_k(x,\bar{b},\bar{c})
42/// \f]
43///
44/// A minimal implementation of a RooAbsAnaConvPdf physics model consists of
45///
46/// - A constructor that declares the required basis functions using the declareBasis() method.
47/// The declareBasis() function assigns a unique identifier code to each declare basis
48///
49/// - An implementation of `coefficient(Int_t code)` returning the coefficient value for each
50/// declared basis function
51///
52/// Optionally, analytical integrals can be provided for the coefficient functions. The
53/// interface for this is quite similar to that for integrals of regular PDFs. Two functions,
54/// \code{.cpp}
55/// Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
56/// Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const
57/// \endcode
58///
59/// advertise the coefficient integration capabilities and implement them respectively.
60/// Please see RooAbsPdf for additional details. Advertised analytical integrals must be
61/// valid for all coefficients.
62
63#include "RooAbsAnaConvPdf.h"
64
65#include "RooFit.h"
66#include "RooMsgService.h"
67#include "Riostream.h"
68#include "RooResolutionModel.h"
69#include "RooRealVar.h"
70#include "RooFormulaVar.h"
71#include "RooConvGenContext.h"
72#include "RooGenContext.h"
73#include "RooTruthModel.h"
74#include "RooConvCoefVar.h"
75#include "RooNameReg.h"
76
77using namespace std;
78
80
81
82////////////////////////////////////////////////////////////////////////////////
83/// Default constructor, required for persistence
84
86 _isCopy(false),
87 _coefNormMgr(this,10)
88{
89}
90
91
92
93////////////////////////////////////////////////////////////////////////////////
94/// Constructor. The supplied resolution model must be constructed with the same
95/// convoluted variable as this physics model ('convVar')
96
97RooAbsAnaConvPdf::RooAbsAnaConvPdf(const char *name, const char *title,
98 const RooResolutionModel& model, RooRealVar& cVar) :
99 RooAbsPdf(name,title), _isCopy(kFALSE),
100 _model("!model","Original resolution model",this,(RooResolutionModel&)model,kFALSE,kFALSE),
101 _convVar("!convVar","Convolution variable",this,cVar,kFALSE,kFALSE),
102 _convSet("!convSet","Set of resModel X basisFunc convolutions",this),
103 _coefNormMgr(this,10),
104 _codeReg(10)
105{
106 _model.absArg()->setAttribute("NOCacheAndTrack") ;
107}
108
109
110
111////////////////////////////////////////////////////////////////////////////////
112
114 RooAbsPdf(other,name), _isCopy(kTRUE),
115 _model("!model",this,other._model),
116 _convVar("!convVar",this,other._convVar),
117 _convSet("!convSet",this,other._convSet),
118 _coefNormMgr(other._coefNormMgr,this),
119 _codeReg(other._codeReg)
120{
121 // Copy constructor
122 if (_model.absArg()) {
123 _model.absArg()->setAttribute("NOCacheAndTrack") ;
124 }
126}
127
128
129
130////////////////////////////////////////////////////////////////////////////////
131/// Destructor
132
134{
135 if (!_isCopy) {
136 std::vector<RooAbsArg*> tmp(_convSet.begin(), _convSet.end());
137
138 for (auto arg : tmp) {
139 _convSet.remove(*arg) ;
140 delete arg ;
141 }
142 }
143
144}
145
146
147////////////////////////////////////////////////////////////////////////////////
148/// Declare a basis function for use in this physics model. The string expression
149/// must be a valid RooFormulVar expression representing the basis function, referring
150/// to the convolution variable as '@0', and any additional parameters (supplied in
151/// 'params' as '@1','@2' etc.
152///
153/// The return value is a unique identifier code, that will be passed to coefficient()
154/// to identify the basis function for which the coefficient is requested. If the
155/// resolution model used does not support the declared basis function, code -1 is
156/// returned.
157///
158
159Int_t RooAbsAnaConvPdf::declareBasis(const char* expression, const RooArgList& params)
160{
161 // Sanity check
162 if (_isCopy) {
163 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): ERROR attempt to "
164 << " declare basis functions in a copied RooAbsAnaConvPdf" << endl ;
165 return -1 ;
166 }
167
168 // Resolution model must support declared basis
169 if (!((RooResolutionModel*)_model.absArg())->isBasisSupported(expression)) {
170 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): resolution model "
171 << _model.absArg()->GetName()
172 << " doesn't support basis function " << expression << endl ;
173 return -1 ;
174 }
175
176 // Instantiate basis function
177 RooArgList basisArgs(_convVar.arg()) ;
178 basisArgs.add(params) ;
179
180 TString basisName(expression) ;
181 for (const auto arg : basisArgs) {
182 basisName.Append("_") ;
183 basisName.Append(arg->GetName()) ;
184 }
185
186 RooFormulaVar* basisFunc = new RooFormulaVar(basisName, expression, basisArgs);
187 basisFunc->setAttribute("RooWorkspace::Recycle") ;
188 basisFunc->setAttribute("NOCacheAndTrack") ;
189 basisFunc->setOperMode(operMode()) ;
190 _basisList.addOwned(*basisFunc) ;
191
192 // Instantiate resModel x basisFunc convolution
193 RooAbsReal* conv = ((RooResolutionModel*)_model.absArg())->convolution(basisFunc,this) ;
194 if (!conv) {
195 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): unable to construct convolution with basis function '"
196 << expression << "'" << endl ;
197 return -1 ;
198 }
199 _convSet.add(*conv) ;
200
201 return _convSet.index(conv) ;
202}
203
204
205
206////////////////////////////////////////////////////////////////////////////////
207/// Change the current resolution model to newModel
208
210{
211 RooArgList newConvSet ;
212 Bool_t allOK(kTRUE) ;
213 for (auto convArg : _convSet) {
214 auto conv = static_cast<RooResolutionModel*>(convArg);
215
216 // Build new resolution model
217 RooResolutionModel* newConv = newModel.convolution((RooFormulaVar*)&conv->basis(),this) ;
218 if (!newConvSet.add(*newConv)) {
219 allOK = kFALSE ;
220 break ;
221 }
222 }
223
224 // Check if all convolutions were successfully built
225 if (!allOK) {
226 // Delete new basis functions created sofar
227 std::for_each(newConvSet.begin(), newConvSet.end(), [](RooAbsArg* arg){delete arg;});
228
229 return kTRUE ;
230 }
231
232 // Replace old convolutions with new set
234 _convSet.addOwned(newConvSet) ;
235
236 // Update server link by hand, since _model.setArg() below will not do this
238
239 _model.setArg((RooResolutionModel&)newModel) ;
240 return kFALSE ;
241}
242
243
244
245
246////////////////////////////////////////////////////////////////////////////////
247/// Create a generator context for this p.d.f. If both the p.d.f and the resolution model
248/// support internal generation of the convolution observable on an infinite domain,
249/// deploy a specialized convolution generator context, which generates the physics distribution
250/// and the smearing separately, adding them a posteriori. If this is not possible return
251/// a (slower) generic generation context that uses accept/reject sampling
252
254 const RooArgSet* auxProto, Bool_t verbose) const
255{
256 // Check if the resolution model specifies a special context to be used.
257 RooResolutionModel* conv = dynamic_cast<RooResolutionModel*>(_model.absArg());
258 assert(conv);
259
260 RooArgSet* modelDep = _model.absArg()->getObservables(&vars) ;
261 modelDep->remove(*convVar(),kTRUE,kTRUE) ;
262 Int_t numAddDep = modelDep->getSize() ;
263 delete modelDep ;
264
265 // Check if physics PDF and resolution model can both directly generate the convolution variable
266 RooArgSet dummy ;
267 Bool_t pdfCanDir = (getGenerator(*convVar(),dummy) != 0) ;
268 Bool_t resCanDir = conv && (conv->getGenerator(*convVar(),dummy)!=0) && conv->isDirectGenSafe(*convVar()) ;
269
270 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
271 // Any resolution model with more dependents than the convolution variable
272 // or pdf or resmodel do not support direct generation
273 string reason ;
274 if (numAddDep>0) reason += "Resolution model has more observables than the convolution variable. " ;
275 if (!pdfCanDir) reason += "PDF does not support internal generation of convolution observable. " ;
276 if (!resCanDir) reason += "Resolution model does not support internal generation of convolution observable. " ;
277
278 coutI(Generation) << "RooAbsAnaConvPdf::genContext(" << GetName() << ") Using regular accept/reject generator for convolution p.d.f because: " << reason.c_str() << endl ;
279 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
280 }
281
282 RooAbsGenContext* context = conv->modelGenContext(*this, vars, prototype, auxProto, verbose);
283 if (context) return context;
284
285 // Any other resolution model: use specialized generator context
286 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
287}
288
289
290
291////////////////////////////////////////////////////////////////////////////////
292/// Return true if it is safe to generate the convolution observable
293/// from the internal generator (this is the case if the chosen resolution
294/// model is the truth model)
295
297{
298
299 // All direct generation of convolution arg if model is truth model
300 if (!TString(_convVar.absArg()->GetName()).CompareTo(arg.GetName()) &&
301 dynamic_cast<RooTruthModel*>(_model.absArg())) {
302 return kTRUE ;
303 }
304
305 return RooAbsPdf::isDirectGenSafe(arg) ;
306}
307
308
309
310////////////////////////////////////////////////////////////////////////////////
311/// Return a pointer to the convolution variable instance used in the resolution model
312
314{
316 if (!conv) return 0 ;
317 return &conv->convVar() ;
318}
319
320
321
322////////////////////////////////////////////////////////////////////////////////
323/// Calculate the current unnormalized value of the PDF
324///
325/// PDF = sum_k coef_k * [ basis_k (x) ResModel ]
326///
327
329{
330 Double_t result(0) ;
331
332 Int_t index(0) ;
333 for (auto convArg : _convSet) {
334 auto conv = static_cast<RooAbsPdf*>(convArg);
335 Double_t coef = coefficient(index++) ;
336 if (coef!=0.) {
337 Double_t c = conv->getVal(0) ;
338 Double_t r = coef ;
339 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") val += coef*conv [" << index-1 << "/"
340 << _convSet.getSize() << "] coef = " << r << " conv = " << c << endl ;
341 result += conv->getVal(0)*coef ;
342 } else {
343 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") [" << index-1 << "/" << _convSet.getSize() << "] coef = 0" << endl ;
344 }
345 }
346
347 return result ;
348}
349
350
351
352////////////////////////////////////////////////////////////////////////////////
353/// Advertise capability to perform (analytical) integrals
354/// internally. For a given integration request over allVars while
355/// normalized over normSet2 and in range 'rangeName', returns
356/// largest subset that can be performed internally in analVars
357/// Return code is unique integer code identifying integration scenario
358/// to be passed to analyticalIntegralWN() to calculate requeste integral
359///
360/// Class RooAbsAnaConv defers analytical integration request to
361/// resolution model and/or coefficient implementations and
362/// aggregates results into composite configuration with a unique
363/// code assigned by RooAICRegistry
364
366 RooArgSet& analVars, const RooArgSet* normSet2, const char* /*rangeName*/) const
367{
368 // Handle trivial no-integration scenario
369 if (allVars.getSize()==0) return 0 ;
370
371 if (_forceNumInt) return 0 ;
372
373 // Select subset of allVars that are actual dependents
374 RooArgSet* allDeps = getObservables(allVars) ;
375 RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
376
377 RooAbsArg *arg ;
378 RooResolutionModel *conv ;
379
380 RooArgSet* intSetAll = new RooArgSet(*allDeps,"intSetAll") ;
381
382 // Split intSetAll in coef/conv parts
383 RooArgSet* intCoefSet = new RooArgSet("intCoefSet") ;
384 RooArgSet* intConvSet = new RooArgSet("intConvSet") ;
385 TIterator* varIter = intSetAll->createIterator() ;
386 TIterator* convIter = _convSet.createIterator() ;
387
388 while(((arg=(RooAbsArg*) varIter->Next()))) {
389 Bool_t ok(kTRUE) ;
390 convIter->Reset() ;
391 while(((conv=(RooResolutionModel*) convIter->Next()))) {
392 if (conv->dependsOn(*arg)) ok=kFALSE ;
393 }
394
395 if (ok) {
396 intCoefSet->add(*arg) ;
397 } else {
398 intConvSet->add(*arg) ;
399 }
400
401 }
402 delete varIter ;
403
404
405 // Split normSetAll in coef/conv parts
406 RooArgSet* normCoefSet = new RooArgSet("normCoefSet") ;
407 RooArgSet* normConvSet = new RooArgSet("normConvSet") ;
408 RooArgSet* normSetAll = normSet ? (new RooArgSet(*normSet,"normSetAll")) : 0 ;
409 if (normSetAll) {
410 varIter = normSetAll->createIterator() ;
411 while(((arg=(RooAbsArg*) varIter->Next()))) {
412 Bool_t ok(kTRUE) ;
413 convIter->Reset() ;
414 while(((conv=(RooResolutionModel*) convIter->Next()))) {
415 if (conv->dependsOn(*arg)) ok=kFALSE ;
416 }
417
418 if (ok) {
419 normCoefSet->add(*arg) ;
420 } else {
421 normConvSet->add(*arg) ;
422 }
423
424 }
425 delete varIter ;
426 }
427 delete convIter ;
428
429 if (intCoefSet->getSize()==0) {
430 delete intCoefSet ; intCoefSet=0 ;
431 }
432 if (intConvSet->getSize()==0) {
433 delete intConvSet ; intConvSet=0 ;
434 }
435 if (normCoefSet->getSize()==0) {
436 delete normCoefSet ; normCoefSet=0 ;
437 }
438 if (normConvSet->getSize()==0) {
439 delete normConvSet ; normConvSet=0 ;
440 }
441
442
443
444 // Store integration configuration in registry
445 Int_t masterCode(0) ;
446 std::vector<Int_t> tmp(1, 0) ;
447
448 masterCode = _codeReg.store(tmp, intCoefSet, intConvSet, normCoefSet, normConvSet) + 1 ; // takes ownership of all sets
449
450 analVars.add(*allDeps) ;
451 delete allDeps ;
452 if (normSet) delete normSet ;
453 if (normSetAll) delete normSetAll ;
454 delete intSetAll ;
455
456// cout << this << "---> masterCode = " << masterCode << endl ;
457
458 return masterCode ;
459}
460
461
462
463
464////////////////////////////////////////////////////////////////////////////////
465/// Return analytical integral defined by given code, which is returned
466/// by getAnalyticalIntegralWN()
467///
468/// For unnormalized integrals the returned value is
469/// \f[
470/// \mathrm{PDF} = \sum_k \int \mathrm{coef}_k \; \mathrm{d}\bar{x}
471/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}\bar{y},
472/// \f]
473/// where \f$ \bar{x} \f$ is the set of coefficient dependents to be integrated,
474/// and \f$ \bar{y} \f$ the set of basis function dependents to be integrated.
475///
476/// For normalized integrals this becomes
477/// \f[
478/// \mathrm{PDF} = \frac{\sum_k \int \mathrm{coef}_k \; \mathrm{d}x
479/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}y}
480/// {\sum_k \int \mathrm{coef}_k \; \mathrm{d}v
481/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}w},
482/// \f]
483/// where
484/// * \f$ x \f$ is the set of coefficient dependents to be integrated,
485/// * \f$ y \f$ the set of basis function dependents to be integrated,
486/// * \f$ v \f$ is the set of coefficient dependents over which is normalized and
487/// * \f$ w \f$ is the set of basis function dependents over which is normalized.
488///
489/// Set \f$ x \f$ must be contained in \f$ v \f$ and set \f$ y \f$ must be contained in \f$ w \f$.
490///
491
492Double_t RooAbsAnaConvPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
493{
494 // WVE needs adaptation to handle new rangeName feature
495
496 // Handle trivial passthrough scenario
497 if (code==0) return getVal(normSet) ;
498
499 // Unpack master code
500 RooArgSet *intCoefSet, *intConvSet, *normCoefSet, *normConvSet ;
501 _codeReg.retrieve(code-1,intCoefSet,intConvSet,normCoefSet,normConvSet) ;
502
503 Int_t index(0) ;
504 Double_t answer(0) ;
505
506 if (normCoefSet==0&&normConvSet==0) {
507
508 // Integral over unnormalized function
509 Double_t integral(0) ;
510 const TNamed *_rangeName = RooNameReg::ptr(rangeName);
511 for (auto convArg : _convSet) {
512 auto conv = static_cast<RooResolutionModel*>(convArg);
513 Double_t coef = getCoefNorm(index++,intCoefSet,_rangeName) ;
514 //cout << "coefInt[" << index << "] = " << coef << " " ; intCoefSet->Print("1") ;
515 if (coef!=0) {
516 integral += coef* conv->getNormObj(0,intConvSet,_rangeName)->getVal();
517 cxcoutD(Eval) << "RooAbsAnaConv::aiWN(" << GetName() << ") [" << index-1 << "] integral += " << conv->getNorm(intConvSet) << endl ;
518 }
519
520 }
521 answer = integral ;
522
523 } else {
524
525 // Integral over normalized function
526 Double_t integral(0) ;
527 Double_t norm(0) ;
528 const TNamed *_rangeName = RooNameReg::ptr(rangeName);
529 for (auto convArg : _convSet) {
530 auto conv = static_cast<RooResolutionModel*>(convArg);
531
532 Double_t coefInt = getCoefNorm(index,intCoefSet,_rangeName) ;
533 //cout << "coefInt[" << index << "] = " << coefInt << "*" << term << " " << (intCoefSet?*intCoefSet:RooArgSet()) << endl ;
534 if (coefInt!=0) {
535 Double_t term = conv->getNormObj(0,intConvSet,_rangeName)->getVal();
536 integral += coefInt*term ;
537 }
538
539 Double_t coefNorm = getCoefNorm(index,normCoefSet) ;
540 //cout << "coefNorm[" << index << "] = " << coefNorm << "*" << term << " " << (normCoefSet?*normCoefSet:RooArgSet()) << endl ;
541 if (coefNorm!=0) {
542 Double_t term = conv->getNormObj(0,normConvSet)->getVal();
543 norm += coefNorm*term ;
544 }
545
546 index++ ;
547 }
548 answer = integral/norm ;
549 }
550
551 return answer ;
552}
553
554
555
556////////////////////////////////////////////////////////////////////////////////
557/// Default implementation of function advertising integration capabilities. The interface is
558/// similar to that of getAnalyticalIntegral except that an integer code is added that
559/// designates the coefficient number for which the integration capabilities are requested
560///
561/// This default implementation advertises that no internal integrals are supported.
562
563Int_t RooAbsAnaConvPdf::getCoefAnalyticalIntegral(Int_t /* coef*/, RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
564{
565 return 0 ;
566}
567
568
569
570////////////////////////////////////////////////////////////////////////////////
571/// Default implementation of function implementing advertised integrals. Only
572/// the pass-through scenario (no integration) is implemented.
573
574Double_t RooAbsAnaConvPdf::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* /*rangeName*/) const
575{
576 if (code==0) return coefficient(coef) ;
577 coutE(InputArguments) << "RooAbsAnaConvPdf::coefAnalyticalIntegral(" << GetName() << ") ERROR: unrecognized integration code: " << code << endl ;
578 assert(0) ;
579 return 1 ;
580}
581
582
583
584////////////////////////////////////////////////////////////////////////////////
585/// This function forces RooRealIntegral to offer all integration dependents
586/// to RooAbsAnaConvPdf::getAnalyticalIntegralWN() for consideration for
587/// internal integration, if RooRealIntegral considers this to be unsafe (e.g. due
588/// to hidden Jacobian terms).
589///
590/// RooAbsAnaConvPdf will not attempt to actually integrate all these dependents
591/// but feed them to the resolution models integration interface, which will
592/// make the final determination on how to integrate these dependents.
593
595{
596 return kTRUE ;
597}
598
599
600
601////////////////////////////////////////////////////////////////////////////////
602/// Returns the normalization integral value of the coefficient with number coefIdx over normalization
603/// set nset in range rangeName
604
605Double_t RooAbsAnaConvPdf::getCoefNorm(Int_t coefIdx, const RooArgSet* nset, const TNamed* rangeName) const
606{
607 if (nset==0) return coefficient(coefIdx) ;
608
609 CacheElem* cache = (CacheElem*) _coefNormMgr.getObj(nset,0,0,rangeName) ;
610 if (!cache) {
611
612 cache = new CacheElem ;
613
614 // Make list of coefficient normalizations
615 Int_t i ;
617
618 for (i=0 ; i<cache->_coefVarList.getSize() ; i++) {
619 RooAbsReal* coefInt = static_cast<RooAbsReal&>(*cache->_coefVarList.at(i)).createIntegral(*nset,RooNameReg::str(rangeName)) ;
620 cache->_normList.addOwned(*coefInt) ;
621 }
622
623 _coefNormMgr.setObj(nset,0,cache,rangeName) ;
624 }
625
626 return ((RooAbsReal*)cache->_normList.at(coefIdx))->getVal() ;
627}
628
629
630
631////////////////////////////////////////////////////////////////////////////////
632/// Build complete list of coefficient variables
633
635{
636 // Instantate a coefficient variables
637 for (Int_t i=0 ; i<_convSet.getSize() ; i++) {
638 RooArgSet* cvars = coefVars(i) ;
639 RooAbsReal* coefVar = new RooConvCoefVar(Form("%s_coefVar_%d",GetName(),i),"coefVar",*this,i,cvars) ;
640 varList.addOwned(*coefVar) ;
641 delete cvars ;
642 }
643
644}
645
646
647////////////////////////////////////////////////////////////////////////////////
648/// Return set of parameters with are used exclusively by the coefficient functions
649
651{
652 RooArgSet* cVars = getParameters((RooArgSet*)0) ;
653 std::vector<RooAbsArg*> tmp;
654 for (auto arg : *cVars) {
655 for (auto convSetArg : _convSet) {
656 if (convSetArg->dependsOn(*arg)) {
657 tmp.push_back(arg);
658 }
659 }
660 }
661
662 cVars->remove(tmp.begin(), tmp.end(), true, true);
663
664 return cVars ;
665}
666
667
668
669
670////////////////////////////////////////////////////////////////////////////////
671/// Print info about this object to the specified stream. In addition to the info
672/// from RooAbsPdf::printStream() we add:
673///
674/// Verbose : detailed information on convolution integrals
675
676void RooAbsAnaConvPdf::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
677{
678 RooAbsPdf::printMultiline(os,contents,verbose,indent);
679
680 os << indent << "--- RooAbsAnaConvPdf ---" << endl;
681 TIter iter = _convSet.createIterator() ;
682 RooResolutionModel* conv ;
683 while (((conv=(RooResolutionModel*)iter.Next()))) {
684 conv->printMultiline(os,contents,verbose,indent) ;
685 }
686}
687
688
689///////////////////////////////////////////////////////////////////////////////
690/// Label OK'ed components with cache-and-track
692{
693 RooFIter citer = _convSet.fwdIterator() ;
694 RooAbsArg* carg ;
695 while ((carg=citer.next())) {
696 if (carg->canNodeBeCached()==Always) {
697 trackNodes.add(*carg) ;
698 //cout << "tracking node RooAddPdf component " << carg->IsA()->GetName() << "::" << carg->GetName() << endl ;
699 }
700 }
701}
ROOT::R::TRInterface & r
Definition Object.C:4
#define c(i)
Definition RSha256.hxx:101
#define coutI(a)
#define cxcoutD(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
const std::vector< Int_t > & retrieve(Int_t masterCode) const
Retrieve the array of integer codes associated with the given master code.
Int_t store(const std::vector< Int_t > &codeList, RooArgSet *set1=0, RooArgSet *set2=0, RooArgSet *set3=0, RooArgSet *set4=0)
Store given arrays of integer codes, and up to four RooArgSets in the registry (each setX pointer may...
List of created basis functions.
RooAbsAnaConvPdf is the base class for PDFs that represent a physics model that can be analytically c...
virtual Bool_t forceAnalyticalInt(const RooAbsArg &dep) const
This function forces RooRealIntegral to offer all integration dependents to RooAbsAnaConvPdf::getAnal...
Double_t getCoefNorm(Int_t coefIdx, const RooArgSet &nset, const char *rangeName) const
friend class RooConvGenContext
virtual Bool_t changeModel(const RooResolutionModel &newModel)
Change the current resolution model to newModel.
virtual Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise capability to perform (analytical) integrals internally.
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
virtual void printMultiline(std::ostream &stream, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print info about this object to the specified stream.
RooObjCacheManager _coefNormMgr
void makeCoefVarList(RooArgList &) const
Build complete list of coefficient variables.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Create a generator context for this p.d.f.
RooAICRegistry _codeReg
Coefficient normalization manager.
virtual ~RooAbsAnaConvPdf()
Destructor.
virtual Double_t evaluate() const
Calculate the current unnormalized value of the PDF.
virtual void setCacheAndTrackHints(RooArgSet &)
Label OK'ed components with cache-and-track.
RooAbsRealLValue * convVar()
Retrieve the convolution variable.
virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Return analytical integral defined by given code, which is returned by getAnalyticalIntegralWN()
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
RooRealProxy _model
virtual RooArgSet * coefVars(Int_t coefIdx) const
Return set of parameters with are used exclusively by the coefficient functions.
RooRealProxy _convVar
RooAbsAnaConvPdf()
Default constructor, required for persistence.
virtual Double_t coefficient(Int_t basisIndex) const =0
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
RooListProxy _convSet
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Return true if it is safe to generate the convolution observable from the internal generator (this is...
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
virtual CacheMode canNodeBeCached() const
Definition RooAbsArg.h:427
friend class RooArgSet
Definition RooAbsArg.h:642
void setAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
void replaceServer(RooAbsArg &oldServer, RooAbsArg &newServer, Bool_t valueProp, Bool_t shapeProp)
Replace 'oldServer' with 'newServer'.
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...
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:499
Int_t getSize() const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
RooFIter fwdIterator() const
One-time forward iterator.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
const_iterator end() const
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add an argument and transfer the ownership to the collection.
const_iterator begin() const
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print multi line detailed information of this RooAbsPdf.
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
Bool_t _forceNumInt
Definition RooAbsReal.h:487
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
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
Definition RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooConvCoefVar is an auxilary class that represents the coefficient of a RooAbsAnaConvPdf implementat...
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:36
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::addOwned()
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
virtual void removeAll() override
Reimplementation of standard RooArgList::removeAll()
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE) override
Reimplementation of standard RooArgList::remove()
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
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.
virtual RooAbsGenContext * modelGenContext(const RooAbsAnaConvPdf &, const RooArgSet &, const RooDataSet *, const RooArgSet *, Bool_t) const
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Print info about this object to the specified stream.
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
bool setArg(T &newRef)
Change object held in proxy into newRef.
const T & arg() const
Return reference to object held in proxy.
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
TObject * Next()
Iterator abstract base class.
Definition TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:136
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
Definition TString.cxx:442
TString & Append(const char *cs)
Definition TString.h:564