Logo ROOT   6.16/01
Reference Guide
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//
19// RooAbsAnaConvPdf is the base class of for PDFs that represents a
20// physics model that can be analytically convolved with a resolution model
21//
22// To achieve factorization between the physics model and the resolution
23// model, each physics model must be able to be written in the form
24// _ _ _ _
25// Phys(x,a,b) = Sum_k coef_k(a) * basis_k(x,b)
26//
27// where basis_k are a limited number of functions in terms of the variable
28// to be convoluted and coef_k are coefficients independent of the convolution
29// variable.
30//
31// Classes derived from RooResolutionModel implement
32// _ _ _ _
33// R_k(x,b,c) = Int(dx') basis_k(x',b) * resModel(x-x',c)
34//
35// which RooAbsAnaConvPdf uses to construct the pdf for [ Phys (x) R ] :
36// _ _ _ _ _ _
37// PDF(x,a,b,c) = Sum_k coef_k(a) * R_k(x,b,c)
38//
39// A minimal implementation of a RooAbsAnaConvPdf physics model consists of
40//
41// - A constructor that declares the required basis functions using the declareBasis() method.
42// The declareBasis() function assigns a unique identifier code to each declare basis
43//
44// - An implementation of coefficient(Int_t code) returning the coefficient value for each
45// declared basis function
46//
47// Optionally, analytical integrals can be provided for the coefficient functions. The
48// interface for this is quite similar to that for integrals of regular PDFs. Two functions,
49//
50// Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
51// Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const
52//
53// advertise the coefficient integration capabilities and implement them respectively.
54// Please see RooAbsPdf for additional details. Advertised analytical integrals must be
55// valid for all coefficients.
56
57
58#include "RooFit.h"
59#include "RooMsgService.h"
60
61#include "Riostream.h"
62#include "Riostream.h"
63#include "RooAbsAnaConvPdf.h"
64#include "RooResolutionModel.h"
65#include "RooRealVar.h"
66#include "RooFormulaVar.h"
67#include "RooConvGenContext.h"
68#include "RooGenContext.h"
69#include "RooTruthModel.h"
70#include "RooConvCoefVar.h"
71#include "RooNameReg.h"
72
73using namespace std;
74
76;
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// Default constructor, required for persistence
81
83 _isCopy(kFALSE),
84 _convNormSet(0),
85 _convSetIter(_convSet.createIterator())
86{
87}
88
89
90
91////////////////////////////////////////////////////////////////////////////////
92/// Constructor. The supplied resolution model must be constructed with the same
93/// convoluted variable as this physics model ('convVar')
94
95RooAbsAnaConvPdf::RooAbsAnaConvPdf(const char *name, const char *title,
96 const RooResolutionModel& model, RooRealVar& cVar) :
97 RooAbsPdf(name,title), _isCopy(kFALSE),
98 _model("!model","Original resolution model",this,(RooResolutionModel&)model,kFALSE,kFALSE),
99 _convVar("!convVar","Convolution variable",this,cVar,kFALSE,kFALSE),
100 _convSet("!convSet","Set of resModel X basisFunc convolutions",this),
101 _convNormSet(0), _convSetIter(_convSet.createIterator()),
102 _coefNormMgr(this,10),
103 _codeReg(10)
104{
105 _convNormSet = new RooArgSet(cVar,"convNormSet") ;
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 // _basisList(other._basisList),
119 _convNormSet(other._convNormSet? new RooArgSet(*other._convNormSet) : new RooArgSet() ),
120 _convSetIter(_convSet.createIterator()),
121 _coefNormMgr(other._coefNormMgr,this),
122 _codeReg(other._codeReg)
123{
124 // Copy constructor
125 if (_model.absArg()) {
126 _model.absArg()->setAttribute("NOCacheAndTrack") ;
127 }
129}
130
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Destructor
135
137{
138 if (_convNormSet) {
139 delete _convNormSet ;
140 }
141
142 delete _convSetIter ;
143
144 if (!_isCopy) {
146 RooAbsArg* arg ;
147 while (((arg = (RooAbsArg*)iter->Next()))) {
148 _convSet.remove(*arg) ;
149 delete arg ;
150 }
151 delete iter ;
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
169Int_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
223{
225 RooResolutionModel* conv ;
226 RooArgList newConvSet ;
227 Bool_t allOK(kTRUE) ;
228 while(((conv=(RooResolutionModel*)cIter->Next()))) {
229
230 // Build new resolution model
231 RooResolutionModel* newConv = newModel.convolution((RooFormulaVar*)&conv->basis(),this) ;
232 if (!newConvSet.add(*newConv)) {
233 allOK = kFALSE ;
234 break ;
235 }
236 }
237 delete cIter ;
238
239 // Check if all convolutions were successfully built
240 if (!allOK) {
241 // Delete new basis functions created sofar
242 TIterator* iter = newConvSet.createIterator() ;
243 while(((conv=(RooResolutionModel*)iter->Next()))) delete conv ;
244 delete iter ;
245
246 return kTRUE ;
247 }
248
249 // Replace old convolutions with new set
251 _convSet.addOwned(newConvSet) ;
252
253 _model.setArg((RooResolutionModel&)newModel) ;
254 return kFALSE ;
255}
256
257
258
259
260////////////////////////////////////////////////////////////////////////////////
261/// Create a generator context for this p.d.f. If both the p.d.f and the resolution model
262/// support internal generation of the convolution observable on an infinite domain,
263/// deploy a specialized convolution generator context, which generates the physics distribution
264/// and the smearing separately, adding them a posteriori. If this is not possible return
265/// a (slower) generic generation context that uses accept/reject sampling
266
268 const RooArgSet* auxProto, Bool_t verbose) const
269{
270 // Check if the resolution model specifies a special context to be used.
271 RooResolutionModel* conv = dynamic_cast<RooResolutionModel*>(_model.absArg());
272 assert(conv);
273
274 RooArgSet* modelDep = _model.absArg()->getObservables(&vars) ;
275 modelDep->remove(*convVar(),kTRUE,kTRUE) ;
276 Int_t numAddDep = modelDep->getSize() ;
277 delete modelDep ;
278
279 // Check if physics PDF and resolution model can both directly generate the convolution variable
281 Bool_t pdfCanDir = (getGenerator(*convVar(),dummy) != 0) ;
282 Bool_t resCanDir = conv && (conv->getGenerator(*convVar(),dummy)!=0) && conv->isDirectGenSafe(*convVar()) ;
283
284 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
285 // Any resolution model with more dependents than the convolution variable
286 // or pdf or resmodel do not support direct generation
287 string reason ;
288 if (numAddDep>0) reason += "Resolution model has more onservables that the convolution variable. " ;
289 if (!pdfCanDir) reason += "PDF does not support internal generation of convolution observable. " ;
290 if (!resCanDir) reason += "Resolution model does not support internal generation of convolution observable. " ;
291
292 coutI(Generation) << "RooAbsAnaConvPdf::genContext(" << GetName() << ") Using regular accept/reject generator for convolution p.d.f because: " << reason.c_str() << endl ;
293 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
294 }
295
296 RooAbsGenContext* context = conv->modelGenContext(*this, vars, prototype, auxProto, verbose);
297 if (context) return context;
298
299 // Any other resolution model: use specialized generator context
300 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
301}
302
303
304
305////////////////////////////////////////////////////////////////////////////////
306/// Return true if it is safe to generate the convolution observable
307/// from the internal generator (this is the case if the chosen resolution
308/// model is the truth model)
309
311{
312
313 // All direct generation of convolution arg if model is truth model
314 if (!TString(_convVar.absArg()->GetName()).CompareTo(arg.GetName()) &&
315 dynamic_cast<RooTruthModel*>(_model.absArg())) {
316 return kTRUE ;
317 }
318
319 return RooAbsPdf::isDirectGenSafe(arg) ;
320}
321
322
323
324////////////////////////////////////////////////////////////////////////////////
325/// Return a pointer to the convolution variable instance used in the resolution model
326
328{
330 if (!conv) return 0 ;
331 return &conv->convVar() ;
332}
333
334
335
336////////////////////////////////////////////////////////////////////////////////
337/// Calculate the current unnormalized value of the PDF
338///
339/// PDF = sum_k coef_k * [ basis_k (x) ResModel ]
340///
341
343{
344 Double_t result(0) ;
345
347 RooAbsPdf* conv ;
348 Int_t index(0) ;
349 while(((conv=(RooAbsPdf*)_convSetIter->Next()))) {
350 Double_t coef = coefficient(index++) ;
351 if (coef!=0.) {
352 Double_t c = conv->getVal(0) ;
353 Double_t r = coef ;
354 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") val += coef*conv [" << index-1 << "/"
355 << _convSet.getSize() << "] coef = " << r << " conv = " << c << endl ;
356 result += conv->getVal(0)*coef ;
357 } else {
358 cxcoutD(Eval) << "RooAbsAnaConvPdf::evaluate(" << GetName() << ") [" << index-1 << "/" << _convSet.getSize() << "] coef = 0" << endl ;
359 }
360 }
361
362 return result ;
363}
364
365
366
367////////////////////////////////////////////////////////////////////////////////
368/// Advertise capability to perform (analytical) integrals
369/// internally. For a given integration request over allVars while
370/// normalized over normSet2 and in range 'rangeName', returns
371/// largest subset that can be performed internally in analVars
372/// Return code is unique integer code identifying integration scenario
373/// to be passed to analyticalIntegralWN() to calculate requeste integral
374///
375/// Class RooAbsAnaConv defers analytical integration request to
376/// resolution model and/or coefficient implementations and
377/// aggregates results into composite configuration with a unique
378/// code assigned by RooAICRegistry
379
381 RooArgSet& analVars, const RooArgSet* normSet2, const char* /*rangeName*/) const
382{
383 // Handle trivial no-integration scenario
384 if (allVars.getSize()==0) return 0 ;
385
386 if (_forceNumInt) return 0 ;
387
388 // Select subset of allVars that are actual dependents
389 RooArgSet* allDeps = getObservables(allVars) ;
390 RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
391
392 RooAbsArg *arg ;
393 RooResolutionModel *conv ;
394
395 RooArgSet* intSetAll = new RooArgSet(*allDeps,"intSetAll") ;
396
397 // Split intSetAll in coef/conv parts
398 RooArgSet* intCoefSet = new RooArgSet("intCoefSet") ;
399 RooArgSet* intConvSet = new RooArgSet("intConvSet") ;
400 TIterator* varIter = intSetAll->createIterator() ;
401 TIterator* convIter = _convSet.createIterator() ;
402
403 while(((arg=(RooAbsArg*) varIter->Next()))) {
404 Bool_t ok(kTRUE) ;
405 convIter->Reset() ;
406 while(((conv=(RooResolutionModel*) convIter->Next()))) {
407 if (conv->dependsOn(*arg)) ok=kFALSE ;
408 }
409
410 if (ok) {
411 intCoefSet->add(*arg) ;
412 } else {
413 intConvSet->add(*arg) ;
414 }
415
416 }
417 delete varIter ;
418
419
420 // Split normSetAll in coef/conv parts
421 RooArgSet* normCoefSet = new RooArgSet("normCoefSet") ;
422 RooArgSet* normConvSet = new RooArgSet("normConvSet") ;
423 RooArgSet* normSetAll = normSet ? (new RooArgSet(*normSet,"normSetAll")) : 0 ;
424 if (normSetAll) {
425 varIter = normSetAll->createIterator() ;
426 while(((arg=(RooAbsArg*) varIter->Next()))) {
427 Bool_t ok(kTRUE) ;
428 convIter->Reset() ;
429 while(((conv=(RooResolutionModel*) convIter->Next()))) {
430 if (conv->dependsOn(*arg)) ok=kFALSE ;
431 }
432
433 if (ok) {
434 normCoefSet->add(*arg) ;
435 } else {
436 normConvSet->add(*arg) ;
437 }
438
439 }
440 delete varIter ;
441 }
442 delete convIter ;
443
444 if (intCoefSet->getSize()==0) {
445 delete intCoefSet ; intCoefSet=0 ;
446 }
447 if (intConvSet->getSize()==0) {
448 delete intConvSet ; intConvSet=0 ;
449 }
450 if (normCoefSet->getSize()==0) {
451 delete normCoefSet ; normCoefSet=0 ;
452 }
453 if (normConvSet->getSize()==0) {
454 delete normConvSet ; normConvSet=0 ;
455 }
456
457
458
459 // Store integration configuration in registry
460 Int_t masterCode(0) ;
461 std::vector<Int_t> tmp(1, 0) ;
462
463 masterCode = _codeReg.store(tmp, intCoefSet, intConvSet, normCoefSet, normConvSet) + 1 ; // takes ownership of all sets
464
465 analVars.add(*allDeps) ;
466 delete allDeps ;
467 if (normSet) delete normSet ;
468 if (normSetAll) delete normSetAll ;
469 delete intSetAll ;
470
471// cout << this << "---> masterCode = " << masterCode << endl ;
472
473 return masterCode ;
474}
475
476
477
478
479////////////////////////////////////////////////////////////////////////////////
480/// Return analytical integral defined by given code, which is returned
481/// by getAnalyticalIntegralWN()
482///
483/// For unnormalized integrals the returned value is
484/// _ _
485/// PDF = sum_k Int(dx) coef_k * Int(dy) [ basis_k (x) ResModel ].
486/// _
487/// where x is the set of coefficient dependents to be integrated
488/// and y the set of basis function dependents to be integrated.
489///
490/// For normalized integrals this becomes
491///
492/// sum_k Int(dx) coef_k * Int(dy) [ basis_k (x) ResModel ].
493/// PDF = --------------------------------------------------------
494/// sum_k Int(dv) coef_k * Int(dw) [ basis_k (x) ResModel ].
495///
496/// where x is the set of coefficient dependents to be integrated,
497/// y the set of basis function dependents to be integrated,
498/// v is the set of coefficient dependents over which is normalized and
499/// w is the set of basis function dependents over which is normalized.
500///
501/// Set x must be contained in v and set y must be contained in w.
502///
503
504Double_t RooAbsAnaConvPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
505{
506 // WVE needs adaptation to handle new rangeName feature
507
508 // Handle trivial passthrough scenario
509 if (code==0) return getVal(normSet) ;
510
511 // Unpack master code
512 RooArgSet *intCoefSet, *intConvSet, *normCoefSet, *normConvSet ;
513 _codeReg.retrieve(code-1,intCoefSet,intConvSet,normCoefSet,normConvSet) ;
514
515 RooResolutionModel* conv ;
516 Int_t index(0) ;
517 Double_t answer(0) ;
519
520 if (normCoefSet==0&&normConvSet==0) {
521
522 // Integral over unnormalized function
523 Double_t integral(0) ;
524 const TNamed *_rangeName = RooNameReg::ptr(rangeName);
525 while(((conv=(RooResolutionModel*)_convSetIter->Next()))) {
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 while(((conv=(RooResolutionModel*)_convSetIter->Next()))) {
543
544 Double_t coefInt = getCoefNorm(index,intCoefSet,_rangeName) ;
545 //cout << "coefInt[" << index << "] = " << coefInt << "*" << term << " " << (intCoefSet?*intCoefSet:RooArgSet()) << endl ;
546 if (coefInt!=0) {
547 Double_t term = (_rangeName ? conv->getNormObj(0,intConvSet,_rangeName)->getVal() : conv->getNorm(intConvSet) ) ;
548 integral += coefInt*term ;
549 }
550
551 Double_t coefNorm = getCoefNorm(index,normCoefSet) ;
552 //cout << "coefNorm[" << index << "] = " << coefNorm << "*" << term << " " << (normCoefSet?*normCoefSet:RooArgSet()) << endl ;
553 if (coefNorm!=0) {
554 Double_t term = conv->getNorm(normConvSet) ;
555 norm += coefNorm*term ;
556 }
557
558 index++ ;
559 }
560 answer = integral/norm ;
561 }
562
563 return answer ;
564}
565
566
567
568////////////////////////////////////////////////////////////////////////////////
569/// Default implementation of function advertising integration capabilities. The interface is
570/// similar to that of getAnalyticalIntegral except that an integer code is added that
571/// designates the coefficient number for which the integration capabilities are requested
572///
573/// This default implementation advertises that no internal integrals are supported.
574
575Int_t RooAbsAnaConvPdf::getCoefAnalyticalIntegral(Int_t /* coef*/, RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
576{
577 return 0 ;
578}
579
580
581
582////////////////////////////////////////////////////////////////////////////////
583/// Default implementation of function implementing advertised integrals. Only
584/// the pass-through scenario (no integration) is implemented.
585
586Double_t RooAbsAnaConvPdf::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* /*rangeName*/) const
587{
588 if (code==0) return coefficient(coef) ;
589 coutE(InputArguments) << "RooAbsAnaConvPdf::coefAnalyticalIntegral(" << GetName() << ") ERROR: unrecognized integration code: " << code << endl ;
590 assert(0) ;
591 return 1 ;
592}
593
594
595
596////////////////////////////////////////////////////////////////////////////////
597/// This function forces RooRealIntegral to offer all integration dependents
598/// to RooAbsAnaConvPdf::getAnalyticalIntegralWN() for consideration for
599/// internal integration, if RooRealIntegral considers this to be unsafe (e.g. due
600/// to hidden Jacobian terms).
601///
602/// RooAbsAnaConvPdf will not attempt to actually integrate all these dependents
603/// but feed them to the resolution models integration interface, which will
604/// make the final determination on how to integrate these dependents.
605
607{
608 return kTRUE ;
609}
610
611
612
613////////////////////////////////////////////////////////////////////////////////
614/// Returns the normalization integral value of the coefficient with number coefIdx over normalization
615/// set nset in range rangeName
616
617Double_t RooAbsAnaConvPdf::getCoefNorm(Int_t coefIdx, const RooArgSet* nset, const TNamed* rangeName) const
618{
619 if (nset==0) return coefficient(coefIdx) ;
620
621 CacheElem* cache = (CacheElem*) _coefNormMgr.getObj(nset,0,0,rangeName) ;
622 if (!cache) {
623
624 cache = new CacheElem ;
625
626 // Make list of coefficient normalizations
627 Int_t i ;
629
630 for (i=0 ; i<cache->_coefVarList.getSize() ; i++) {
631 RooAbsReal* coefInt = static_cast<RooAbsReal&>(*cache->_coefVarList.at(i)).createIntegral(*nset,RooNameReg::str(rangeName)) ;
632 cache->_normList.addOwned(*coefInt) ;
633 }
634
635 _coefNormMgr.setObj(nset,0,cache,rangeName) ;
636 }
637
638 return ((RooAbsReal*)cache->_normList.at(coefIdx))->getVal() ;
639}
640
641
642
643////////////////////////////////////////////////////////////////////////////////
644/// Build complete list of coefficient variables
645
647{
648 // Instantate a coefficient variables
649 for (Int_t i=0 ; i<_convSet.getSize() ; i++) {
650 RooArgSet* cvars = coefVars(i) ;
651 RooAbsReal* coefVar = new RooConvCoefVar(Form("%s_coefVar_%d",GetName(),i),"coefVar",*this,i,cvars) ;
652 varList.addOwned(*coefVar) ;
653 delete cvars ;
654 }
655
656}
657
658
659////////////////////////////////////////////////////////////////////////////////
660/// Return set of parameters with are used exclusively by the coefficient functions
661
663{
664 RooArgSet* cVars = getParameters((RooArgSet*)0) ;
665 TIterator* iter = cVars->createIterator() ;
666 RooAbsArg* arg ;
667 Int_t i ;
668 while(((arg=(RooAbsArg*)iter->Next()))) {
669 for (i=0 ; i<_convSet.getSize() ; i++) {
670 if (_convSet.at(i)->dependsOn(*arg)) {
671 cVars->remove(*arg,kTRUE) ;
672 }
673 }
674 }
675 delete iter ;
676 return cVars ;
677}
678
679
680
681
682////////////////////////////////////////////////////////////////////////////////
683/// Print info about this object to the specified stream. In addition to the info
684/// from RooAbsPdf::printStream() we add:
685///
686/// Verbose : detailed information on convolution integrals
687
688void RooAbsAnaConvPdf::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
689{
690 RooAbsPdf::printMultiline(os,contents,verbose,indent);
691
692 os << indent << "--- RooAbsAnaConvPdf ---" << endl;
694 RooResolutionModel* conv ;
695 while (((conv=(RooResolutionModel*)iter->Next()))) {
696 conv->printMultiline(os,contents,verbose,indent) ;
697 }
698}
699
700
701///////////////////////////////////////////////////////////////////////////////
702/// Label OK'ed components with cache-and-track
704{
705 RooFIter citer = _convSet.fwdIterator() ;
706 RooAbsArg* carg ;
707 while ((carg=citer.next())) {
708 if (carg->canNodeBeCached()==Always) {
709 trackNodes.add(*carg) ;
710 //cout << "tracking node RooAddPdf component " << carg->IsA()->GetName() << "::" << carg->GetName() << endl ;
711 }
712 }
713}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define c(i)
Definition: RSha256.hxx:101
static RooMathCoreReg dummy
#define coutI(a)
Definition: RooMsgService.h:31
#define cxcoutD(a)
Definition: RooMsgService.h:79
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
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...
Iterator over _convNormSet.
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.
const RooRealVar * convVar() const
Return a pointer to the convolution variable instance used in the resolution model.
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
TIterator * _convSetIter
Subset of last normalization that applies to convolutions.
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.
RooArgList _basisList
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
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.
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.
RooArgSet * _convNormSet
List of created basis functions.
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 (of arbitrary type) an...
Definition: RooAbsArg.h:66
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
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.
Definition: RooAbsArg.cxx:736
virtual CacheMode canNodeBeCached() const
Definition: RooAbsArg.h:317
friend class RooArgSet
Definition: RooAbsArg.h:471
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:533
void setAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
Definition: RooAbsArg.cxx:241
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Change cache operation mode to given mode.
Definition: RooAbsArg.cxx:1746
OperMode operMode() const
Definition: RooAbsArg.h:399
Int_t getSize() const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
RooFIter fwdIterator() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
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
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
friend class CacheElem
Definition: RooAbsPdf.h:328
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...
Definition: RooAbsPdf.cxx:2047
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print multi line detailed information of this RooAbsPdf.
Definition: RooAbsPdf.cxx:1639
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=0) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
Definition: RooAbsPdf.cxx:401
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Definition: RooAbsPdf.cxx:2082
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Bool_t _forceNumInt
Definition: RooAbsReal.h:394
Double_t getVal(const RooArgSet *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
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 ...
Definition: RooAbsReal.cxx:502
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
Int_t index(const RooAbsArg *arg) const
Definition: RooArgList.h:76
RooAbsArg * absArg() const
Definition: RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
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:31
RooAbsArg * next()
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
Definition: RooGenContext.h:30
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Reimplementation of standard RooArgList::remove()
virtual void removeAll()
Reimplementation of standard RooArgList::removeAll()
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::addOwned()
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition: RooNameReg.cxx:135
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:125
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
virtual Bool_t setArg(RooAbsReal &newRef)
Change object held in proxy into newRef.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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.
RooRealVar & convVar() const
Return the convolution variable of the resolution model.
Double_t getNorm(const RooArgSet *nset=0) const
Return the integral of this PDF over all elements of 'nset'.
const RooFormulaVar & basis() const
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
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
@ Generation
Definition: RooGlobalFunc.h:57
@ InputArguments
Definition: RooGlobalFunc.h:58
STL namespace.