Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
64 #include "RooFit.h"
65 #include "RooMsgService.h"
66 
67 #include "Riostream.h"
68 #include "Riostream.h"
69 #include "RooAbsAnaConvPdf.h"
70 #include "RooResolutionModel.h"
71 #include "RooRealVar.h"
72 #include "RooFormulaVar.h"
73 #include "RooConvGenContext.h"
74 #include "RooGenContext.h"
75 #include "RooTruthModel.h"
76 #include "RooConvCoefVar.h"
77 #include "RooNameReg.h"
78 
79 using namespace std;
80 
81 ClassImp(RooAbsAnaConvPdf);
82 
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 /// Default constructor, required for persistence
86 
87 RooAbsAnaConvPdf::RooAbsAnaConvPdf() :
88  _isCopy(kFALSE),
89  _convNormSet(nullptr)
90 {
91 }
92 
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Constructor. The supplied resolution model must be constructed with the same
97 /// convoluted variable as this physics model ('convVar')
98 
99 RooAbsAnaConvPdf::RooAbsAnaConvPdf(const char *name, const char *title,
100  const RooResolutionModel& model, RooRealVar& cVar) :
101  RooAbsPdf(name,title), _isCopy(kFALSE),
102  _model("!model","Original resolution model",this,(RooResolutionModel&)model,kFALSE,kFALSE),
103  _convVar("!convVar","Convolution variable",this,cVar,kFALSE,kFALSE),
104  _convSet("!convSet","Set of resModel X basisFunc convolutions",this),
105  _convNormSet(nullptr),
106  _coefNormMgr(this,10),
107  _codeReg(10)
108 {
109  _convNormSet = new RooArgSet(cVar,"convNormSet") ;
110  _model.absArg()->setAttribute("NOCacheAndTrack") ;
111 }
112 
113 
114 
115 ////////////////////////////////////////////////////////////////////////////////
116 
117 RooAbsAnaConvPdf::RooAbsAnaConvPdf(const RooAbsAnaConvPdf& other, const char* name) :
118  RooAbsPdf(other,name), _isCopy(kTRUE),
119  _model("!model",this,other._model),
120  _convVar("!convVar",this,other._convVar),
121  _convSet("!convSet",this,other._convSet),
122  // _basisList(other._basisList),
123  _convNormSet(other._convNormSet? new RooArgSet(*other._convNormSet) : new RooArgSet() ),
124  _coefNormMgr(other._coefNormMgr,this),
125  _codeReg(other._codeReg)
126 {
127  // Copy constructor
128  if (_model.absArg()) {
129  _model.absArg()->setAttribute("NOCacheAndTrack") ;
130  }
131  other._basisList.snapshot(_basisList);
132 }
133 
134 
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Destructor
138 
139 RooAbsAnaConvPdf::~RooAbsAnaConvPdf()
140 {
141  if (_convNormSet) {
142  delete _convNormSet ;
143  }
144 
145  if (!_isCopy) {
146  std::vector<RooAbsArg*> tmp(_convSet.begin(), _convSet.end());
147 
148  for (auto arg : tmp) {
149  _convSet.remove(*arg) ;
150  delete arg ;
151  }
152  }
153 
154 }
155 
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Declare a basis function for use in this physics model. The string expression
159 /// must be a valid RooFormulVar expression representing the basis function, referring
160 /// to the convolution variable as '@0', and any additional parameters (supplied in
161 /// 'params' as '@1','@2' etc.
162 ///
163 /// The return value is a unique identifier code, that will be passed to coefficient()
164 /// to identify the basis function for which the coefficient is requested. If the
165 /// resolution model used does not support the declared basis function, code -1 is
166 /// returned.
167 ///
168 
169 Int_t RooAbsAnaConvPdf::declareBasis(const char* expression, const RooArgList& params)
170 {
171  // Sanity check
172  if (_isCopy) {
173  coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): ERROR attempt to "
174  << " declare basis functions in a copied RooAbsAnaConvPdf" << endl ;
175  return -1 ;
176  }
177 
178  // Resolution model must support declared basis
179  if (!((RooResolutionModel*)_model.absArg())->isBasisSupported(expression)) {
180  coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): resolution model "
181  << _model.absArg()->GetName()
182  << " doesn't support basis function " << expression << endl ;
183  return -1 ;
184  }
185 
186  // Instantiate basis function
187  RooArgList basisArgs(_convVar.arg()) ;
188  basisArgs.add(params) ;
189 
190  TString basisName(expression) ;
191  TIterator* iter = basisArgs.createIterator() ;
192  RooAbsArg* arg ;
193  while(((arg=(RooAbsArg*)iter->Next()))) {
194  basisName.Append("_") ;
195  basisName.Append(arg->GetName()) ;
196  }
197  delete iter ;
198 
199  RooFormulaVar* basisFunc = new RooFormulaVar(basisName,expression,basisArgs) ;
200  basisFunc->setAttribute("RooWorkspace::Recycle") ;
201  basisFunc->setAttribute("NOCacheAndTrack") ;
202  basisFunc->setOperMode(operMode()) ;
203  _basisList.addOwned(*basisFunc) ;
204 
205  // Instantiate resModel x basisFunc convolution
206  RooAbsReal* conv = ((RooResolutionModel*)_model.absArg())->convolution(basisFunc,this) ;
207  if (!conv) {
208  coutE(InputArguments) << "RooAbsAnaConvPdf::declareBasis(" << GetName() << "): unable to construct convolution with basis function '"
209  << expression << "'" << endl ;
210  return -1 ;
211  }
212  _convSet.add(*conv) ;
213 
214  return _convSet.index(conv) ;
215 }
216 
217 
218 
219 ////////////////////////////////////////////////////////////////////////////////
220 /// Change the current resolution model to newModel
221 
222 Bool_t RooAbsAnaConvPdf::changeModel(const RooResolutionModel& newModel)
223 {
224  RooArgList newConvSet ;
225  Bool_t allOK(kTRUE) ;
226  for (auto convArg : _convSet) {
227  auto conv = static_cast<RooResolutionModel*>(convArg);
228 
229  // Build new resolution model
230  RooResolutionModel* newConv = newModel.convolution((RooFormulaVar*)&conv->basis(),this) ;
231  if (!newConvSet.add(*newConv)) {
232  allOK = kFALSE ;
233  break ;
234  }
235  }
236 
237  // Check if all convolutions were successfully built
238  if (!allOK) {
239  // Delete new basis functions created sofar
240  std::for_each(newConvSet.begin(), newConvSet.end(), [](RooAbsArg* arg){delete arg;});
241 
242  return kTRUE ;
243  }
244 
245  // Replace old convolutions with new set
246  _convSet.removeAll() ;
247  _convSet.addOwned(newConvSet) ;
248 
249  // Update server link by hand, since _model.setArg() below will not do this
250  replaceServer((RooAbsArg&)_model.arg(),(RooAbsArg&)newModel,kFALSE,kFALSE) ;
251 
252  _model.setArg((RooResolutionModel&)newModel) ;
253  return kFALSE ;
254 }
255 
256 
257 
258 
259 ////////////////////////////////////////////////////////////////////////////////
260 /// Create a generator context for this p.d.f. If both the p.d.f and the resolution model
261 /// support internal generation of the convolution observable on an infinite domain,
262 /// deploy a specialized convolution generator context, which generates the physics distribution
263 /// and the smearing separately, adding them a posteriori. If this is not possible return
264 /// a (slower) generic generation context that uses accept/reject sampling
265 
266 RooAbsGenContext* RooAbsAnaConvPdf::genContext(const RooArgSet &vars, const RooDataSet *prototype,
267  const RooArgSet* auxProto, Bool_t verbose) const
268 {
269  // Check if the resolution model specifies a special context to be used.
270  RooResolutionModel* conv = dynamic_cast<RooResolutionModel*>(_model.absArg());
271  assert(conv);
272 
273  RooArgSet* modelDep = _model.absArg()->getObservables(&vars) ;
274  modelDep->remove(*convVar(),kTRUE,kTRUE) ;
275  Int_t numAddDep = modelDep->getSize() ;
276  delete modelDep ;
277 
278  // Check if physics PDF and resolution model can both directly generate the convolution variable
279  RooArgSet dummy ;
280  Bool_t pdfCanDir = (getGenerator(*convVar(),dummy) != 0) ;
281  Bool_t resCanDir = conv && (conv->getGenerator(*convVar(),dummy)!=0) && conv->isDirectGenSafe(*convVar()) ;
282 
283  if (numAddDep>0 || !pdfCanDir || !resCanDir) {
284  // Any resolution model with more dependents than the convolution variable
285  // or pdf or resmodel do not support direct generation
286  string reason ;
287  if (numAddDep>0) reason += "Resolution model has more observables than the convolution variable. " ;
288  if (!pdfCanDir) reason += "PDF does not support internal generation of convolution observable. " ;
289  if (!resCanDir) reason += "Resolution model does not support internal generation of convolution observable. " ;
290 
291  coutI(Generation) << "RooAbsAnaConvPdf::genContext(" << GetName() << ") Using regular accept/reject generator for convolution p.d.f because: " << reason.c_str() << endl ;
292  return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
293  }
294 
295  RooAbsGenContext* context = conv->modelGenContext(*this, vars, prototype, auxProto, verbose);
296  if (context) return context;
297 
298  // Any other resolution model: use specialized generator context
299  return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
300 }
301 
302 
303 
304 ////////////////////////////////////////////////////////////////////////////////
305 /// Return true if it is safe to generate the convolution observable
306 /// from the internal generator (this is the case if the chosen resolution
307 /// model is the truth model)
308 
309 Bool_t RooAbsAnaConvPdf::isDirectGenSafe(const RooAbsArg& arg) const
310 {
311 
312  // All direct generation of convolution arg if model is truth model
313  if (!TString(_convVar.absArg()->GetName()).CompareTo(arg.GetName()) &&
314  dynamic_cast<RooTruthModel*>(_model.absArg())) {
315  return kTRUE ;
316  }
317 
318  return RooAbsPdf::isDirectGenSafe(arg) ;
319 }
320 
321 
322 
323 ////////////////////////////////////////////////////////////////////////////////
324 /// Return a pointer to the convolution variable instance used in the resolution model
325 
326 const RooRealVar* RooAbsAnaConvPdf::convVar() const
327 {
328  RooResolutionModel* conv = (RooResolutionModel*) _convSet.at(0) ;
329  if (!conv) return 0 ;
330  return &conv->convVar() ;
331 }
332 
333 
334 
335 ////////////////////////////////////////////////////////////////////////////////
336 /// Calculate the current unnormalized value of the PDF
337 ///
338 /// PDF = sum_k coef_k * [ basis_k (x) ResModel ]
339 ///
340 
341 Double_t RooAbsAnaConvPdf::evaluate() const
342 {
343  Double_t result(0) ;
344 
345  Int_t index(0) ;
346  for (auto convArg : _convSet) {
347  auto conv = static_cast<RooAbsPdf*>(convArg);
348  Double_t coef = coefficient(index++) ;
349  if (coef!=0.) {
350  Double_t c = conv->getVal(0) ;
351  Double_t r = coef ;
352  cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") val += coef*conv [" << index-1 << "/"
353  << _convSet.getSize() << "] coef = " << r << " conv = " << c << endl ;
354  result += conv->getVal(0)*coef ;
355  } else {
356  cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") [" << index-1 << "/" << _convSet.getSize() << "] coef = 0" << endl ;
357  }
358  }
359 
360  return result ;
361 }
362 
363 
364 
365 ////////////////////////////////////////////////////////////////////////////////
366 /// Advertise capability to perform (analytical) integrals
367 /// internally. For a given integration request over allVars while
368 /// normalized over normSet2 and in range 'rangeName', returns
369 /// largest subset that can be performed internally in analVars
370 /// Return code is unique integer code identifying integration scenario
371 /// to be passed to analyticalIntegralWN() to calculate requeste integral
372 ///
373 /// Class RooAbsAnaConv defers analytical integration request to
374 /// resolution model and/or coefficient implementations and
375 /// aggregates results into composite configuration with a unique
376 /// code assigned by RooAICRegistry
377 
378 Int_t RooAbsAnaConvPdf::getAnalyticalIntegralWN(RooArgSet& allVars,
379  RooArgSet& analVars, const RooArgSet* normSet2, const char* /*rangeName*/) const
380 {
381  // Handle trivial no-integration scenario
382  if (allVars.getSize()==0) return 0 ;
383 
384  if (_forceNumInt) return 0 ;
385 
386  // Select subset of allVars that are actual dependents
387  RooArgSet* allDeps = getObservables(allVars) ;
388  RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
389 
390  RooAbsArg *arg ;
391  RooResolutionModel *conv ;
392 
393  RooArgSet* intSetAll = new RooArgSet(*allDeps,"intSetAll") ;
394 
395  // Split intSetAll in coef/conv parts
396  RooArgSet* intCoefSet = new RooArgSet("intCoefSet") ;
397  RooArgSet* intConvSet = new RooArgSet("intConvSet") ;
398  TIterator* varIter = intSetAll->createIterator() ;
399  TIterator* convIter = _convSet.createIterator() ;
400 
401  while(((arg=(RooAbsArg*) varIter->Next()))) {
402  Bool_t ok(kTRUE) ;
403  convIter->Reset() ;
404  while(((conv=(RooResolutionModel*) convIter->Next()))) {
405  if (conv->dependsOn(*arg)) ok=kFALSE ;
406  }
407 
408  if (ok) {
409  intCoefSet->add(*arg) ;
410  } else {
411  intConvSet->add(*arg) ;
412  }
413 
414  }
415  delete varIter ;
416 
417 
418  // Split normSetAll in coef/conv parts
419  RooArgSet* normCoefSet = new RooArgSet("normCoefSet") ;
420  RooArgSet* normConvSet = new RooArgSet("normConvSet") ;
421  RooArgSet* normSetAll = normSet ? (new RooArgSet(*normSet,"normSetAll")) : 0 ;
422  if (normSetAll) {
423  varIter = normSetAll->createIterator() ;
424  while(((arg=(RooAbsArg*) varIter->Next()))) {
425  Bool_t ok(kTRUE) ;
426  convIter->Reset() ;
427  while(((conv=(RooResolutionModel*) convIter->Next()))) {
428  if (conv->dependsOn(*arg)) ok=kFALSE ;
429  }
430 
431  if (ok) {
432  normCoefSet->add(*arg) ;
433  } else {
434  normConvSet->add(*arg) ;
435  }
436 
437  }
438  delete varIter ;
439  }
440  delete convIter ;
441 
442  if (intCoefSet->getSize()==0) {
443  delete intCoefSet ; intCoefSet=0 ;
444  }
445  if (intConvSet->getSize()==0) {
446  delete intConvSet ; intConvSet=0 ;
447  }
448  if (normCoefSet->getSize()==0) {
449  delete normCoefSet ; normCoefSet=0 ;
450  }
451  if (normConvSet->getSize()==0) {
452  delete normConvSet ; normConvSet=0 ;
453  }
454 
455 
456 
457  // Store integration configuration in registry
458  Int_t masterCode(0) ;
459  std::vector<Int_t> tmp(1, 0) ;
460 
461  masterCode = _codeReg.store(tmp, intCoefSet, intConvSet, normCoefSet, normConvSet) + 1 ; // takes ownership of all sets
462 
463  analVars.add(*allDeps) ;
464  delete allDeps ;
465  if (normSet) delete normSet ;
466  if (normSetAll) delete normSetAll ;
467  delete intSetAll ;
468 
469 // cout << this << "---> masterCode = " << masterCode << endl ;
470 
471  return masterCode ;
472 }
473 
474 
475 
476 
477 ////////////////////////////////////////////////////////////////////////////////
478 /// Return analytical integral defined by given code, which is returned
479 /// by getAnalyticalIntegralWN()
480 ///
481 /// For unnormalized integrals the returned value is
482 /// \f[
483 /// \mathrm{PDF} = \sum_k \int \mathrm{coef}_k \; \mathrm{d}\bar{x}
484 /// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}\bar{y},
485 /// \f]
486 /// where \f$ \bar{x} \f$ is the set of coefficient dependents to be integrated,
487 /// and \f$ \bar{y} \f$ the set of basis function dependents to be integrated.
488 ///
489 /// For normalized integrals this becomes
490 /// \f[
491 /// \mathrm{PDF} = \frac{\sum_k \int \mathrm{coef}_k \; \mathrm{d}x
492 /// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}y}
493 /// {\sum_k \int \mathrm{coef}_k \; \mathrm{d}v
494 /// \cdot \int \mathrm{basis}_k (x) \mathrm{ResModel} \; \mathrm{d}w},
495 /// \f]
496 /// where
497 /// * \f$ x \f$ is the set of coefficient dependents to be integrated,
498 /// * \f$ y \f$ the set of basis function dependents to be integrated,
499 /// * \f$ v \f$ is the set of coefficient dependents over which is normalized and
500 /// * \f$ w \f$ is the set of basis function dependents over which is normalized.
501 ///
502 /// Set \f$ x \f$ must be contained in \f$ v \f$ and set \f$ y \f$ must be contained in \f$ w \f$.
503 ///
504 
505 Double_t RooAbsAnaConvPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
506 {
507  // WVE needs adaptation to handle new rangeName feature
508 
509  // Handle trivial passthrough scenario
510  if (code==0) return getVal(normSet) ;
511 
512  // Unpack master code
513  RooArgSet *intCoefSet, *intConvSet, *normCoefSet, *normConvSet ;
514  _codeReg.retrieve(code-1,intCoefSet,intConvSet,normCoefSet,normConvSet) ;
515 
516  Int_t index(0) ;
517  Double_t answer(0) ;
518 
519  if (normCoefSet==0&&normConvSet==0) {
520 
521  // Integral over unnormalized function
522  Double_t integral(0) ;
523  const TNamed *_rangeName = RooNameReg::ptr(rangeName);
524  for (auto convArg : _convSet) {
525  auto conv = static_cast<RooResolutionModel*>(convArg);
526  Double_t coef = getCoefNorm(index++,intCoefSet,_rangeName) ;
527  //cout << "coefInt[" << index << "] = " << coef << " " ; intCoefSet->Print("1") ;
528  if (coef!=0) {
529  integral += coef*(_rangeName ? conv->getNormObj(0,intConvSet,_rangeName)->getVal() : conv->getNorm(intConvSet) ) ;
530  cxcoutD(Eval) << "RooAbsAnaConv::aiWN(" << GetName() << ") [" << index-1 << "] integral += " << conv->getNorm(intConvSet) << endl ;
531  }
532 
533  }
534  answer = integral ;
535 
536  } else {
537 
538  // Integral over normalized function
539  Double_t integral(0) ;
540  Double_t norm(0) ;
541  const TNamed *_rangeName = RooNameReg::ptr(rangeName);
542  for (auto convArg : _convSet) {
543  auto conv = static_cast<RooResolutionModel*>(convArg);
544 
545  Double_t coefInt = getCoefNorm(index,intCoefSet,_rangeName) ;
546  //cout << "coefInt[" << index << "] = " << coefInt << "*" << term << " " << (intCoefSet?*intCoefSet:RooArgSet()) << endl ;
547  if (coefInt!=0) {
548  Double_t term = (_rangeName ? conv->getNormObj(0,intConvSet,_rangeName)->getVal() : conv->getNorm(intConvSet) ) ;
549  integral += coefInt*term ;
550  }
551 
552  Double_t coefNorm = getCoefNorm(index,normCoefSet) ;
553  //cout << "coefNorm[" << index << "] = " << coefNorm << "*" << term << " " << (normCoefSet?*normCoefSet:RooArgSet()) << endl ;
554  if (coefNorm!=0) {
555  Double_t term = conv->getNorm(normConvSet) ;
556  norm += coefNorm*term ;
557  }
558 
559  index++ ;
560  }
561  answer = integral/norm ;
562  }
563 
564  return answer ;
565 }
566 
567 
568 
569 ////////////////////////////////////////////////////////////////////////////////
570 /// Default implementation of function advertising integration capabilities. The interface is
571 /// similar to that of getAnalyticalIntegral except that an integer code is added that
572 /// designates the coefficient number for which the integration capabilities are requested
573 ///
574 /// This default implementation advertises that no internal integrals are supported.
575 
576 Int_t RooAbsAnaConvPdf::getCoefAnalyticalIntegral(Int_t /* coef*/, RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
577 {
578  return 0 ;
579 }
580 
581 
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// Default implementation of function implementing advertised integrals. Only
585 /// the pass-through scenario (no integration) is implemented.
586 
587 Double_t RooAbsAnaConvPdf::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* /*rangeName*/) const
588 {
589  if (code==0) return coefficient(coef) ;
590  coutE(InputArguments) << "RooAbsAnaConvPdf::coefAnalyticalIntegral(" << GetName() << ") ERROR: unrecognized integration code: " << code << endl ;
591  assert(0) ;
592  return 1 ;
593 }
594 
595 
596 
597 ////////////////////////////////////////////////////////////////////////////////
598 /// This function forces RooRealIntegral to offer all integration dependents
599 /// to RooAbsAnaConvPdf::getAnalyticalIntegralWN() for consideration for
600 /// internal integration, if RooRealIntegral considers this to be unsafe (e.g. due
601 /// to hidden Jacobian terms).
602 ///
603 /// RooAbsAnaConvPdf will not attempt to actually integrate all these dependents
604 /// but feed them to the resolution models integration interface, which will
605 /// make the final determination on how to integrate these dependents.
606 
607 Bool_t RooAbsAnaConvPdf::forceAnalyticalInt(const RooAbsArg& /*dep*/) const
608 {
609  return kTRUE ;
610 }
611 
612 
613 
614 ////////////////////////////////////////////////////////////////////////////////
615 /// Returns the normalization integral value of the coefficient with number coefIdx over normalization
616 /// set nset in range rangeName
617 
618 Double_t RooAbsAnaConvPdf::getCoefNorm(Int_t coefIdx, const RooArgSet* nset, const TNamed* rangeName) const
619 {
620  if (nset==0) return coefficient(coefIdx) ;
621 
622  CacheElem* cache = (CacheElem*) _coefNormMgr.getObj(nset,0,0,rangeName) ;
623  if (!cache) {
624 
625  cache = new CacheElem ;
626 
627  // Make list of coefficient normalizations
628  Int_t i ;
629  makeCoefVarList(cache->_coefVarList) ;
630 
631  for (i=0 ; i<cache->_coefVarList.getSize() ; i++) {
632  RooAbsReal* coefInt = static_cast<RooAbsReal&>(*cache->_coefVarList.at(i)).createIntegral(*nset,RooNameReg::str(rangeName)) ;
633  cache->_normList.addOwned(*coefInt) ;
634  }
635 
636  _coefNormMgr.setObj(nset,0,cache,rangeName) ;
637  }
638 
639  return ((RooAbsReal*)cache->_normList.at(coefIdx))->getVal() ;
640 }
641 
642 
643 
644 ////////////////////////////////////////////////////////////////////////////////
645 /// Build complete list of coefficient variables
646 
647 void RooAbsAnaConvPdf::makeCoefVarList(RooArgList& varList) const
648 {
649  // Instantate a coefficient variables
650  for (Int_t i=0 ; i<_convSet.getSize() ; i++) {
651  RooArgSet* cvars = coefVars(i) ;
652  RooAbsReal* coefVar = new RooConvCoefVar(Form("%s_coefVar_%d",GetName(),i),"coefVar",*this,i,cvars) ;
653  varList.addOwned(*coefVar) ;
654  delete cvars ;
655  }
656 
657 }
658 
659 
660 ////////////////////////////////////////////////////////////////////////////////
661 /// Return set of parameters with are used exclusively by the coefficient functions
662 
663 RooArgSet* RooAbsAnaConvPdf::coefVars(Int_t /*coefIdx*/) const
664 {
665  RooArgSet* cVars = getParameters((RooArgSet*)0) ;
666  std::vector<RooAbsArg*> tmp;
667  for (auto arg : *cVars) {
668  for (auto convSetArg : _convSet) {
669  if (convSetArg->dependsOn(*arg)) {
670  tmp.push_back(arg);
671  }
672  }
673  }
674 
675  cVars->remove(tmp.begin(), tmp.end(), true, true);
676 
677  return cVars ;
678 }
679 
680 
681 
682 
683 ////////////////////////////////////////////////////////////////////////////////
684 /// Print info about this object to the specified stream. In addition to the info
685 /// from RooAbsPdf::printStream() we add:
686 ///
687 /// Verbose : detailed information on convolution integrals
688 
689 void RooAbsAnaConvPdf::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
690 {
691  RooAbsPdf::printMultiline(os,contents,verbose,indent);
692 
693  os << indent << "--- RooAbsAnaConvPdf ---" << endl;
694  TIterator* iter = _convSet.createIterator() ;
695  RooResolutionModel* conv ;
696  while (((conv=(RooResolutionModel*)iter->Next()))) {
697  conv->printMultiline(os,contents,verbose,indent) ;
698  }
699 }
700 
701 
702 ///////////////////////////////////////////////////////////////////////////////
703 /// Label OK'ed components with cache-and-track
704 void RooAbsAnaConvPdf::setCacheAndTrackHints(RooArgSet& trackNodes)
705 {
706  RooFIter citer = _convSet.fwdIterator() ;
707  RooAbsArg* carg ;
708  while ((carg=citer.next())) {
709  if (carg->canNodeBeCached()==Always) {
710  trackNodes.add(*carg) ;
711  //cout << "tracking node RooAddPdf component " << carg->IsA()->GetName() << "::" << carg->GetName() << endl ;
712  }
713  }
714 }