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/// 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 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
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 std::endl, std::string, std::ostream;
78
79
80
81////////////////////////////////////////////////////////////////////////////////
82/// Default constructor, required for persistence
83
85 _isCopy(false),
86 _coefNormMgr(this,10)
87{
88}
89
90
91
92////////////////////////////////////////////////////////////////////////////////
93/// Constructor. The supplied resolution model must be constructed with the same
94/// convoluted variable as this physics model ('convVar')
95
96RooAbsAnaConvPdf::RooAbsAnaConvPdf(const char *name, const char *title,
97 const RooResolutionModel& model, RooRealVar& cVar) :
98 RooAbsPdf(name,title), _isCopy(false),
99 _model("!model","Original resolution model",this,(RooResolutionModel&)model,false,false),
100 _convVar("!convVar","Convolution variable",this,cVar,false,false),
101 _convSet("!convSet","Set of resModel X basisFunc convolutions",this),
102 _coefNormMgr(this,10),
103 _codeReg(10)
104{
105 _model.absArg()->setAttribute("NOCacheAndTrack") ;
106}
107
108
109
110////////////////////////////////////////////////////////////////////////////////
111
113 RooAbsPdf(other,name), _isCopy(true),
114 _model("!model",this,other._model),
115 _convVar("!convVar",this,other._convVar),
116 _convSet("!convSet",this,other._convSet),
117 _coefNormMgr(other._coefNormMgr,this),
118 _codeReg(other._codeReg)
119{
120 // Copy constructor
121 if (_model.absArg()) {
122 _model.absArg()->setAttribute("NOCacheAndTrack") ;
123 }
124 other._basisList.snapshot(_basisList);
125}
126
127
128
129////////////////////////////////////////////////////////////////////////////////
130/// Destructor
131
133{
134 if (!_isCopy) {
135 std::vector<RooAbsArg*> tmp(_convSet.begin(), _convSet.end());
136
137 for (auto arg : tmp) {
138 _convSet.remove(*arg) ;
139 delete arg ;
140 }
141 }
142
143}
144
145
146////////////////////////////////////////////////////////////////////////////////
147/// Declare a basis function for use in this physics model. The string expression
148/// must be a valid RooFormulVar expression representing the basis function, referring
149/// to the convolution variable as '@0', and any additional parameters (supplied in
150/// 'params' as '@1','@2' etc.
151///
152/// The return value is a unique identifier code, that will be passed to coefficient()
153/// to identify the basis function for which the coefficient is requested. If the
154/// resolution model used does not support the declared basis function, code -1 is
155/// returned.
156///
157
158Int_t RooAbsAnaConvPdf::declareBasis(const char* expression, const RooArgList& params)
159{
160 // Sanity check
161 if (_isCopy) {
162 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): ERROR attempt to "
163 << " declare basis functions in a copied RooAbsAnaConvPdf" << std::endl ;
164 return -1 ;
165 }
166
167 // Resolution model must support declared basis
168 if (!(static_cast<RooResolutionModel*>(_model.absArg()))->isBasisSupported(expression)) {
169 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): resolution model "
170 << _model.absArg()->GetName()
171 << " doesn't support basis function " << expression << std::endl ;
172 return -1 ;
173 }
174
175 // Instantiate basis function
177 basisArgs.add(params) ;
178
179 TString basisName(expression) ;
180 for (const auto arg : basisArgs) {
181 basisName.Append("_") ;
182 basisName.Append(arg->GetName()) ;
183 }
184
185 auto basisFunc = std::make_unique<RooFormulaVar>(basisName, expression, basisArgs);
186 basisFunc->setAttribute("RooWorkspace::Recycle") ;
187 basisFunc->setAttribute("NOCacheAndTrack") ;
188 basisFunc->setOperMode(operMode()) ;
189
190 // Instantiate resModel x basisFunc convolution
191 RooAbsReal* conv = static_cast<RooResolutionModel*>(_model.absArg())->convolution(basisFunc.get(),this);
192 _basisList.addOwned(std::move(basisFunc));
193 if (!conv) {
194 coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): unable to construct convolution with basis function '"
195 << expression << "'" << std::endl ;
196 return -1 ;
197 }
198 _convSet.add(*conv) ;
199
200 return _convSet.index(conv) ;
201}
202
203
204
205////////////////////////////////////////////////////////////////////////////////
206/// Change the current resolution model to newModel
207
209{
211 bool allOK(true) ;
212 for (auto convArg : _convSet) {
213 auto conv = static_cast<RooResolutionModel*>(convArg);
214
215 // Build new resolution model
216 std::unique_ptr<RooResolutionModel> newConv{newModel.convolution(const_cast<RooFormulaVar*>(&conv->basis()),this)};
217 if (!newConvSet.addOwned(std::move(newConv))) {
218 allOK = false ;
219 break ;
220 }
221 }
222
223 // Check if all convolutions were successfully built
224 if (!allOK) {
225 return true ;
226 }
227
228 // Replace old convolutions with new set
230 _convSet.addOwned(std::move(newConvSet));
231
232 const std::string attrib = std::string("ORIGNAME:") + _model->GetName();
233 const bool oldAttrib = newModel.getAttribute(attrib.c_str());
234 const_cast<RooResolutionModel&>(newModel).setAttribute(attrib.c_str());
235
236 redirectServers(RooArgSet{newModel}, false, true);
237
238 // reset temporary attribute for server redirection
240
241 return false ;
242}
243
244
245
246
247////////////////////////////////////////////////////////////////////////////////
248/// Create a generator context for this p.d.f. If both the p.d.f and the resolution model
249/// support internal generation of the convolution observable on an infinite domain,
250/// deploy a specialized convolution generator context, which generates the physics distribution
251/// and the smearing separately, adding them a posteriori. If this is not possible return
252/// a (slower) generic generation context that uses accept/reject sampling
253
255 const RooArgSet* auxProto, bool verbose) const
256{
257 // Check if the resolution model specifies a special context to be used.
258 RooResolutionModel* conv = dynamic_cast<RooResolutionModel*>(_model.absArg());
259 assert(conv);
260
261 std::unique_ptr<RooArgSet> modelDep {_model->getObservables(&vars)};
262 modelDep->remove(*convVar(),true,true) ;
263 Int_t numAddDep = modelDep->size() ;
264
265 // Check if physics PDF and resolution model can both directly generate the convolution variable
266 RooArgSet dummy ;
267 bool pdfCanDir = (getGenerator(*convVar(),dummy) != 0) ;
268 bool 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() << std::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 true ;
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{
315 auto* conv = static_cast<RooResolutionModel*>(_convSet.at(0));
316 if (!conv) return nullptr;
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 result(0) ;
331
332 Int_t index(0) ;
333 for (auto convArg : _convSet) {
334 auto conv = static_cast<RooAbsPdf*>(convArg);
335 double coef = coefficient(index++) ;
336 if (coef!=0.) {
337 const double c = conv->getVal(nullptr);
338 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") val += coef*conv [" << index-1 << "/"
339 << _convSet.size() << "] coef = " << coef << " conv = " << c << std::endl ;
340 result += c * coef;
341 } else {
342 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") [" << index-1 << "/" << _convSet.size() << "] coef = 0" << std::endl ;
343 }
344 }
345
346 return result ;
347}
348
349
350
351////////////////////////////////////////////////////////////////////////////////
352/// Advertise capability to perform (analytical) integrals
353/// internally. For a given integration request over allVars while
354/// normalized over normSet2 and in range 'rangeName', returns
355/// largest subset that can be performed internally in analVars
356/// Return code is unique integer code identifying integration scenario
357/// to be passed to analyticalIntegralWN() to calculate requeste integral
358///
359/// Class RooAbsAnaConv defers analytical integration request to
360/// resolution model and/or coefficient implementations and
361/// aggregates results into composite configuration with a unique
362/// code assigned by RooAICRegistry
363
365 RooArgSet& analVars, const RooArgSet* normSet2, const char* /*rangeName*/) const
366{
367 // Handle trivial no-integration scenario
368 if (allVars.empty()) return 0 ;
369
370 if (_forceNumInt) return 0 ;
371
372 // Select subset of allVars that are actual dependents
374 getObservables(&allVars, allDeps);
375 std::unique_ptr<RooArgSet> normSet{normSet2 ? getObservables(normSet2) : nullptr};
376
377 RooArgSet intSetAll{allDeps,"intSetAll"};
378
379 // Split intSetAll in coef/conv parts
380 auto intCoefSet = std::make_unique<RooArgSet>("intCoefSet");
381 auto intConvSet = std::make_unique<RooArgSet>("intConvSet");
382
383 for (RooAbsArg * arg : intSetAll) {
384 bool ok(true) ;
385 for (RooAbsArg * conv : _convSet) {
386 if (conv->dependsOn(*arg)) ok=false ;
387 }
388
389 if (ok) {
390 intCoefSet->add(*arg) ;
391 } else {
392 intConvSet->add(*arg) ;
393 }
394
395 }
396
397 // Split normSetAll in coef/conv parts
398 auto normCoefSet = std::make_unique<RooArgSet>("normCoefSet");
399 auto normConvSet = std::make_unique<RooArgSet>("normConvSet");
400 if (normSet) {
401 for (RooAbsArg * arg : *normSet) {
402 bool ok(true) ;
403 for (RooAbsArg * conv : _convSet) {
404 if (conv->dependsOn(*arg)) ok=false ;
405 }
406
407 if (ok) {
408 normCoefSet->add(*arg) ;
409 } else {
410 normConvSet->add(*arg) ;
411 }
412
413 }
414 }
415
416 if (intCoefSet->empty()) intCoefSet.reset();
417 if (intConvSet->empty()) intConvSet.reset();
418 if (normCoefSet->empty()) normCoefSet.reset();
419 if (normConvSet->empty()) normConvSet.reset();
420
421
422 // Store integration configuration in registry
423 Int_t masterCode(0) ;
424 std::vector<Int_t> tmp(1, 0) ;
425
426 // takes ownership of all sets
428 intCoefSet.release(),
429 intConvSet.release(),
430 normCoefSet.release(),
431 normConvSet.release()) + 1;
432
433 analVars.add(allDeps) ;
434
435 return masterCode ;
436}
437
438
439
440
441////////////////////////////////////////////////////////////////////////////////
442/// Return analytical integral defined by given code, which is returned
443/// by getAnalyticalIntegralWN()
444///
445/// For unnormalized integrals the returned value is
446/// \f[
447/// \mathrm{PDF} = \sum_k \int \mathrm{coef}_k \; \mathrm{d}\bar{x}
448/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}\bar{y},
449/// \f]
450/// where \f$ \bar{x} \f$ is the set of coefficient dependents to be integrated,
451/// and \f$ \bar{y} \f$ the set of basis function dependents to be integrated.
452///
453/// For normalized integrals this becomes
454/// \f[
455/// \mathrm{PDF} = \frac{\sum_k \int \mathrm{coef}_k \; \mathrm{d}x
456/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}y}
457/// {\sum_k \int \mathrm{coef}_k \; \mathrm{d}v
458/// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}w},
459/// \f]
460/// where
461/// * \f$ x \f$ is the set of coefficient dependents to be integrated,
462/// * \f$ y \f$ the set of basis function dependents to be integrated,
463/// * \f$ v \f$ is the set of coefficient dependents over which is normalized and
464/// * \f$ w \f$ is the set of basis function dependents over which is normalized.
465///
466/// Set \f$ x \f$ must be contained in \f$ v \f$ and set \f$ y \f$ must be contained in \f$ w \f$.
467///
468
470{
471 // WVE needs adaptation to handle new rangeName feature
472
473 // Handle trivial passthrough scenario
474 if (code == 0)
475 return getVal(normSet);
476
477 // Unpack master code
483
484 Int_t index(0);
485
486 if (normCoefSet == nullptr && normConvSet == nullptr) {
487 // Integral over unnormalized function
488 double integral(0);
490 for (auto *conv : static_range_cast<RooAbsPdf *>(_convSet)) {
491 double coef = getCoefNorm(index++, intCoefSet, rangeNamePtr);
492 if (coef != 0) {
493 const double term = coef * conv->getNormObj(nullptr, intConvSet, rangeNamePtr)->getVal();
494 integral += term;
495 cxcoutD(Eval) << "RooAbsAnaConv::aiWN(" << GetName() << ") [" << index - 1 << "] integral += " << term
496 << std::endl;
497 }
498 }
499 return integral;
500 }
501
502 // Integral over normalized function
503 double integral(0);
504 double norm(0);
506 for (auto *conv : static_range_cast<RooAbsPdf *>(_convSet)) {
507
509 if (coefInt != 0) {
510 double term = conv->getNormObj(nullptr, intConvSet, rangeNamePtr)->getVal();
511 integral += coefInt * term;
512 }
513
515 if (coefNorm != 0) {
516 double term = conv->getNormObj(nullptr, normConvSet)->getVal();
517 norm += coefNorm * term;
518 }
519
520 index++;
521 }
522 return integral / norm;
523}
524
525
526
527////////////////////////////////////////////////////////////////////////////////
528/// Default implementation of function advertising integration capabilities. The interface is
529/// similar to that of getAnalyticalIntegral except that an integer code is added that
530/// designates the coefficient number for which the integration capabilities are requested
531///
532/// This default implementation advertises that no internal integrals are supported.
533
534Int_t RooAbsAnaConvPdf::getCoefAnalyticalIntegral(Int_t /* coef*/, RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
535{
536 return 0 ;
537}
538
539
540
541////////////////////////////////////////////////////////////////////////////////
542/// Default implementation of function implementing advertised integrals. Only
543/// the pass-through scenario (no integration) is implemented.
544
545double RooAbsAnaConvPdf::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* /*rangeName*/) const
546{
547 if (code==0) return coefficient(coef) ;
548 coutE(InputArguments) << "RooAbsAnaConvPdf::coefAnalyticalIntegral(" << GetName() << ") ERROR: unrecognized integration code: " << code << std::endl ;
549 assert(0) ;
550 return 1 ;
551}
552
553
554
555////////////////////////////////////////////////////////////////////////////////
556/// This function forces RooRealIntegral to offer all integration dependents
557/// to RooAbsAnaConvPdf::getAnalyticalIntegralWN() for consideration for
558/// internal integration, if RooRealIntegral considers this to be unsafe (e.g. due
559/// to hidden Jacobian terms).
560///
561/// RooAbsAnaConvPdf will not attempt to actually integrate all these dependents
562/// but feed them to the resolution models integration interface, which will
563/// make the final determination on how to integrate these dependents.
564
566{
567 return true ;
568}
569
570
571
572////////////////////////////////////////////////////////////////////////////////
573/// Returns the normalization integral value of the coefficient with number coefIdx over normalization
574/// set nset in range rangeName
575
577{
578 if (nset==nullptr) return coefficient(coefIdx) ;
579
580 CacheElem* cache = static_cast<CacheElem*>(_coefNormMgr.getObj(nset,nullptr,nullptr,rangeName)) ;
581 if (!cache) {
582
583 cache = new CacheElem ;
584
585 // Make list of coefficient normalizations
587
588 for (std::size_t i=0 ; i<cache->_coefVarList.size() ; i++) {
589 cache->_normList.addOwned(std::unique_ptr<RooAbsReal>{static_cast<RooAbsReal&>(*cache->_coefVarList.at(i)).createIntegral(*nset,RooNameReg::str(rangeName))});
590 }
591
592 _coefNormMgr.setObj(nset,nullptr,cache,rangeName) ;
593 }
594
595 return (static_cast<RooAbsReal*>(cache->_normList.at(coefIdx)))->getVal() ;
596}
597
598
599
600////////////////////////////////////////////////////////////////////////////////
601/// Build complete list of coefficient variables
602
604{
605 // Instantiate a coefficient variables
606 for (std::size_t i=0 ; i<_convSet.size() ; i++) {
607 auto cvars = coefVars(i);
608 std::string name = std::string{GetName()} + "_coefVar_" + std::to_string(i);
609 varList.addOwned(std::make_unique<RooConvCoefVar>(name.c_str(),"coefVar",*this,i,&*cvars));
610 }
611
612}
613
614
615////////////////////////////////////////////////////////////////////////////////
616/// Return set of parameters with are used exclusively by the coefficient functions
617
619{
620 std::unique_ptr<RooArgSet> cVars{getParameters(static_cast<RooArgSet*>(nullptr))};
621 std::vector<RooAbsArg*> tmp;
622 for (auto arg : *cVars) {
623 for (auto convSetArg : _convSet) {
624 if (convSetArg->dependsOn(*arg)) {
625 tmp.push_back(arg);
626 }
627 }
628 }
629
630 cVars->remove(tmp.begin(), tmp.end(), true, true);
631
632 return RooFit::makeOwningPtr(std::move(cVars));
633}
634
635
636
637
638////////////////////////////////////////////////////////////////////////////////
639/// Print info about this object to the specified stream. In addition to the info
640/// from RooAbsPdf::printStream() we add:
641///
642/// Verbose : detailed information on convolution integrals
643
644void RooAbsAnaConvPdf::printMultiline(ostream& os, Int_t contents, bool verbose, TString indent) const
645{
646 RooAbsPdf::printMultiline(os,contents,verbose,indent);
647
648 os << indent << "--- RooAbsAnaConvPdf ---" << std::endl;
649 for (RooAbsArg * conv : _convSet) {
650 conv->printMultiline(os,contents,verbose,indent) ;
651 }
652}
653
654
655///////////////////////////////////////////////////////////////////////////////
656/// Label OK'ed components with cache-and-track
658{
659 for (auto const* carg : static_range_cast<RooAbsArg*>(_convSet)) {
660 if (carg->canNodeBeCached()==Always) {
661 trackNodes.add(*carg) ;
662 //cout << "tracking node RooAddPdf component " << carg->ClassName() << "::" << carg->GetName() << std::endl ;
663 }
664 }
665}
666
667std::unique_ptr<RooAbsArg>
669{
670 // If there is only one component in the linear sum of convolutions, we can
671 // just return that one, normalized.
672 if(_convSet.size() == 1) {
673 if (normSet.empty()) {
674 return _convSet[0].compileForNormSet(normSet, ctx);
675 }
676 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(_convSet[0].Clone()));
678
679 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, normSet);
680
681 // The direct servers are this pdf and the normalization integral, which
682 // don't need to be compiled further.
683 for (RooAbsArg *server : newArg->servers()) {
684 server->setAttribute("_COMPILED");
685 }
686 newArg->setAttribute("_COMPILED");
687 newArg->addOwnedComponents(std::move(pdfClone));
688 return newArg;
689 }
690
691 // Here, we can't use directly the function from the RooAbsPdf base class,
692 // because the convolution argument servers need to be evaluated
693 // unnormalized, even if they are pdfs.
694
695 if (normSet.empty()) {
697 }
698 std::unique_ptr<RooAbsAnaConvPdf> pdfClone(static_cast<RooAbsAnaConvPdf *>(this->Clone()));
699
700 // The actual resolution model is not serving the RooAbsAnaConvPdf
701 // in the evaluation. It was only used get the convolutions with a given
702 // basis. We can remove it for the compiled model.
703 pdfClone->removeServer(const_cast<RooAbsReal &>(pdfClone->_model.arg()), true);
704
705 // The other servers will be compiled with the original normSet, but the
706 // _convSet has to be evaluated unnormalized.
708 for (RooAbsArg *convArg : _convSet) {
709 if (auto convArgClone = ctx.compile(*convArg, *pdfClone, {})) {
711 }
712 }
713 pdfClone->redirectServers(convArgClones, false, true);
714
715 // Compile remaining servers that are evaluated normalized
717
718 // Finally, this RooAbsAnaConvPdf needs to be normalized
719 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, normSet);
720
721 // The direct servers are this pdf and the normalization integral, which
722 // don't need to be compiled further.
723 for (RooAbsArg *server : newArg->servers()) {
724 server->setAttribute("_COMPILED");
725 }
726 newArg->setAttribute("_COMPILED");
727 newArg->addOwnedComponents(std::move(pdfClone));
728 return newArg;
729}
#define c(i)
Definition RSha256.hxx:101
#define coutI(a)
#define cxcoutD(a)
#define coutE(a)
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
const_iterator begin() const
const_iterator end() const
const std::vector< Int_t > & retrieve(Int_t masterCode) const
Retrieve the array of integer codes associated with the given master code.
Int_t store(const std::vector< Int_t > &codeList, RooArgSet *set1=nullptr, RooArgSet *set2=nullptr, RooArgSet *set3=nullptr, RooArgSet *set4=nullptr)
Store given arrays of integer codes, and up to four RooArgSets in the registry (each setX pointer may...
Base class for PDFs that represent a physics model that can be analytically convolved with a resoluti...
friend class RooConvGenContext
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise capability to perform (analytical) integrals internally.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral defined by given code, which is returned by getAnalyticalIntegralWN()
virtual double coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=nullptr) const
Default implementation of function implementing advertised integrals.
virtual bool changeModel(const RooResolutionModel &newModel)
Change the current resolution model to newModel.
double getCoefNorm(Int_t coefIdx, const RooArgSet &nset, const char *rangeName) const
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components with cache-and-track.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
This function forces RooRealIntegral to offer all integration dependents to RooAbsAnaConvPdf::getAnal...
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Create a generator context for this p.d.f.
virtual double coefficient(Int_t basisIndex) const =0
RooArgList _basisList
! List of created basis functions
RooObjCacheManager _coefNormMgr
! Coefficient normalization manager
void makeCoefVarList(RooArgList &) const
Build complete list of coefficient variables.
RooAICRegistry _codeReg
! Registry of analytical integration codes
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const
Default implementation of function advertising integration capabilities.
~RooAbsAnaConvPdf() override
Destructor.
bool isDirectGenSafe(const RooAbsArg &arg) const override
Return true if it is safe to generate the convolution observable from the internal generator (this is...
RooAbsRealLValue * convVar()
Retrieve the convolution variable.
double evaluate() const override
Calculate the current unnormalized value of the PDF.
RooRealProxy _model
Original model.
void printMultiline(std::ostream &stream, Int_t contents, bool verbose=false, TString indent="") const override
Print info about this object to the specified stream.
RooRealProxy _convVar
Convolution variable.
RooAbsAnaConvPdf()
Default constructor, required for persistence.
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
RooListProxy _convSet
Set of (resModel (x) basisFunc) convolution objects.
virtual RooFit::OwningPtr< RooArgSet > coefVars(Int_t coefIdx) const
Return set of parameters with are used exclusively by the coefficient functions.
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
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.
bool redirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool isRecursionStep=false)
Replace all direct servers of this object with the new servers in newServerList.
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:89
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:425
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
const_iterator end() const
Storage_t::size_type size() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
const_iterator begin() const
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print multi line detailed information of this RooAbsPdf.
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
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:538
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
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
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
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 addOwned(RooAbsArg &var, bool silent=false) override
Overloaded RooCollection_t::addOwned() method insert object into owning set and registers object as s...
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...
bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false) override
Remove object 'var' from set and deregister 'var' as server to owner.
Container class to hold unbinned data.
Definition RooDataSet.h:34
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
T * compile(T &arg, RooAbsArg &owner, RooArgSet const &normSet)
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Implements a universal generator context for all RooAbsPdf classes that do not have or need a special...
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.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
virtual RooAbsGenContext * modelGenContext(const RooAbsAnaConvPdf &, const RooArgSet &, const RooDataSet *, const RooArgSet *, bool) const
RooAbsRealLValue & convVar() const
Return the convolution variable of the resolution model.
const T & arg() const
Return reference to object held in proxy.
Implements a RooResolution model that corresponds to a delta function.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Basic string class.
Definition TString.h:139
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
Definition TString.cxx:457
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40