Logo ROOT   6.16/01
Reference Guide
RooProdPdf.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\file RooProdPdf.cxx
19\class RooProdPdf
20\ingroup Roofitcore
21
22RooProdPdf is an efficient implementation of a product of PDFs of the form
23
24\f[PDF_1 * PDF_2 * ... * PDF_N \f]
25
26PDFs may share observables. If that is the case any irreducable subset
27of PDFS that share observables will be normalized with explicit numeric
28integration as any built-in normalization will no longer be valid.
29Alternatively, products using conditional PDFs can be defined, e.g.
30
31\f[F(x|y) * G(y) \f]
32
33meaning a pdf F(x) _given_ y and a PDF G(y). In this contruction F is only
34normalized w.r.t x and G is normalized w.r.t y. The product in this construction
35is properly normalized.
36If exactly one of the component PDFs supports extended likelihood fits, the
37product will also be usable in extended mode, returning the number of expected
38events from the extendable component PDF. The extendable component does not
39have to appear in any specific place in the list.
40
41**/
42
43#include "RooProdPdf.h"
44#include "RooRealProxy.h"
45#include "RooProdGenContext.h"
46#include "RooGenProdProj.h"
47#include "RooProduct.h"
48#include "RooNameReg.h"
49#include "RooMsgService.h"
50#include "RooFormulaVar.h"
51#include "RooRealVar.h"
52#include "RooAddition.h"
53#include "RooGlobalFunc.h"
54#include "RooConstVar.h"
55#include "RooWorkspace.h"
56#include "RooRangeBoolean.h"
57#include "RooCustomizer.h"
58#include "RooRealIntegral.h"
59#include "RooTrace.h"
60
61#include <cstring>
62#include <sstream>
63#include <algorithm>
64
65#ifndef _WIN32
66#include <strings.h>
67#endif
68
69
70#include "TSystem.h"
71
72using namespace std;
73
75;
76
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// Default constructor
81
83 _curNormSet(0),
84 _cutOff(0),
85 _extendedIndex(-1),
86 _useDefaultGen(kFALSE),
87 _refRangeName(0),
88 _selfNorm(kTRUE)
89{
90 // Default constructor
92}
93
94
95
96////////////////////////////////////////////////////////////////////////////////
97/// Dummy constructor
98
99RooProdPdf::RooProdPdf(const char *name, const char *title, Double_t cutOff) :
100 RooAbsPdf(name,title),
101 _cacheMgr(this,10),
102 _genCode(10),
103 _cutOff(cutOff),
104 _pdfList("!pdfs","List of PDFs",this),
105 _extendedIndex(-1),
106 _useDefaultGen(kFALSE),
107 _refRangeName(0),
108 _selfNorm(kTRUE)
109{
111}
112
113
114
115////////////////////////////////////////////////////////////////////////////////
116/// Constructor with 2 PDFs (most frequent use case).
117///
118/// The optional cutOff parameter can be used as a speed optimization if
119/// one or more of the PDF have sizable regions with very small values,
120/// which would pull the entire product of PDFs to zero in those regions.
121///
122/// After each PDF multiplication, the running product is compared with
123/// the cutOff parameter. If the running product is smaller than the
124/// cutOff value, the product series is terminated and remaining PDFs
125/// are not evaluated.
126///
127/// There is no magic value of the cutOff, the user should experiment
128/// to find the appropriate balance between speed and precision.
129/// If a cutoff is specified, the PDFs most likely to be small should
130/// be put first in the product. The default cutOff value is zero.
131///
132
133RooProdPdf::RooProdPdf(const char *name, const char *title,
134 RooAbsPdf& pdf1, RooAbsPdf& pdf2, Double_t cutOff) :
135 RooAbsPdf(name,title),
136 _cacheMgr(this,10),
137 _genCode(10),
138 _cutOff(cutOff),
139 _pdfList("!pdfs","List of PDFs",this),
140 _extendedIndex(-1),
141 _useDefaultGen(kFALSE),
142 _refRangeName(0),
143 _selfNorm(kTRUE)
144{
145 _pdfList.add(pdf1) ;
146 RooArgSet* nset1 = new RooArgSet("nset") ;
147 _pdfNSetList.Add(nset1) ;
148 if (pdf1.canBeExtended()) {
149 _extendedIndex = _pdfList.index(&pdf1) ;
150 }
151
152 _pdfList.add(pdf2) ;
153 RooArgSet* nset2 = new RooArgSet("nset") ;
154 _pdfNSetList.Add(nset2) ;
155
156 if (pdf2.canBeExtended()) {
157 if (_extendedIndex>=0) {
158 // Protect against multiple extended terms
159 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName()
160 << ") multiple components with extended terms detected,"
161 << " product will not be extendible." << endl ;
162 _extendedIndex=-1 ;
163 } else {
165 }
166 }
168}
169
170
171
172////////////////////////////////////////////////////////////////////////////////
173/// Constructor from a list of PDFs
174///
175/// The optional cutOff parameter can be used as a speed optimization if
176/// one or more of the PDF have sizable regions with very small values,
177/// which would pull the entire product of PDFs to zero in those regions.
178///
179/// After each PDF multiplication, the running product is compared with
180/// the cutOff parameter. If the running product is smaller than the
181/// cutOff value, the product series is terminated and remaining PDFs
182/// are not evaluated.
183///
184/// There is no magic value of the cutOff, the user should experiment
185/// to find the appropriate balance between speed and precision.
186/// If a cutoff is specified, the PDFs most likely to be small should
187/// be put first in the product. The default cutOff value is zero.
188
189RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgList& inPdfList, Double_t cutOff) :
190 RooAbsPdf(name,title),
191 _cacheMgr(this,10),
192 _genCode(10),
193 _cutOff(cutOff),
194 _pdfList("!pdfs","List of PDFs",this),
195 _extendedIndex(-1),
196 _useDefaultGen(kFALSE),
197 _refRangeName(0),
198 _selfNorm(kTRUE)
199{
200 RooFIter iter = inPdfList.fwdIterator();
201 RooAbsArg* arg ;
202 Int_t numExtended(0) ;
203 while((arg=(RooAbsArg*)iter.next())) {
204 RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg) ;
205 if (!pdf) {
206 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName() << ") list arg "
207 << arg->GetName() << " is not a PDF, ignored" << endl ;
208 continue ;
209 }
210 _pdfList.add(*pdf) ;
211
212 RooArgSet* nset = new RooArgSet("nset") ;
213 _pdfNSetList.Add(nset) ;
214
215 if (pdf->canBeExtended()) {
217 numExtended++ ;
218 }
219 }
220
221 // Protect against multiple extended terms
222 if (numExtended>1) {
223 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName()
224 << ") WARNING: multiple components with extended terms detected,"
225 << " product will not be extendible." << endl ;
226 _extendedIndex = -1 ;
227 }
228
230}
231
232
233
234////////////////////////////////////////////////////////////////////////////////
235/// Constructor from named argument list
236///
237/// fullPdf -- Set of 'regular' PDFS that are normalized over all their observables
238/// Argument | Description
239/// ---------------------------|------------
240/// Conditional(pdfSet,depSet) | Add PDF to product with condition that it<br>only be normalized over specified observables<br>any remaining observables will be conditional<br>observables
241///
242///
243/// For example, given a PDF F(x,y) and G(y)
244///
245/// RooProdPdf("P","P",G,Conditional(F,x)) will construct a 2-dimensional PDF as follows:
246///
247/// \f[P(x,y) = \frac{G(y)}{\int_y G(y)} * \frac{F(x,y)}{\int_x G(x,y)} \f]
248///
249/// which is a well normalized and properly defined PDF, but different from the
250///
251/// \f[P'(x,y) = \frac{F(x,y)}{\int_y G(y)} / \int_x\int_y F(x,y)*G(y) \f]
252///
253/// In the former case the y distribution of P is identical to that of G, while
254/// F only is used to determine the correlation between X and Y. In the latter
255/// case the Y distribution is defined by the product of F and G.
256///
257/// This P(x,y) construction is analoguous to generating events from F(x,y) with
258/// a prototype dataset sampled from G(y)
259
260RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgSet& fullPdfSet,
261 const RooCmdArg& arg1, const RooCmdArg& arg2,
262 const RooCmdArg& arg3, const RooCmdArg& arg4,
263 const RooCmdArg& arg5, const RooCmdArg& arg6,
264 const RooCmdArg& arg7, const RooCmdArg& arg8) :
265 RooAbsPdf(name,title),
266 _cacheMgr(this,10),
267 _genCode(10),
268 _cutOff(0),
269 _pdfList("!pdfs","List of PDFs",this),
270 _extendedIndex(-1),
271 _useDefaultGen(kFALSE),
272 _refRangeName(0),
273 _selfNorm(kTRUE)
274{
276 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
277 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
278 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
279 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
280
281 initializeFromCmdArgList(fullPdfSet,l) ;
283}
284
285
286
287////////////////////////////////////////////////////////////////////////////////
288/// Constructor from named argument list
289
290RooProdPdf::RooProdPdf(const char* name, const char* title,
291 const RooCmdArg& arg1, const RooCmdArg& arg2,
292 const RooCmdArg& arg3, const RooCmdArg& arg4,
293 const RooCmdArg& arg5, const RooCmdArg& arg6,
294 const RooCmdArg& arg7, const RooCmdArg& arg8) :
295 RooAbsPdf(name,title),
296 _cacheMgr(this,10),
297 _genCode(10),
298 _cutOff(0),
299 _pdfList("!pdfList","List of PDFs",this),
300 _extendedIndex(-1),
301 _useDefaultGen(kFALSE),
302 _refRangeName(0),
303 _selfNorm(kTRUE)
304{
306 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
307 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
308 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
309 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
310
313}
314
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// Internal constructor from list of named arguments
319
320RooProdPdf::RooProdPdf(const char* name, const char* title, const RooArgSet& fullPdfSet, const RooLinkedList& cmdArgList) :
321 RooAbsPdf(name,title),
322 _cacheMgr(this,10),
323 _genCode(10),
324 _cutOff(0),
325 _pdfList("!pdfs","List of PDFs",this),
326 _extendedIndex(-1),
327 _useDefaultGen(kFALSE),
328 _refRangeName(0),
329 _selfNorm(kTRUE)
330{
331 initializeFromCmdArgList(fullPdfSet, cmdArgList) ;
333}
334
335
336
337////////////////////////////////////////////////////////////////////////////////
338/// Copy constructor
339
340RooProdPdf::RooProdPdf(const RooProdPdf& other, const char* name) :
341 RooAbsPdf(other,name),
342 _cacheMgr(other._cacheMgr,this),
343 _genCode(other._genCode),
344 _cutOff(other._cutOff),
345 _pdfList("!pdfs",this,other._pdfList),
346 _extendedIndex(other._extendedIndex),
347 _useDefaultGen(other._useDefaultGen),
348 _refRangeName(other._refRangeName),
349 _selfNorm(other._selfNorm),
350 _defNormSet(other._defNormSet)
351{
352 // Clone contents of normalizarion set list
353 RooFIter iter = other._pdfNSetList.fwdIterator();
354 RooArgSet* nset ;
355 while((nset=(RooArgSet*)iter.next())) {
356 RooArgSet* tmp = (RooArgSet*) nset->snapshot() ;
357 tmp->setName(nset->GetName()) ;
358 _pdfNSetList.Add(tmp) ;
359 }
361}
362
363
364
365////////////////////////////////////////////////////////////////////////////////
366/// Initialize RooProdPdf configuration from given list of RooCmdArg configuration arguments
367/// and set of 'regular' p.d.f.s in product
368
370{
371 Int_t numExtended(0) ;
372
373 // Process set of full PDFS
374 RooFIter siter = fullPdfSet.fwdIterator() ;
375 RooAbsPdf* pdf ;
376 while((pdf=(RooAbsPdf*)siter.next())) {
377 _pdfList.add(*pdf) ;
378 RooArgSet* nset1 = new RooArgSet("nset") ;
379 _pdfNSetList.Add(nset1) ;
380
381 if (pdf->canBeExtended()) {
383 numExtended++ ;
384 }
385
386 }
387
388 // Process list of conditional PDFs
389 RooFIter iter = l.fwdIterator();
390 RooCmdArg* carg ;
391 while((carg=(RooCmdArg*)iter.next())) {
392
393 if (0 == strcmp(carg->GetName(), "Conditional")) {
394
395 Int_t argType = carg->getInt(0) ;
396 RooArgSet* pdfSet = (RooArgSet*) carg->getSet(0) ;
397 RooArgSet* normSet = (RooArgSet*) carg->getSet(1) ;
398
399 RooFIter siter2 = pdfSet->fwdIterator() ;
400 RooAbsPdf* thePdf ;
401 while ((thePdf=(RooAbsPdf*)siter2.next())) {
402 _pdfList.add(*thePdf) ;
403
404 RooArgSet* tmp = (RooArgSet*) normSet->snapshot() ;
405 tmp->setName(0 == argType ? "nset" : "cset") ;
406 _pdfNSetList.Add(tmp) ;
407
408 if (thePdf->canBeExtended()) {
409 _extendedIndex = _pdfList.index(thePdf) ;
410 numExtended++ ;
411 }
412
413 }
414
415 } else if (0 != strlen(carg->GetName())) {
416 coutW(InputArguments) << "Unknown arg: " << carg->GetName() << endl ;
417 }
418 }
419
420 // Protect against multiple extended terms
421 if (numExtended>1) {
422 coutW(InputArguments) << "RooProdPdf::RooProdPdf(" << GetName()
423 << ") WARNING: multiple components with extended terms detected,"
424 << " product will not be extendible." << endl ;
425 _extendedIndex = -1 ;
426 }
427
428
429}
430
431
432
433////////////////////////////////////////////////////////////////////////////////
434/// Destructor
435
437{
440}
441
442
443
444////////////////////////////////////////////////////////////////////////////////
445/// Overload getVal() to intercept normalization set for use in evaluate()
446
448{
449 _curNormSet = (RooArgSet*)set ;
450 return RooAbsPdf::getValV(set) ;
451}
452
453
454
455////////////////////////////////////////////////////////////////////////////////
456/// Calculate current value of object
457
459{
460 Int_t code ;
461 CacheElem* cache = (CacheElem*) _cacheMgr.getObj(_curNormSet,0,&code) ;
462
463 // If cache doesn't have our configuration, recalculate here
464 if (!cache) {
465 RooArgList *plist(0) ;
466 RooLinkedList *nlist(0) ;
467 getPartIntList(_curNormSet,0,plist,nlist,code) ;
468 cache = (CacheElem*) _cacheMgr.getObj(_curNormSet,0,&code) ;
469 }
470
471
472 return calculate(*cache) ;
473}
474
475
476
477////////////////////////////////////////////////////////////////////////////////
478/// Calculate running product of pdfs terms, using the supplied
479/// normalization set in 'normSetList' for each component
480
481Double_t RooProdPdf::calculate(const RooArgList* partIntList, const RooLinkedList* normSetList) const
482{
483 RooAbsReal* partInt;
484 RooArgSet* normSet;
485 Double_t value = 1.0;
486 RooFIter plIter = partIntList->fwdIterator(), nlIter = normSetList->fwdIterator();
487
488 for (partInt = (RooAbsReal*) plIter.next(),
489 normSet = (RooArgSet*) nlIter.next(); partInt && normSet;
490 partInt = (RooAbsReal*) plIter.next(),
491 normSet = (RooArgSet*) nlIter.next()) {
492 const Double_t piVal = partInt->getVal(normSet->getSize() > 0 ? normSet : 0);
493// cout << "RooProdPdf::calculate(" << GetName() << ") partInt(" << partInt->GetName() << ") = " << piVal << " normSet = " << normSet << " " << (normSet->getSize()>0 ? *normSet : RooArgSet()) << endl ;
494 value *= piVal;
495 if (value <= _cutOff) { break; }
496 }
497
498 return value ;
499}
500
501
502
503////////////////////////////////////////////////////////////////////////////////
504/// Calculate running product of pdfs terms, using the supplied
505/// normalization set in 'normSetList' for each component
506
508{
509 //cout << "RooProdPdf::calculate from cache" << endl ;
510
511 if (cache._isRearranged) {
512 if (dologD(Eval)) {
513 cxcoutD(Eval) << "RooProdPdf::calculate(" << GetName() << ") rearranged product calculation"
514 << " calculate: num = " << cache._rearrangedNum->GetName() << " = " << cache._rearrangedNum->getVal() << endl ;
515// cache._rearrangedNum->printComponentTree("",0,5) ;
516 cxcoutD(Eval) << "calculate: den = " << cache._rearrangedDen->GetName() << " = " << cache._rearrangedDen->getVal() << endl ;
517// cache._rearrangedDen->printComponentTree("",0,5) ;
518 }
519
520 return cache._rearrangedNum->getVal() / cache._rearrangedDen->getVal();
521 } else {
522 RooAbsReal* partInt;
523 RooArgSet* normSet;
524 RooFIter plIter = cache._partList.fwdIterator();
525 RooFIter nlIter = cache._normList.fwdIterator();
526 Double_t value = 1.0;
527 for (partInt = (RooAbsReal*) plIter.next(),
528 normSet = (RooArgSet*) nlIter.next(); partInt && normSet;
529 partInt = (RooAbsReal*) plIter.next(),
530 normSet = (RooArgSet*) nlIter.next()) {
531 const Double_t piVal = partInt->getVal(normSet->getSize() > 0 ? normSet : 0);
532 value *= piVal ;
533 if (value <= _cutOff) break;
534 }
535// cout << "return value = " << value << endl ;
536 return value ;
537 }
538}
539
540
541
542////////////////////////////////////////////////////////////////////////////////
543/// Factorize product in irreducible terms for given choice of integration/normalization
544
545void RooProdPdf::factorizeProduct(const RooArgSet& normSet, const RooArgSet& intSet,
546 RooLinkedList& termList, RooLinkedList& normList,
547 RooLinkedList& impDepList, RooLinkedList& crossDepList,
548 RooLinkedList& intList) const
549{
550 // List of all term dependents: normalization and imported
551 RooLinkedList depAllList;
552 RooLinkedList depIntNoNormList;
553
554 // Setup lists for factorization terms and their dependents
555 RooArgSet* term(0);
556 RooArgSet* termNormDeps(0);
557 RooArgSet* termAllDeps(0);
558 RooArgSet* termIntDeps(0);
559 RooArgSet* termIntNoNormDeps(0);
560
561 // Loop over the PDFs
562 RooAbsPdf* pdf;
563 RooArgSet* pdfNSetOrig;
564 for (RooFIter pdfIter = _pdfList.fwdIterator(),
565 nIter = _pdfNSetList.fwdIterator();
566 (pdfNSetOrig = (RooArgSet*) nIter.next(),
567 pdf = (RooAbsPdf*) pdfIter.next()); ) {
568 RooArgSet* pdfNSet, *pdfCSet;
569
570 // Reduce pdfNSet to actual dependents
571 if (0 == strcmp("nset", pdfNSetOrig->GetName())) {
572 pdfNSet = pdf->getObservables(*pdfNSetOrig);
573 pdfCSet = new RooArgSet;
574 } else if (0 == strcmp("cset", pdfNSetOrig->GetName())) {
575 RooArgSet* tmp = pdf->getObservables(normSet);
576 tmp->remove(*pdfNSetOrig, kTRUE, kTRUE);
577 pdfCSet = pdfNSetOrig;
578 pdfNSet = tmp;
579 } else {
580 // Legacy mode. Interpret at NSet for backward compatibility
581 pdfNSet = pdf->getObservables(*pdfNSetOrig);
582 pdfCSet = new RooArgSet;
583 }
584
585
586 RooArgSet pdfNormDeps; // Dependents to be normalized for the PDF
587 RooArgSet pdfAllDeps; // All dependents of this PDF
588
589 // Make list of all dependents of this PDF
590 RooArgSet* tmp = pdf->getObservables(normSet);
591 pdfAllDeps.add(*tmp);
592 delete tmp;
593
594
595// cout << GetName() << ": pdf = " << pdf->GetName() << " pdfAllDeps = " << pdfAllDeps << " pdfNSet = " << *pdfNSet << " pdfCSet = " << *pdfCSet << endl;
596
597 // Make list of normalization dependents for this PDF;
598 if (pdfNSet->getSize() > 0) {
599 // PDF is conditional
600 RooArgSet* tmp2 = (RooArgSet*) pdfAllDeps.selectCommon(*pdfNSet);
601 pdfNormDeps.add(*tmp2);
602 delete tmp2;
603 } else {
604 // PDF is regular
605 pdfNormDeps.add(pdfAllDeps);
606 }
607
608// cout << GetName() << ": pdfNormDeps for " << pdf->GetName() << " = " << pdfNormDeps << endl;
609
610 RooArgSet* pdfIntSet = pdf->getObservables(intSet) ;
611
612 // WVE if we have no norm deps, conditional observables should be taken out of pdfIntSet
613 if (0 == pdfNormDeps.getSize() && pdfCSet->getSize() > 0) {
614 pdfIntSet->remove(*pdfCSet, kTRUE, kTRUE);
615// cout << GetName() << ": have no norm deps, removing conditional observables from intset" << endl;
616 }
617
618 RooArgSet pdfIntNoNormDeps(*pdfIntSet);
619 pdfIntNoNormDeps.remove(pdfNormDeps, kTRUE, kTRUE);
620
621// cout << GetName() << ": pdf = " << pdf->GetName() << " intset = " << *pdfIntSet << " pdfIntNoNormDeps = " << pdfIntNoNormDeps << endl;
622
623 // Check if this PDF has dependents overlapping with one of the existing terms
624 Bool_t done(kFALSE);
625 for (RooFIter lIter = termList.fwdIterator(),
626 ldIter = normList.fwdIterator(),
627 laIter = depAllList.fwdIterator();
628 (termNormDeps = (RooArgSet*) ldIter.next(),
629 termAllDeps = (RooArgSet*) laIter.next(),
630 term = (RooArgSet*) lIter.next()); ) {
631 // PDF should be added to existing term if
632 // 1) It has overlapping normalization dependents with any other PDF in existing term
633 // 2) It has overlapping dependents of any class for which integration is requested
634 // 3) If normalization happens over multiple ranges, and those ranges are both defined
635 // in either observable
636
637 Bool_t normOverlap = pdfNormDeps.overlaps(*termNormDeps);
638 //Bool_t intOverlap = pdfIntSet->overlaps(*termAllDeps);
639
640 if (normOverlap) {
641// cout << GetName() << ": this term overlaps with term " << (*term) << " in normalization observables" << endl;
642
643 term->add(*pdf);
644 termNormDeps->add(pdfNormDeps, kFALSE);
645 termAllDeps->add(pdfAllDeps, kFALSE);
646 if (!termIntDeps) {
647 termIntDeps = new RooArgSet("termIntDeps");
648 }
649 termIntDeps->add(*pdfIntSet, kFALSE);
650 if (!termIntNoNormDeps) {
651 termIntNoNormDeps = new RooArgSet("termIntNoNormDeps");
652 }
653 termIntNoNormDeps->add(pdfIntNoNormDeps, kFALSE);
654 done = kTRUE;
655 }
656 }
657
658 // If not, create a new term
659 if (!done) {
660 if (!(0 == pdfNormDeps.getSize() && 0 == pdfAllDeps.getSize() &&
661 0 == pdfIntSet->getSize()) || 0 == normSet.getSize()) {
662// cout << GetName() << ": creating new term" << endl;
663 term = new RooArgSet("term");
664 termNormDeps = new RooArgSet("termNormDeps");
665 termAllDeps = new RooArgSet("termAllDeps");
666 termIntDeps = new RooArgSet("termIntDeps");
667 termIntNoNormDeps = new RooArgSet("termIntNoNormDeps");
668
669 term->add(*pdf);
670 termNormDeps->add(pdfNormDeps, kFALSE);
671 termAllDeps->add(pdfAllDeps, kFALSE);
672 termIntDeps->add(*pdfIntSet, kFALSE);
673 termIntNoNormDeps->add(pdfIntNoNormDeps, kFALSE);
674
675 termList.Add(term);
676 normList.Add(termNormDeps);
677 depAllList.Add(termAllDeps);
678 intList.Add(termIntDeps);
679 depIntNoNormList.Add(termIntNoNormDeps);
680 }
681 }
682
683 // We own the reduced version of pdfNSet
684 delete pdfNSet;
685 delete pdfIntSet;
686 if (pdfCSet != pdfNSetOrig) {
687 delete pdfCSet;
688 }
689 }
690
691 // Loop over list of terms again to determine 'imported' observables
692 RooArgSet *normDeps, *allDeps, *intNoNormDeps;
693 for (RooFIter lIter = termList.fwdIterator(),
694 ldIter = normList.fwdIterator(),
695 laIter = depAllList.fwdIterator(),
696 innIter = depIntNoNormList.fwdIterator();
697 (normDeps = (RooArgSet*) ldIter.next(),
698 allDeps = (RooArgSet*) laIter.next(),
699 intNoNormDeps = (RooArgSet*) innIter.next(),
700 term=(RooArgSet*)lIter.next()); ) {
701 // Make list of wholly imported dependents
702 RooArgSet impDeps(*allDeps);
703 impDeps.remove(*normDeps, kTRUE, kTRUE);
704 impDepList.Add(impDeps.snapshot());
705// cout << GetName() << ": list of imported dependents for term " << (*term) << " set to " << impDeps << endl ;
706
707 // Make list of cross dependents (term is self contained for these dependents,
708 // but components import dependents from other components)
709 RooArgSet* crossDeps = (RooArgSet*) intNoNormDeps->selectCommon(*normDeps);
710 crossDepList.Add(crossDeps->snapshot());
711// cout << GetName() << ": list of cross dependents for term " << (*term) << " set to " << *crossDeps << endl ;
712 delete crossDeps;
713 }
714
715 depAllList.Delete();
716 depIntNoNormList.Delete();
717
718 return;
719}
720
721
722
723
724////////////////////////////////////////////////////////////////////////////////
725/// Return list of (partial) integrals of product terms for integration
726/// of p.d.f over observables iset while normalization over observables nset.
727/// Also return list of normalization sets to be used to evaluate
728/// each component in the list correctly.
729
730void RooProdPdf::getPartIntList(const RooArgSet* nset, const RooArgSet* iset,
731 pRooArgList& partList, pRooLinkedList& nsetList,
732 Int_t& code, const char* isetRangeName) const
733{
734// cout << " FOLKERT::RooProdPdf::getPartIntList(" << GetName() <<") nset = " << (nset?*nset:RooArgSet()) << endl
735// << " _normRange = " << _normRange << endl
736// << " iset = " << (iset?*iset:RooArgSet()) << endl
737// << " isetRangeName = " << (isetRangeName?isetRangeName:"<null>") << endl ;
738
739 // Check if this configuration was created before
740 Int_t sterileIdx(-1);
741
742 CacheElem* cache = (CacheElem*) _cacheMgr.getObj(nset,iset,&sterileIdx,isetRangeName);
743 if (cache) {
744 code = _cacheMgr.lastIndex();
745 partList = &cache->_partList;
746 nsetList = &cache->_normList;
747 return;
748 }
749
750 // Create containers for partial integral components to be generated
751 cache = new CacheElem;
752
753 // Factorize the product in irreducible terms for this nset
754 RooLinkedList terms, norms, imp, ints, cross;
755// cout << "RooProdPdf::getPIL -- now calling factorizeProduct()" << endl ;
756
757
758 // Normalization set used for factorization
759 RooArgSet factNset(nset ? (*nset) : _defNormSet);
760// cout << GetName() << "factNset = " << factNset << endl ;
761
762 factorizeProduct(factNset, iset ? (*iset) : RooArgSet(), terms, norms, imp, cross, ints);
763
764 RooArgSet *norm, *integ, *xdeps, *imps;
765
766 // Group irriducible terms that need to be (partially) integrated together
767 RooLinkedList groupedList;
768 RooArgSet outerIntDeps;
769// cout << "RooProdPdf::getPIL -- now calling groupProductTerms()" << endl;
770 groupProductTerms(groupedList, outerIntDeps, terms, norms, imp, ints, cross);
771 RooFIter gIter = groupedList.fwdIterator();
773
774 // Loop over groups
775// cout<<"FK: pdf("<<GetName()<<") Starting selecting F(x|y)!"<<endl;
776 // Find groups of type F(x|y), i.e. termImpSet!=0, construct ratio object
777 map<string, RooArgSet> ratioTerms;
778 while ((group = (RooLinkedList*) gIter.next())) {
779 if (1 == group->GetSize()) {
780// cout<<"FK: Starting Single Term"<<endl;
781
782 RooArgSet* term = (RooArgSet*) group->At(0);
783
784 Int_t termIdx = terms.IndexOf(term);
785 norm=(RooArgSet*) norms.At(termIdx);
786 imps=(RooArgSet*)imp.At(termIdx);
787 RooArgSet termNSet(*norm), termImpSet(*imps);
788
789// cout<<"FK: termImpSet.getSize() = "<<termImpSet.getSize()<< " " << termImpSet << endl;
790// cout<<"FK: _refRangeName = "<<_refRangeName<<endl;
791
792 if (termImpSet.getSize() > 0 && 0 != _refRangeName) {
793
794// cout << "WVE now here" << endl;
795
796 // WVE we can skip this if the ref range is equal to the normalization range
797 Bool_t rangeIdentical(kTRUE);
798 RooFIter niter = termNSet.fwdIterator();
799 RooRealVar* normObs;
800// cout << "_normRange = " << _normRange << " _refRangeName = " << RooNameReg::str(_refRangeName) << endl ;
801 while ((normObs = (RooRealVar*) niter.next())) {
802 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
803 if (_normRange.Length() > 0) {
804 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
805 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
806 }
807 else{
808 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
809 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
810 }
811 }
812// cout<<"FK: rangeIdentical Single = "<<(rangeIdentical ? 'T':'F')<<endl;
813 // coverity[CONSTANT_EXPRESSION_RESULT]
814 if (!rangeIdentical || 1) {
815// cout << "PREPARING RATIO HERE (SINGLE TERM)" << endl ;
816 RooAbsReal* ratio = makeCondPdfRatioCorr(*(RooAbsReal*)term->first(), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
817 ostringstream str; termImpSet.printValue(str);
818// cout << GetName() << "inserting ratio term" << endl;
819 ratioTerms[str.str()].add(*ratio);
820 }
821 }
822
823 } else {
824// cout<<"FK: Starting Composite Term"<<endl;
825
826 RooArgSet compTermSet, compTermNorm;
827 RooFIter tIter = group->fwdIterator();
828 RooArgSet* term;
829 while ((term = (RooArgSet*) tIter.next())) {
830
831 Int_t termIdx = terms.IndexOf(term);
832 norm=(RooArgSet*) norms.At(termIdx);
833 imps=(RooArgSet*)imp.At(termIdx);
834 RooArgSet termNSet(*norm), termImpSet(*imps);
835
836 if (termImpSet.getSize() > 0 && 0 != _refRangeName) {
837
838 // WVE we can skip this if the ref range is equal to the normalization range
839 Bool_t rangeIdentical(kTRUE);
840 RooFIter niter = termNSet.fwdIterator();
841 RooRealVar* normObs;
842 //FK: Here the refRange should be compared to _normRange, if it's set, and to the normObs range if it's not set
843 if(_normRange.Length() > 0) {
844 while ((normObs = (RooRealVar*) niter.next())) {
845 if (normObs->getMin(_normRange.Data()) != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
846 if (normObs->getMax(_normRange.Data()) != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
847 }
848 } else {
849 while ((normObs = (RooRealVar*) niter.next())) {
850 if (normObs->getMin() != normObs->getMin(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
851 if (normObs->getMax() != normObs->getMax(RooNameReg::str(_refRangeName))) rangeIdentical = kFALSE;
852 }
853 }
854// cout<<"FK: rangeIdentical Composite = "<<(rangeIdentical ? 'T':'F') <<endl;
855 if (!rangeIdentical || 1) {
856// cout << "PREPARING RATIO HERE (COMPOSITE TERM)" << endl ;
857 RooAbsReal* ratio = makeCondPdfRatioCorr(*(RooAbsReal*)term->first(), termNSet, termImpSet, normRange(), RooNameReg::str(_refRangeName));
858 ostringstream str; termImpSet.printValue(str);
859 ratioTerms[str.str()].add(*ratio);
860 }
861 }
862 }
863 }
864
865 }
866
867 // Find groups with y as termNSet
868 // Replace G(y) with (G(y),ratio)
869 gIter = groupedList.fwdIterator();
870 while ((group = (RooLinkedList*) gIter.next())) {
871 if (1 == group->GetSize()) {
872 RooArgSet* term = (RooArgSet*) group->At(0);
873
874 Int_t termIdx = terms.IndexOf(term);
875 norm = (RooArgSet*) norms.At(termIdx);
876 imps = (RooArgSet*) imp.At(termIdx);
877 RooArgSet termNSet(*norm), termImpSet(*imps);
878
879 // If termNset matches index of ratioTerms, insert ratio here
880 ostringstream str; termNSet.printValue(str);
881 if (ratioTerms[str.str()].getSize() > 0) {
882// cout << "MUST INSERT RATIO OBJECT IN TERM (SINGLE) " << *term << endl;
883 term->add(ratioTerms[str.str()]);
884 }
885 } else {
886 RooArgSet compTermSet, compTermNorm;
887 RooFIter tIter = group->fwdIterator();
888 RooArgSet* term;
889 while ((term = (RooArgSet*) tIter.next())) {
890 Int_t termIdx = terms.IndexOf(term);
891 norm = (RooArgSet*) norms.At(termIdx);
892 imps = (RooArgSet*) imp.At(termIdx);
893 RooArgSet termNSet(*norm), termImpSet(*imps);
894
895 // If termNset matches index of ratioTerms, insert ratio here
896 ostringstream str; termNSet.printValue(str);
897 if (ratioTerms[str.str()].getSize() > 0) {
898// cout << "MUST INSERT RATIO OBJECT IN TERM (COMPOSITE)" << *term << endl;
899 term->add(ratioTerms[str.str()]);
900 }
901 }
902 }
903 }
904
905 gIter = groupedList.fwdIterator();
906 while ((group = (RooLinkedList*) gIter.next())) {
907// cout << GetName() << ":now processing group" << endl;
908// group->Print("1");
909
910 if (1 == group->GetSize()) {
911// cout << "processing atomic item" << endl;
912 RooArgSet* term = (RooArgSet*) group->At(0);
913
914 Int_t termIdx = terms.IndexOf(term);
915 norm = (RooArgSet*) norms.At(termIdx);
916 integ = (RooArgSet*) ints.At(termIdx);
917 xdeps = (RooArgSet*) cross.At(termIdx);
918 imps = (RooArgSet*) imp.At(termIdx);
919
920 RooArgSet termNSet, termISet, termXSet, termImpSet;
921
922 // Take list of normalization, integrated dependents from factorization algorithm
923 termISet.add(*integ);
924 termNSet.add(*norm);
925
926 // Cross-imported integrated dependents
927 termXSet.add(*xdeps);
928 termImpSet.add(*imps);
929
930// cout << GetName() << ": termISet = " << termISet << endl;
931// cout << GetName() << ": termNSet = " << termNSet << endl;
932// cout << GetName() << ": termXSet = " << termXSet << endl;
933// cout << GetName() << ": termImpSet = " << termImpSet << endl;
934
935 // Add prefab term to partIntList.
936 Bool_t isOwned(kFALSE);
937 vector<RooAbsReal*> func = processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned);
938 if (func[0]) {
939 cache->_partList.add(*func[0]);
940 if (isOwned) cache->_ownedList.addOwned(*func[0]);
941
942 cache->_normList.Add(norm->snapshot(kFALSE));
943
944 cache->_numList.addOwned(*func[1]);
945 cache->_denList.addOwned(*func[2]);
946// cout << "func[0]=" << func[0]->IsA()->GetName() << "::" << func[0]->GetName() << endl;
947// cout << "func[1]=" << func[1]->IsA()->GetName() << "::" << func[1]->GetName() << endl;
948// cout << "func[2]=" << func[2]->IsA()->GetName() << "::" << func[2]->GetName() << endl;
949 }
950 } else {
951// cout << "processing composite item" << endl;
952 RooArgSet compTermSet, compTermNorm, compTermNum, compTermDen;
953 RooFIter tIter = group->fwdIterator();
954 RooArgSet* term;
955 while ((term = (RooArgSet*) tIter.next())) {
956// cout << GetName() << ": processing term " << (*term) << " of composite item" << endl ;
957 Int_t termIdx = terms.IndexOf(term);
958 norm = (RooArgSet*) norms.At(termIdx);
959 integ = (RooArgSet*) ints.At(termIdx);
960 xdeps = (RooArgSet*) cross.At(termIdx);
961 imps = (RooArgSet*) imp.At(termIdx);
962
963 RooArgSet termNSet, termISet, termXSet, termImpSet;
964 termISet.add(*integ);
965 termNSet.add(*norm);
966 termXSet.add(*xdeps);
967 termImpSet.add(*imps);
968
969 // Remove outer integration dependents from termISet
970 termISet.remove(outerIntDeps, kTRUE, kTRUE);
971// cout << "termISet = "; termISet.Print("1");
972
973// cout << GetName() << ": termISet = " << termISet << endl;
974// cout << GetName() << ": termNSet = " << termNSet << endl;
975// cout << GetName() << ": termXSet = " << termXSet << endl;
976// cout << GetName() << ": termImpSet = " << termImpSet << endl;
977
978 Bool_t isOwned;
979 vector<RooAbsReal*> func = processProductTerm(nset, iset, isetRangeName, term, termNSet, termISet, isOwned, kTRUE);
980// cout << GetName() << ": created composite term component " << func[0]->GetName() << endl;
981 if (func[0]) {
982 compTermSet.add(*func[0]);
983 if (isOwned) cache->_ownedList.addOwned(*func[0]);
984 compTermNorm.add(*norm, kFALSE);
985
986 compTermNum.add(*func[1]);
987 compTermDen.add(*func[2]);
988 //cache->_numList.add(*func[1]);
989 //cache->_denList.add(*func[2]);
990
991// cout << "func[0]=" << func[0]->IsA()->GetName() << "::" << func[0]->GetName() << endl;
992// cout << "func[1]=" << func[1]->IsA()->GetName() << "::" << func[1]->GetName() << endl;
993// cout << "func[2]=" << func[2]->IsA()->GetName() << "::" << func[2]->GetName() << endl;
994 }
995 }
996
997// cout << GetName() << ": constructing special composite product" << endl;
998// cout << GetName() << ": compTermSet = " ; compTermSet.Print("1");
999
1000 // WVE THIS NEEDS TO BE REARRANGED
1001
1002 // compTermset is set van partial integrals to be multiplied
1003 // prodtmp = product (compTermSet)
1004 // inttmp = int ( prodtmp ) d (outerIntDeps) _range_isetRangeName
1005
1006 const std::string prodname = makeRGPPName("SPECPROD", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
1007 RooProduct* prodtmp = new RooProduct(prodname.c_str(), prodname.c_str(), compTermSet);
1008 cache->_ownedList.addOwned(*prodtmp);
1009
1010 const std::string intname = makeRGPPName("SPECINT", compTermSet, outerIntDeps, RooArgSet(), isetRangeName);
1011 RooRealIntegral* inttmp = new RooRealIntegral(intname.c_str(), intname.c_str(), *prodtmp, outerIntDeps, 0, 0, isetRangeName);
1012 inttmp->setStringAttribute("PROD_TERM_TYPE", "SPECINT");
1013
1014 cache->_ownedList.addOwned(*inttmp);
1015 cache->_partList.add(*inttmp);
1016
1017 // Product of numerator terms
1018 const string prodname_num = makeRGPPName("SPECPROD_NUM", compTermNum, RooArgSet(), RooArgSet(), 0);
1019 RooProduct* prodtmp_num = new RooProduct(prodname_num.c_str(), prodname_num.c_str(), compTermNum);
1020 prodtmp_num->addOwnedComponents(compTermNum);
1021 cache->_ownedList.addOwned(*prodtmp_num);
1022
1023 // Product of denominator terms
1024 const string prodname_den = makeRGPPName("SPECPROD_DEN", compTermDen, RooArgSet(), RooArgSet(), 0);
1025 RooProduct* prodtmp_den = new RooProduct(prodname_den.c_str(), prodname_den.c_str(), compTermDen);
1026 prodtmp_den->addOwnedComponents(compTermDen);
1027 cache->_ownedList.addOwned(*prodtmp_den);
1028
1029 // Ratio
1030 string name = Form("SPEC_RATIO(%s,%s)", prodname_num.c_str(), prodname_den.c_str());
1031 RooFormulaVar* ndr = new RooFormulaVar(name.c_str(), "@0/@1", RooArgList(*prodtmp_num, *prodtmp_den));
1032
1033 // Integral of ratio
1034 RooAbsReal* numtmp = ndr->createIntegral(outerIntDeps,isetRangeName);
1035 numtmp->addOwnedComponents(*ndr);
1036
1037 cache->_numList.addOwned(*numtmp);
1039 cache->_normList.Add(compTermNorm.snapshot(kFALSE));
1040 }
1041 }
1042
1043 // Store the partial integral list and return the assigned code
1044 code = _cacheMgr.setObj(nset, iset, (RooAbsCacheElement*)cache, RooNameReg::ptr(isetRangeName));
1045
1046 // Fill references to be returned
1047 partList = &cache->_partList;
1048 nsetList = &cache->_normList;
1049
1050 // WVE DEBUG PRINTING
1051// cout << "RooProdPdf::getPartIntList(" << GetName() << ") made cache " << cache << " with the following nset pointers ";
1052// TIterator* nliter = nsetList->MakeIterator();
1053// RooArgSet* ns;
1054// while((ns=(RooArgSet*)nliter->Next())) {
1055// cout << ns << " ";
1056// }
1057// cout << endl;
1058// delete nliter;
1059
1060// cout << " FOLKERT::RooProdPdf::getPartIntList END(" << GetName() <<") nset = " << (nset?*nset:RooArgSet()) << endl
1061// << " _normRange = " << _normRange << endl
1062// << " iset = " << (iset?*iset:RooArgSet()) << endl
1063// << " partList = ";
1064// if(partListPointer) partListPointer->Print();
1065// cout << " nsetList = ";
1066// if(nsetListPointer) nsetListPointer->Print("");
1067// cout << " code = " << code << endl
1068// << " isetRangeName = " << (isetRangeName?isetRangeName:"<null>") << endl;
1069
1070
1071 // Need to rearrange product in case of multiple ranges
1072 if (_normRange.Contains(",")) {
1073 rearrangeProduct(*cache);
1074 }
1075
1076 // We own contents of all lists filled by factorizeProduct()
1077 groupedList.Delete();
1078 terms.Delete();
1079 ints.Delete();
1080 imp.Delete();
1081 norms.Delete();
1082 cross.Delete();
1083}
1084
1085
1086
1087
1088////////////////////////////////////////////////////////////////////////////////
1089/// For single normalization ranges
1090
1091RooAbsReal* RooProdPdf::makeCondPdfRatioCorr(RooAbsReal& pdf, const RooArgSet& termNset, const RooArgSet& /*termImpSet*/, const char* normRangeTmp, const char* refRange) const
1092{
1093 RooAbsReal* ratio_num = pdf.createIntegral(termNset,normRangeTmp) ;
1094 RooAbsReal* ratio_den = pdf.createIntegral(termNset,refRange) ;
1095 RooFormulaVar* ratio = new RooFormulaVar(Form("ratio(%s,%s)",ratio_num->GetName(),ratio_den->GetName()),"@0/@1",
1096 RooArgList(*ratio_num,*ratio_den)) ;
1097
1098 ratio->addOwnedComponents(RooArgSet(*ratio_num,*ratio_den)) ;
1099 ratio->setAttribute("RATIO_TERM") ;
1100 return ratio ;
1101}
1102
1103
1104
1105
1106////////////////////////////////////////////////////////////////////////////////
1107
1109{
1110 RooAbsReal* part, *num, *den ;
1111 RooArgSet nomList ;
1112
1113 list<string> rangeComps ;
1114 {
1115 char* buf = new char[strlen(_normRange.Data()) + 1] ;
1116 strcpy(buf,_normRange.Data()) ;
1117 char* save(0) ;
1118 char* token = R__STRTOK_R(buf,",",&save) ;
1119 while(token) {
1120 rangeComps.push_back(token) ;
1121 token = R__STRTOK_R(0,",",&save) ;
1122 }
1123 delete[] buf;
1124 }
1125
1126
1127 map<string,RooArgSet> denListList ;
1128 RooArgSet specIntDeps ;
1129 string specIntRange ;
1130
1131// cout << "THIS IS REARRANGEPRODUCT" << endl ;
1132
1133 RooFIter iterp = cache._partList.fwdIterator() ;
1134 RooFIter iter1 = cache._numList.fwdIterator() ;
1135 RooFIter iter2 = cache._denList.fwdIterator() ;
1136 RooFIter itern = cache._normList.fwdIterator() ;
1137 while((part=(RooAbsReal*)iterp.next())) {
1138
1139 itern.next() ;
1140 num = (RooAbsReal*) iter1.next() ;
1141 den = (RooAbsReal*) iter2.next() ;
1142
1143// cout << "now processing part " << part->GetName() << " of type " << part->getStringAttribute("PROD_TERM_TYPE") << endl ;
1144// cout << "corresponding numerator = " << num->GetName() << endl ;
1145// cout << "corresponding denominator = " << den->GetName() << endl ;
1146
1147
1148 RooFormulaVar* ratio(0) ;
1149 RooArgSet origNumTerm ;
1150
1151 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1152
1153 RooRealIntegral* orig = (RooRealIntegral*) num;
1154 RooFormulaVar* specratio = (RooFormulaVar*) &orig->integrand() ;
1155 RooProduct* func = (RooProduct*) specratio->getParameter(0) ;
1156
1157 RooArgSet* comps = orig->getComponents() ;
1158 RooFIter iter = comps->fwdIterator() ;
1159 RooAbsArg* carg ;
1160 while((carg=(RooAbsArg*)iter.next())) {
1161 if (carg->getAttribute("RATIO_TERM")) {
1162 ratio = (RooFormulaVar*)carg ;
1163 break ;
1164 }
1165 }
1166 delete comps ;
1167
1168 if (ratio) {
1169 RooCustomizer cust(*func,"blah") ;
1170 cust.replaceArg(*ratio,RooFit::RooConst(1)) ;
1171 RooAbsArg* funcCust = cust.build() ;
1172// cout << "customized function = " << endl ;
1173// funcCust->printComponentTree() ;
1174 nomList.add(*funcCust) ;
1175 } else {
1176 nomList.add(*func) ;
1177 }
1178
1179
1180 } else {
1181
1182 // Find the ratio term
1183 RooAbsReal* func = num;
1184 // If top level object is integral, navigate to integrand
1185 if (func->InheritsFrom(RooRealIntegral::Class())) {
1186 func = (RooAbsReal*) &((RooRealIntegral*)(func))->integrand();
1187 }
1188 if (func->InheritsFrom(RooProduct::Class())) {
1189// cout << "product term found: " ; func->Print() ;
1190 RooArgSet comps(((RooProduct*)(func))->components()) ;
1191 RooFIter iter = comps.fwdIterator() ;
1192 RooAbsArg* arg ;
1193 while((arg=(RooAbsArg*)iter.next())) {
1194 if (arg->getAttribute("RATIO_TERM")) {
1195 ratio = (RooFormulaVar*)(arg) ;
1196 } else {
1197 origNumTerm.add(*arg) ;
1198 }
1199 }
1200 }
1201
1202 if (ratio) {
1203// cout << "Found ratio term in numerator: " << ratio->GetName() << endl ;
1204// cout << "Adding only original term to numerator: " << origNumTerm << endl ;
1205 nomList.add(origNumTerm) ;
1206 } else {
1207 nomList.add(*num) ;
1208 }
1209
1210 }
1211
1212 for (list<string>::iterator iter = rangeComps.begin() ; iter != rangeComps.end() ; ++iter) {
1213 // If denominator is an integral, make a clone with the integration range adjusted to
1214 // the selected component of the normalization integral
1215// cout << "NOW PROCESSING DENOMINATOR " << den->IsA()->GetName() << "::" << den->GetName() << endl ;
1216
1217 if (string("SPECINT")==part->getStringAttribute("PROD_TERM_TYPE")) {
1218
1219// cout << "create integral: SPECINT case" << endl ;
1220 RooRealIntegral* orig = (RooRealIntegral*) num;
1221 RooFormulaVar* specRatio = (RooFormulaVar*) &orig->integrand() ;
1222 specIntDeps.add(orig->intVars()) ;
1223 if (orig->intRange()) {
1224 specIntRange = orig->intRange() ;
1225 }
1226 //RooProduct* numtmp = (RooProduct*) specRatio->getParameter(0) ;
1227 RooProduct* dentmp = (RooProduct*) specRatio->getParameter(1) ;
1228
1229// cout << "numtmp = " << numtmp->IsA()->GetName() << "::" << numtmp->GetName() << endl ;
1230// cout << "dentmp = " << dentmp->IsA()->GetName() << "::" << dentmp->GetName() << endl ;
1231
1232// cout << "denominator components are " << dentmp->components() << endl ;
1233 RooArgSet comps(dentmp->components()) ;
1234 RooFIter piter = comps.fwdIterator() ;
1235 RooAbsReal* parg ;
1236 while((parg=(RooAbsReal*)piter.next())) {
1237// cout << "now processing denominator component " << parg->IsA()->GetName() << "::" << parg->GetName() << endl ;
1238
1239 if (ratio && parg->dependsOn(*ratio)) {
1240// cout << "depends in value of ratio" << endl ;
1241
1242 // Make specialize ratio instance
1243 RooAbsReal* specializedRatio = specializeRatio(*(RooFormulaVar*)ratio,iter->c_str()) ;
1244
1245// cout << "specRatio = " << endl ;
1246// specializedRatio->printComponentTree() ;
1247
1248 // Replace generic ratio with specialized ratio
1249 RooAbsArg *partCust(0) ;
1250 if (parg->InheritsFrom(RooAddition::Class())) {
1251
1252
1253
1254 RooAddition* tmpadd = (RooAddition*)(parg) ;
1255
1256 RooCustomizer cust(*tmpadd->list1().first(),Form("blah_%s",iter->c_str())) ;
1257 cust.replaceArg(*ratio,*specializedRatio) ;
1258 partCust = cust.build() ;
1259
1260 } else {
1261 RooCustomizer cust(*parg,Form("blah_%s",iter->c_str())) ;
1262 cust.replaceArg(*ratio,*specializedRatio) ;
1263 partCust = cust.build() ;
1264 }
1265
1266 // Print customized denominator
1267// cout << "customized function = " << endl ;
1268// partCust->printComponentTree() ;
1269
1270 RooAbsReal* specializedPartCust = specializeIntegral(*(RooAbsReal*)partCust,iter->c_str()) ;
1271
1272 // Finally divide again by ratio
1273 string name = Form("%s_divided_by_ratio",specializedPartCust->GetName()) ;
1274 RooFormulaVar* specIntFinal = new RooFormulaVar(name.c_str(),"@0/@1",RooArgList(*specializedPartCust,*specializedRatio)) ;
1275
1276 denListList[*iter].add(*specIntFinal) ;
1277 } else {
1278
1279// cout << "does NOT depend on value of ratio" << endl ;
1280// parg->Print("t") ;
1281
1282 denListList[*iter].add(*specializeIntegral(*parg,iter->c_str())) ;
1283
1284 }
1285 }
1286// cout << "end iteration over denominator components" << endl ;
1287 } else {
1288
1289 if (ratio) {
1290
1291 RooAbsReal* specRatio = specializeRatio(*(RooFormulaVar*)ratio,iter->c_str()) ;
1292
1293 // If integral is 'Int r(y)*g(y) dy ' then divide a posteriori by r(y)
1294// cout << "have ratio, orig den = " << den->GetName() << endl ;
1295
1296 RooArgSet tmp(origNumTerm) ;
1297 tmp.add(*specRatio) ;
1298 const string pname = makeRGPPName("PROD",tmp,RooArgSet(),RooArgSet(),0) ;
1299 RooProduct* specDenProd = new RooProduct(pname.c_str(),pname.c_str(),tmp) ;
1300 RooAbsReal* specInt(0) ;
1301
1303 specInt = specDenProd->createIntegral(((RooRealIntegral*)den)->intVars(),iter->c_str()) ;
1304 } else if (den->InheritsFrom(RooAddition::Class())) {
1305 RooAddition* orig = (RooAddition*)den ;
1306 RooRealIntegral* origInt = (RooRealIntegral*) orig->list1().first() ;
1307 specInt = specDenProd->createIntegral(origInt->intVars(),iter->c_str()) ;
1308 } else {
1309 throw string("this should not happen") ;
1310 }
1311
1312 //RooAbsReal* specInt = specializeIntegral(*den,iter->c_str()) ;
1313 string name = Form("%s_divided_by_ratio",specInt->GetName()) ;
1314 RooFormulaVar* specIntFinal = new RooFormulaVar(name.c_str(),"@0/@1",RooArgList(*specInt,*specRatio)) ;
1315 denListList[*iter].add(*specIntFinal) ;
1316 } else {
1317 denListList[*iter].add(*specializeIntegral(*den,iter->c_str())) ;
1318 }
1319
1320 }
1321 }
1322
1323 }
1324
1325 // Do not rearrage terms if numerator and denominator are effectively empty
1326 if (nomList.getSize()==0) {
1327 return ;
1328 }
1329
1330 string name = Form("%s_numerator",GetName()) ;
1331 // WVE FIX THIS (2)
1332
1333 RooAbsReal* numerator = new RooProduct(name.c_str(),name.c_str(),nomList) ;
1334
1335 RooArgSet products ;
1336// cout << "nomList = " << nomList << endl ;
1337 for (map<string,RooArgSet>::iterator iter = denListList.begin() ; iter != denListList.end() ; ++iter) {
1338// cout << "denList[" << iter->first << "] = " << iter->second << endl ;
1339 name = Form("%s_denominator_comp_%s",GetName(),iter->first.c_str()) ;
1340 // WVE FIX THIS (2)
1341 RooProduct* prod_comp = new RooProduct(name.c_str(),name.c_str(),iter->second) ;
1342 products.add(*prod_comp) ;
1343 }
1344 name = Form("%s_denominator_sum",GetName()) ;
1345 RooAbsReal* norm = new RooAddition(name.c_str(),name.c_str(),products) ;
1346 norm->addOwnedComponents(products) ;
1347
1348 if (specIntDeps.getSize()>0) {
1349 // Apply posterior integration required for SPECINT case
1350
1351 string namesr = Form("SPEC_RATIO(%s,%s)",numerator->GetName(),norm->GetName()) ;
1352 RooFormulaVar* ndr = new RooFormulaVar(namesr.c_str(),"@0/@1",RooArgList(*numerator,*norm)) ;
1353
1354 // Integral of ratio
1355 RooAbsReal* numtmp = ndr->createIntegral(specIntDeps,specIntRange.c_str()) ;
1356
1357 numerator = numtmp ;
1358 norm = (RooAbsReal*) RooFit::RooConst(1).Clone() ;
1359 }
1360
1361
1362// cout << "numerator" << endl ;
1363// numerator->printComponentTree("",0,5) ;
1364// cout << "denominator" << endl ;
1365// norm->printComponentTree("",0,5) ;
1366
1367
1368 // WVE DEBUG
1369 //RooMsgService::instance().debugWorkspace()->import(RooArgSet(*numerator,*norm)) ;
1370
1371 cache._rearrangedNum.reset(numerator);
1372 cache._rearrangedDen.reset(norm);
1373 cache._isRearranged = kTRUE ;
1374
1375}
1376
1377
1378////////////////////////////////////////////////////////////////////////////////
1379
1380RooAbsReal* RooProdPdf::specializeRatio(RooFormulaVar& input, const char* targetRangeName) const
1381{
1382 RooRealIntegral* numint = (RooRealIntegral*) input.getParameter(0) ;
1383 RooRealIntegral* denint = (RooRealIntegral*) input.getParameter(1) ;
1384
1385 RooAbsReal* numint_spec = specializeIntegral(*numint,targetRangeName) ;
1386
1387 RooAbsReal* ret = new RooFormulaVar(Form("ratio(%s,%s)",numint_spec->GetName(),denint->GetName()),"@0/@1",RooArgList(*numint_spec,*denint)) ;
1388 ret->addOwnedComponents(*numint_spec) ;
1389
1390 return ret ;
1391}
1392
1393
1394
1395////////////////////////////////////////////////////////////////////////////////
1396
1397RooAbsReal* RooProdPdf::specializeIntegral(RooAbsReal& input, const char* targetRangeName) const
1398{
1399 if (input.InheritsFrom(RooRealIntegral::Class())) {
1400
1401 // If input is integral, recreate integral but override integration range to be targetRangeName
1402 RooRealIntegral* orig = (RooRealIntegral*)&input ;
1403// cout << "creating integral: integrand = " << orig->integrand().GetName() << " vars = " << orig->intVars() << " range = " << targetRangeName << endl ;
1404 return orig->integrand().createIntegral(orig->intVars(),targetRangeName) ;
1405
1406 } else if (input.InheritsFrom(RooAddition::Class())) {
1407
1408 // If input is sum of integrals, recreate integral from first component of set, but override integration range to be targetRangeName
1409 RooAddition* orig = (RooAddition*)&input ;
1410 RooRealIntegral* origInt = (RooRealIntegral*) orig->list1().first() ;
1411// cout << "creating integral from addition: integrand = " << origInt->integrand().GetName() << " vars = " << origInt->intVars() << " range = " << targetRangeName << endl ;
1412 return origInt->integrand().createIntegral(origInt->intVars(),targetRangeName) ;
1413
1414 } else {
1415
1416// cout << "specializeIntegral: unknown input type " << input.IsA()->GetName() << "::" << input.GetName() << endl ;
1417 }
1418
1419 return &input ;
1420}
1421
1422
1423////////////////////////////////////////////////////////////////////////////////
1424/// Group product into terms that can be calculated independently
1425
1427 const RooLinkedList& terms, const RooLinkedList& norms,
1428 const RooLinkedList& imps, const RooLinkedList& ints, const RooLinkedList& /*cross*/) const
1429{
1430 // Start out with each term in its own group
1431 RooFIter tIter = terms.fwdIterator() ;
1432 RooArgSet* term ;
1433 while((term=(RooArgSet*)tIter.next())) {
1435 group->Add(term) ;
1436 groupedTerms.Add(group) ;
1437 }
1438
1439 // Make list of imported dependents that occur in any term
1440 RooArgSet allImpDeps ;
1441 RooFIter iIter = imps.fwdIterator() ;
1442 RooArgSet *impDeps ;
1443 while((impDeps=(RooArgSet*)iIter.next())) {
1444 allImpDeps.add(*impDeps,kFALSE) ;
1445 }
1446
1447 // Make list of integrated dependents that occur in any term
1448 RooArgSet allIntDeps ;
1449 iIter = ints.fwdIterator() ;
1450 RooArgSet *intDeps ;
1451 while((intDeps=(RooArgSet*)iIter.next())) {
1452 allIntDeps.add(*intDeps,kFALSE) ;
1453 }
1454
1455 RooArgSet* tmp = (RooArgSet*) allIntDeps.selectCommon(allImpDeps) ;
1456 outerIntDeps.removeAll() ;
1457 outerIntDeps.add(*tmp) ;
1458 delete tmp ;
1459
1460 // Now iteratively merge groups that should be (partially) integrated together
1461 RooFIter oidIter = outerIntDeps.fwdIterator() ;
1462 RooAbsArg* outerIntDep ;
1463 while ((outerIntDep =(RooAbsArg*)oidIter.next())) {
1464
1465 // Collect groups that feature this dependent
1466 RooLinkedList* newGroup = 0 ;
1467
1468 // Loop over groups
1470 RooFIter glIter = groupedTerms.fwdIterator() ;
1471 Bool_t needMerge = kFALSE ;
1472 while((group=(RooLinkedList*)glIter.next())) {
1473
1474 // See if any term in this group depends in any ay on outerDepInt
1475 RooArgSet* term2 ;
1476 RooFIter tIter2 = group->fwdIterator() ;
1477 while((term2=(RooArgSet*)tIter2.next())) {
1478
1479 Int_t termIdx = terms.IndexOf(term2) ;
1480 RooArgSet* termNormDeps = (RooArgSet*) norms.At(termIdx) ;
1481 RooArgSet* termIntDeps = (RooArgSet*) ints.At(termIdx) ;
1482 RooArgSet* termImpDeps = (RooArgSet*) imps.At(termIdx) ;
1483
1484 if (termNormDeps->contains(*outerIntDep) ||
1485 termIntDeps->contains(*outerIntDep) ||
1486 termImpDeps->contains(*outerIntDep)) {
1487 needMerge = kTRUE ;
1488 }
1489
1490 }
1491
1492 if (needMerge) {
1493 // Create composite group if not yet existing
1494 if (newGroup==0) {
1495 newGroup = new RooLinkedList ;
1496 }
1497
1498 // Add terms of this group to new term
1499 tIter2 = group->fwdIterator() ;
1500 while((term2=(RooArgSet*)tIter2.next())) {
1501 newGroup->Add(term2) ;
1502 }
1503
1504 // Remove this group from list and delete it (but not its contents)
1505 groupedTerms.Remove(group) ;
1506 delete group ;
1507 }
1508 }
1509 // If a new group has been created to merge terms dependent on current outerIntDep, add it to group list
1510 if (newGroup) {
1511 groupedTerms.Add(newGroup) ;
1512 }
1513
1514 }
1515}
1516
1517
1518
1519////////////////////////////////////////////////////////////////////////////////
1520/// Calculate integrals of factorized product terms over observables iset while normalized
1521/// to observables in nset.
1522
1523std::vector<RooAbsReal*> RooProdPdf::processProductTerm(const RooArgSet* nset, const RooArgSet* iset, const char* isetRangeName,
1524 const RooArgSet* term,const RooArgSet& termNSet, const RooArgSet& termISet,
1525 Bool_t& isOwned, Bool_t forceWrap) const
1526{
1527// cout << " FOLKERT::RooProdPdf(" << GetName() <<") processProductTerm nset = " << (nset?*nset:RooArgSet()) << endl
1528// << " _normRange = " << _normRange << endl
1529// << " iset = " << (iset?*iset:RooArgSet()) << endl
1530// << " isetRangeName = " << (isetRangeName?isetRangeName:"<null>") << endl
1531// << " term = " << (term?*term:RooArgSet()) << endl
1532// << " termNSet = " << termNSet << endl
1533// << " termISet = " << termISet << endl
1534// << " isOwned = " << isOwned << endl
1535// << " forceWrap = " << forceWrap << endl ;
1536
1537 vector<RooAbsReal*> ret(3) ; ret[0] = 0 ; ret[1] = 0 ; ret[2] = 0 ;
1538
1539 // CASE I: factorizing term: term is integrated over all normalizing observables
1540 // -----------------------------------------------------------------------------
1541 // Check if all observbales of this term are integrated. If so the term cancels
1542 if (termNSet.getSize()>0 && termNSet.getSize()==termISet.getSize() && isetRangeName==0) {
1543
1544
1545 //cout << "processProductTerm(" << GetName() << ") case I " << endl ;
1546
1547 // Term factorizes
1548 return ret ;
1549 }
1550
1551 // CASE II: Dropped terms: if term is entirely unnormalized, it should be dropped
1552 // ------------------------------------------------------------------------------
1553 if (nset && termNSet.getSize()==0) {
1554
1555 //cout << "processProductTerm(" << GetName() << ") case II " << endl ;
1556
1557 // Drop terms that are not asked to be normalized
1558 return ret ;
1559 }
1560
1561 if (iset && termISet.getSize()>0) {
1562 if (term->getSize()==1) {
1563
1564 // CASE IIIa: Normalized and partially integrated single PDF term
1565 //---------------------------------------------------------------
1566
1567 RooAbsPdf* pdf = (RooAbsPdf*) term->first() ;
1568
1569 RooAbsReal* partInt = pdf->createIntegral(termISet,termNSet,isetRangeName) ;
1570 //partInt->setOperMode(operMode()) ;
1571 partInt->setStringAttribute("PROD_TERM_TYPE","IIIa") ;
1572
1573 isOwned=kTRUE ;
1574
1575 //cout << "processProductTerm(" << GetName() << ") case IIIa func = " << partInt->GetName() << endl ;
1576
1577 ret[0] = partInt ;
1578
1579 // Split mode results
1580 ret[1] = pdf->createIntegral(termISet,isetRangeName) ;
1581 ret[2] = pdf->createIntegral(termNSet,normRange()) ;
1582
1583 return ret ;
1584
1585
1586 } else {
1587
1588 // CASE IIIb: Normalized and partially integrated composite PDF term
1589 //---------------------------------------------------------------
1590
1591 // Use auxiliary class RooGenProdProj to calculate this term
1592 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1593 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName) ;
1594 partInt->setStringAttribute("PROD_TERM_TYPE","IIIb") ;
1595 //partInt->setOperMode(operMode()) ;
1596
1597 //cout << "processProductTerm(" << GetName() << ") case IIIb func = " << partInt->GetName() << endl ;
1598
1599 isOwned=kTRUE ;
1600 ret[0] = partInt ;
1601
1602 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),0) ;
1603
1604 // WVE FIX THIS
1605 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1606
1607 ret[1] = tmp_prod->createIntegral(termISet,isetRangeName) ;
1608 ret[2] = tmp_prod->createIntegral(termNSet,normRange()) ;
1609
1610 return ret ;
1611 }
1612 }
1613
1614 // CASE IVa: Normalized non-integrated composite PDF term
1615 // -------------------------------------------------------
1616 if (nset && nset->getSize()>0 && term->getSize()>1) {
1617 // Composite term needs normalized integration
1618
1619 const std::string name = makeRGPPName("GENPROJ_",*term,termISet,termNSet,isetRangeName) ;
1620 RooAbsReal* partInt = new RooGenProdProj(name.c_str(),name.c_str(),*term,termISet,termNSet,isetRangeName,normRange()) ;
1622
1623 partInt->setStringAttribute("PROD_TERM_TYPE","IVa") ;
1624 //partInt->setOperMode(operMode()) ;
1625
1626 //cout << "processProductTerm(" << GetName() << ") case IVa func = " << partInt->GetName() << endl ;
1627
1628 isOwned=kTRUE ;
1629 ret[0] = partInt ;
1630
1631 const std::string name1 = makeRGPPName("PROD",*term,RooArgSet(),RooArgSet(),0) ;
1632
1633 // WVE FIX THIS
1634 RooProduct* tmp_prod = new RooProduct(name1.c_str(),name1.c_str(),*term) ;
1635
1636 ret[1] = tmp_prod->createIntegral(termISet,isetRangeName) ;
1637 ret[2] = tmp_prod->createIntegral(termNSet,normRange()) ;
1638
1639 return ret ;
1640 }
1641
1642 // CASE IVb: Normalized, non-integrated single PDF term
1643 // -----------------------------------------------------
1644 RooFIter pIter = term->fwdIterator() ;
1645 RooAbsPdf* pdf ;
1646 while((pdf=(RooAbsPdf*)pIter.next())) {
1647
1648 if (forceWrap) {
1649
1650 // Construct representative name of normalization wrapper
1651 TString name(pdf->GetName()) ;
1652 name.Append("_NORM[") ;
1653 RooFIter nIter = termNSet.fwdIterator() ;
1654 RooAbsArg* arg ;
1655 Bool_t first(kTRUE) ;
1656 while((arg=(RooAbsArg*)nIter.next())) {
1657 if (!first) {
1658 name.Append(",") ;
1659 } else {
1660 first=kFALSE ;
1661 }
1662 name.Append(arg->GetName()) ;
1663 }
1664 if (normRange()) {
1665 name.Append("|") ;
1666 name.Append(normRange()) ;
1667 }
1668 name.Append("]") ;
1669
1670 RooAbsReal* partInt = new RooRealIntegral(name.Data(),name.Data(),*pdf,RooArgSet(),&termNSet) ;
1671 partInt->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1672 isOwned=kTRUE ;
1673
1674 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << partInt->GetName() << endl ;
1675
1676 ret[0] = partInt ;
1677
1678 ret[1] = pdf->createIntegral(RooArgSet()) ;
1679 ret[2] = pdf->createIntegral(termNSet,normRange()) ;
1680
1681 return ret ;
1682
1683
1684 } else {
1685 isOwned=kFALSE ;
1686
1687 //cout << "processProductTerm(" << GetName() << ") case IVb func = " << pdf->GetName() << endl ;
1688
1689
1690 pdf->setStringAttribute("PROD_TERM_TYPE","IVb") ;
1691 ret[0] = pdf ;
1692
1693 ret[1] = pdf->createIntegral(RooArgSet()) ;
1694 ret[2] = termNSet.getSize()>0 ? pdf->createIntegral(termNSet,normRange()) : ((RooAbsReal*)RooFit::RooConst(1).clone("1")) ;
1695 return ret ;
1696 }
1697 }
1698
1699 coutE(Eval) << "RooProdPdf::processProductTerm(" << GetName() << ") unidentified term!!!" << endl ;
1700 return ret ;
1701}
1702
1703
1704
1705
1706////////////////////////////////////////////////////////////////////////////////
1707/// Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1708
1709std::string RooProdPdf::makeRGPPName(const char* pfx, const RooArgSet& term, const RooArgSet& iset,
1710 const RooArgSet& nset, const char* isetRangeName) const
1711{
1712 // Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
1713
1714 std::ostringstream os(pfx);
1715 os << "[";
1716
1717 RooFIter pIter = term.fwdIterator() ;
1718 // Encode component names
1719 Bool_t first(kTRUE) ;
1720 RooAbsPdf* pdf ;
1721 while ((pdf=(RooAbsPdf*)pIter.next())) {
1722 if (!first) os << "_X_";
1723 first = kFALSE;
1724 os << pdf->GetName();
1725 }
1726 os << "]" << integralNameSuffix(iset,&nset,isetRangeName,kTRUE);
1727
1728 return os.str();
1729}
1730
1731
1732
1733////////////////////////////////////////////////////////////////////////////////
1734/// Force RooRealIntegral to offer all observables for internal integration
1735
1737{
1738 return kTRUE ;
1739}
1740
1741
1742
1743////////////////////////////////////////////////////////////////////////////////
1744/// Determine which part (if any) of given integral can be performed analytically.
1745/// If any analytical integration is possible, return integration scenario code.
1746///
1747/// RooProdPdf implements two strategies in implementing analytical integrals
1748///
1749/// First, PDF components whose entire set of dependents are requested to be integrated
1750/// can be dropped from the product, as they will integrate out to 1 by construction
1751///
1752/// Second, RooProdPdf queries each remaining component PDF for its analytical integration
1753/// capability of the requested set ('allVars'). It finds the largest common set of variables
1754/// that can be integrated by all remaining components. If such a set exists, it reconfirms that
1755/// each component is capable of analytically integrating the common set, and combines the components
1756/// individual integration codes into a single integration code valid for RooProdPdf.
1757
1759 const RooArgSet* normSet, const char* rangeName) const
1760{
1761 if (_forceNumInt) return 0 ;
1762
1763 // Declare that we can analytically integrate all requested observables
1764 analVars.add(allVars) ;
1765
1766 // Retrieve (or create) the required partial integral list
1767 Int_t code ;
1768 RooArgList *plist(0) ;
1769 RooLinkedList *nlist(0) ;
1770 getPartIntList(normSet,&allVars,plist,nlist,code,rangeName) ;
1771
1772 return code+1 ;
1773}
1774
1775
1776
1777
1778////////////////////////////////////////////////////////////////////////////////
1779/// Return analytical integral defined by given scenario code
1780
1781Double_t RooProdPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
1782{
1783 // No integration scenario
1784 if (code==0) {
1785 return getVal(normSet) ;
1786 }
1787
1788
1789 // WVE needs adaptation for rangename feature
1790
1791 // Partial integration scenarios
1792 CacheElem* cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
1793
1794 RooArgList* partIntList(0) ;
1795 RooLinkedList* normList(0) ;
1796
1797 // If cache has been sterilized, revive this slot
1798 if (cache==0) {
1799 RooArgSet* vars = getParameters(RooArgSet()) ;
1800 RooArgSet* nset = _cacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
1801 RooArgSet* iset = _cacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
1802
1803 Int_t code2(-1) ;
1804 getPartIntList(nset,iset,partIntList,normList,code2,rangeName) ;
1805
1806 delete vars ;
1807
1808 // preceding call to getPartIntList guarantees non-null return
1809 // coverity[NULL_RETURNS]
1810 cache = (CacheElem*) _cacheMgr.getObj(nset,iset,&code2,rangeName) ;
1811
1812 delete nset ;
1813 delete iset ;
1814
1815 } else {
1816
1817 partIntList = &cache->_partList ;
1818 normList = &cache->_normList ;
1819
1820 }
1821
1822// Double_t val = calculate(partIntList,normList) ;
1823
1824 Double_t val = calculate(*cache,kTRUE) ;
1825// cout << "RPP::aIWN(" << GetName() << ") ,code = " << code << ", value = " << val << endl ;
1826
1827 return val ;
1828}
1829
1830
1831
1832////////////////////////////////////////////////////////////////////////////////
1833/// Obsolete
1834
1836{ return kFALSE ; }
1837
1838
1839
1840
1841////////////////////////////////////////////////////////////////////////////////
1842/// If this product contains exactly one extendable p.d.f return the extension abilities of
1843/// that p.d.f, otherwise return CanNotBeExtended
1844
1846{
1848}
1849
1850
1851
1852////////////////////////////////////////////////////////////////////////////////
1853/// Return the expected number of events associated with the extentable input p.d.f
1854/// in the product. If there is no extendable term, return zero and issue and error
1855
1857{
1858 if (_extendedIndex<0) {
1859 coutE(Generation) << "ERROR: Requesting expected number of events from a RooProdPdf that does not contain an extended p.d.f" << endl ;
1860 }
1861 assert(_extendedIndex>=0) ;
1862 return ((RooAbsPdf*)_pdfList.at(_extendedIndex))->expectedEvents(nset) ;
1863}
1864
1865
1866
1867////////////////////////////////////////////////////////////////////////////////
1868/// Return generator context optimized for generating events from product p.d.f.s
1869
1871 const RooArgSet* auxProto, Bool_t verbose) const
1872{
1873 if (_useDefaultGen) return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ;
1874 return new RooProdGenContext(*this,vars,prototype,auxProto,verbose) ;
1875}
1876
1877
1878
1879////////////////////////////////////////////////////////////////////////////////
1880/// Query internal generation capabilities of component p.d.f.s and aggregate capabilities
1881/// into master configuration passed to the generator context
1882
1883Int_t RooProdPdf::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
1884{
1885 if (!_useDefaultGen) return 0 ;
1886
1887 // Find the subset directVars that only depend on a single PDF in the product
1888 RooArgSet directSafe ;
1889 RooFIter dIter = directVars.fwdIterator() ;
1890 RooAbsArg* arg ;
1891 while((arg=(RooAbsArg*)dIter.next())) {
1892 if (isDirectGenSafe(*arg)) directSafe.add(*arg) ;
1893 }
1894
1895
1896 // Now find direct integrator for relevant components ;
1897 RooAbsPdf* pdf ;
1898 std::vector<Int_t> code;
1899 code.reserve(64);
1900 RooFIter pdfIter = _pdfList.fwdIterator();
1901 while((pdf=(RooAbsPdf*)pdfIter.next())) {
1902 RooArgSet pdfDirect ;
1903 Int_t pdfCode = pdf->getGenerator(directSafe,pdfDirect,staticInitOK);
1904 code.push_back(pdfCode);
1905 if (pdfCode != 0) {
1906 generateVars.add(pdfDirect) ;
1907 }
1908 }
1909
1910
1911 if (generateVars.getSize()>0) {
1912 Int_t masterCode = _genCode.store(code) ;
1913 return masterCode+1 ;
1914 } else {
1915 return 0 ;
1916 }
1917}
1918
1919
1920
1921////////////////////////////////////////////////////////////////////////////////
1922/// Forward one-time initialization call to component generation initialization
1923/// methods.
1924
1926{
1927 if (!_useDefaultGen) return ;
1928
1929 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1930 RooAbsPdf* pdf ;
1931 Int_t i(0) ;
1932 RooFIter pdfIter = _pdfList.fwdIterator();
1933 while((pdf=(RooAbsPdf*)pdfIter.next())) {
1934 if (codeList[i]!=0) {
1935 pdf->initGenerator(codeList[i]) ;
1936 }
1937 i++ ;
1938 }
1939}
1940
1941
1942
1943////////////////////////////////////////////////////////////////////////////////
1944/// Generate a single event with configuration specified by 'code'
1945/// Defer internal generation to components as encoded in the _genCode
1946/// registry for given generator code.
1947
1949{
1950 if (!_useDefaultGen) return ;
1951
1952 const std::vector<Int_t>& codeList = _genCode.retrieve(code-1) ;
1953 RooAbsPdf* pdf ;
1954 Int_t i(0) ;
1955 RooFIter pdfIter = _pdfList.fwdIterator();
1956 while((pdf=(RooAbsPdf*)pdfIter.next())) {
1957 if (codeList[i]!=0) {
1958 pdf->generateEvent(codeList[i]) ;
1959 }
1960 i++ ;
1961 }
1962}
1963
1964
1965
1966////////////////////////////////////////////////////////////////////////////////
1967/// Destructor
1968
1970{
1971 _normList.Delete() ; //WVE THIS IS AN INTENTIAL LEAK -- MUST FIX LATER
1972// cout << "RooProdPdf::CacheElem dtor, this = " << this << endl ;
1973}
1974
1975
1976
1977////////////////////////////////////////////////////////////////////////////////
1978/// Return RooAbsArg components contained in the cache
1979
1981{
1982 RooArgList ret ;
1983 ret.add(_partList) ;
1984 ret.add(_numList) ;
1985 ret.add(_denList) ;
1986 if (_rearrangedNum) ret.add(*_rearrangedNum) ;
1987 if (_rearrangedDen) ret.add(*_rearrangedDen) ;
1988 return ret ;
1989
1990}
1991
1992
1993
1994////////////////////////////////////////////////////////////////////////////////
1995/// Hook function to print cache contents in tree printing of RooProdPdf
1996
1997void RooProdPdf::CacheElem::printCompactTreeHook(ostream& os, const char* indent, Int_t curElem, Int_t maxElem)
1998{
1999 if (curElem==0) {
2000 os << indent << "RooProdPdf begin partial integral cache" << endl ;
2001 }
2002
2003 RooFIter iter = _partList.fwdIterator() ;
2004 RooAbsArg* arg ;
2005 TString indent2(indent) ;
2006 indent2 += Form("[%d] ",curElem) ;
2007 while((arg=(RooAbsArg*)iter.next())) {
2008 arg->printCompactTree(os,indent2) ;
2009 }
2010
2011 if (curElem==maxElem) {
2012 os << indent << "RooProdPdf end partial integral cache" << endl ;
2013 }
2014}
2015
2016
2017
2018////////////////////////////////////////////////////////////////////////////////
2019/// Forward determination of safety of internal generator code to
2020/// component p.d.f that would generate the given observable
2021
2023{
2024 // Only override base class behaviour if default generator method is enabled
2025 if (!_useDefaultGen) return RooAbsPdf::isDirectGenSafe(arg) ;
2026
2027 // Argument may appear in only one PDF component
2028 RooAbsPdf* pdf, *thePdf(0) ;
2029 RooFIter pdfIter = _pdfList.fwdIterator();
2030 while((pdf=(RooAbsPdf*)pdfIter.next())) {
2031
2032 if (pdf->dependsOn(arg)) {
2033 // Found PDF depending on arg
2034
2035 // If multiple PDFs depend on arg directGen is not safe
2036 if (thePdf) return kFALSE ;
2037
2038 thePdf = pdf ;
2039 }
2040 }
2041 // Forward call to relevant component PDF
2042 return thePdf?(thePdf->isDirectGenSafe(arg)):kFALSE ;
2043}
2044
2045
2046
2047////////////////////////////////////////////////////////////////////////////////
2048/// Look up user specified normalization set for given input PDF component
2049
2051{
2052 Int_t idx = _pdfList.index(&pdf) ;
2053 if (idx<0) return 0 ;
2054 return (RooArgSet*) _pdfNSetList.At(idx) ;
2055}
2056
2057
2058
2059////////////////////////////////////////////////////////////////////////////////
2060/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams
2061/// The observables set is required to distinguish unambiguously p.d.f in terms
2062/// of observables and parameters, which are not constraints, and p.d.fs in terms
2063/// of parameters only, which can serve as constraints p.d.f.s
2064
2065RooArgSet* RooProdPdf::getConstraints(const RooArgSet& observables, RooArgSet& constrainedParams, Bool_t stripDisconnected) const
2066{
2067 RooArgSet constraints ;
2068 RooArgSet pdfParams, conParams ;
2069
2070 // Loop over p.d.f. components
2071 RooFIter piter = _pdfList.fwdIterator() ;
2072 RooAbsPdf* pdf ;
2073 while((pdf=(RooAbsPdf*)piter.next())) {
2074 // A constraint term is a p.d.f that does not depend on any of the listed observables
2075 // but does depends on any of the parameters that should be constrained
2076 if (!pdf->dependsOnValue(observables) && pdf->dependsOnValue(constrainedParams)) {
2077 constraints.add(*pdf) ;
2078 RooArgSet* tmp = pdf->getParameters(observables) ;
2079 conParams.add(*tmp,kTRUE) ;
2080 delete tmp ;
2081 } else {
2082 RooArgSet* tmp = pdf->getParameters(observables) ;
2083 pdfParams.add(*tmp,kTRUE) ;
2084 delete tmp ;
2085 }
2086 }
2087
2088 // Strip any constraints that are completely decoupled from the other product terms
2089 RooArgSet* finalConstraints = new RooArgSet("constraints") ;
2090 RooFIter citer = constraints.fwdIterator() ;
2091 while((pdf=(RooAbsPdf*)citer.next())) {
2092 if (pdf->dependsOnValue(pdfParams) || !stripDisconnected) {
2093 finalConstraints->add(*pdf) ;
2094 } else {
2095 coutI(Minimization) << "RooProdPdf::getConstraints(" << GetName() << ") omitting term " << pdf->GetName()
2096 << " as constraint term as it does not share any parameters with the other pdfs in product. "
2097 << "To force inclusion in likelihood, add an explicit Constrain() argument for the target parameter" << endl ;
2098 }
2099 }
2100
2101 // Now remove from constrainedParams all parameters that occur exclusively in constraint term and not in regular pdf term
2102
2103 RooArgSet* cexl = (RooArgSet*) conParams.selectCommon(constrainedParams) ;
2104 cexl->remove(pdfParams,kTRUE,kTRUE) ;
2105 constrainedParams.remove(*cexl,kTRUE,kTRUE) ;
2106 delete cexl ;
2107
2108 return finalConstraints ;
2109}
2110
2111
2112
2113
2114////////////////////////////////////////////////////////////////////////////////
2115/// Return all parameter constraint p.d.f.s on parameters listed in constrainedParams
2116/// The observables set is required to distinguish unambiguously p.d.f in terms
2117/// of observables and parameters, which are not constraints, and p.d.fs in terms
2118/// of parameters only, which can serve as constraints p.d.f.s
2119
2121{
2122 RooFIter iter = _pdfList.fwdIterator() ;
2123 RooAbsArg* arg ;
2124 RooArgSet* connectedPars = new RooArgSet("connectedPars") ;
2125 while((arg=iter.next())) {
2126 // Check if term is relevant
2127 if (arg->dependsOn(observables)) {
2128 RooArgSet* tmp = arg->getParameters(observables) ;
2129 connectedPars->add(*tmp) ;
2130 delete tmp ;
2131 } else {
2132 }
2133 }
2134 return connectedPars ;
2135}
2136
2137
2138
2139
2140////////////////////////////////////////////////////////////////////////////////
2141
2142void RooProdPdf::getParametersHook(const RooArgSet* nset, RooArgSet* params, Bool_t stripDisconnected) const
2143{
2144 if (!stripDisconnected) return ;
2145 if (!nset || nset->getSize()==0) return ;
2146
2147 // Get/create appropriate term list for this normalization set
2148 RooArgList *plist(0) ;
2149 RooLinkedList *nlist(0) ;
2150 Int_t code ;
2151 getPartIntList(nset,0,plist,nlist,code) ;
2152
2153 // Strip any terms from params that do not depend on any term
2154 RooFIter piter = params->fwdIterator() ;
2155 RooAbsReal* term, *param ;
2156 RooArgSet tostrip ;
2157 while((param=(RooAbsReal*)piter.next())) {
2158 Bool_t anyDep(kFALSE) ;
2159 RooFIter titer = plist->fwdIterator() ;
2160 while((term=(RooAbsReal*)titer.next())) {
2161 if (term->dependsOnValue(*param)) {
2162 anyDep=kTRUE ;
2163 }
2164 }
2165 if (!anyDep) {
2166 tostrip.add(*param) ;
2167 }
2168 }
2169
2170 if (tostrip.getSize()>0) {
2171 params->remove(tostrip,kTRUE,kTRUE);
2172 }
2173
2174}
2175
2176
2177
2178////////////////////////////////////////////////////////////////////////////////
2179/// Interface function used by test statistics to freeze choice of range
2180/// for interpretation of conditional product terms
2181
2182void RooProdPdf::selectNormalizationRange(const char* rangeName, Bool_t force)
2183{
2184 if (!force && _refRangeName) {
2185 return ;
2186 }
2187
2188 fixRefRange(rangeName) ;
2189}
2190
2191
2192
2193
2194////////////////////////////////////////////////////////////////////////////////
2195
2196void RooProdPdf::fixRefRange(const char* rangeName)
2197{
2198 _refRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
2199}
2200
2201
2202
2203////////////////////////////////////////////////////////////////////////////////
2204/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2205
2206std::list<Double_t>* RooProdPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
2207{
2208 RooAbsPdf* pdf ;
2209 RooFIter pdfIter = _pdfList.fwdIterator();
2210 while((pdf=(RooAbsPdf*)pdfIter.next())) {
2211 list<Double_t>* hint = pdf->plotSamplingHint(obs,xlo,xhi) ;
2212 if (hint) {
2213 return hint ;
2214 }
2215 }
2216
2217 return 0 ;
2218}
2219
2220
2221
2222////////////////////////////////////////////////////////////////////////////////
2223/// If all components that depend on obs are binned that so is the product
2224
2226{
2227 RooAbsPdf* pdf ;
2228 RooFIter pdfIter = _pdfList.fwdIterator();
2229 while((pdf=(RooAbsPdf*)pdfIter.next())) {
2230 if (pdf->dependsOn(obs) && !pdf->isBinnedDistribution(obs)) {
2231 return kFALSE ;
2232 }
2233 }
2234
2235 return kTRUE ;
2236}
2237
2238
2239
2240
2241
2242
2243////////////////////////////////////////////////////////////////////////////////
2244/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
2245
2246std::list<Double_t>* RooProdPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
2247{
2248 RooAbsPdf* pdf ;
2249 RooFIter pdfIter = _pdfList.fwdIterator();
2250 while((pdf=(RooAbsPdf*)pdfIter.next())) {
2251 list<Double_t>* hint = pdf->binBoundaries(obs,xlo,xhi) ;
2252 if (hint) {
2253 return hint ;
2254 }
2255 }
2256
2257 return 0 ;
2258}
2259
2260
2261////////////////////////////////////////////////////////////////////////////////
2262/// Label OK'ed components of a RooProdPdf with cache-and-track, _and_ label all RooProdPdf
2263/// descendants with extra informationa about (conditional) normalization, needed to be able
2264/// to Cache-And-Track them outside the RooprodPdf context
2265
2267{
2268 RooFIter piter = pdfList().fwdIterator() ;
2269 RooAbsArg* parg ;
2270 while ((parg=piter.next())) {
2271
2272 if (parg->canNodeBeCached()==Always) {
2273 trackNodes.add(*parg) ;
2274 //cout << "tracking node RooProdPdf component " << parg->IsA()->GetName() << "::" << parg->GetName() << endl ;
2275
2276 // Additional processing to fix normalization sets in case product defines conditional observables
2277 RooArgSet* pdf_nset = findPdfNSet((RooAbsPdf&)(*parg)) ;
2278 if (pdf_nset) {
2279 // Check if conditional normalization is specified
2280 if (string("nset")==pdf_nset->GetName() && pdf_nset->getSize()>0) {
2281 RooNameSet n(*pdf_nset) ;
2282 parg->setStringAttribute("CATNormSet",n.content()) ;
2283 }
2284 if (string("cset")==pdf_nset->GetName()) {
2285 RooNameSet c(*pdf_nset) ;
2286 parg->setStringAttribute("CATCondSet",c.content()) ;
2287 }
2288 } else {
2289 coutW(Optimization) << "RooProdPdf::setCacheAndTrackHints(" << GetName() << ") WARNING product pdf does not specify a normalization set for component " << parg->GetName() << endl ;
2290 }
2291 }
2292 }
2293}
2294
2295
2296
2297////////////////////////////////////////////////////////////////////////////////
2298/// Customized printing of arguments of a RooProdPdf to more intuitively reflect the contents of the
2299/// product operator construction
2300
2301void RooProdPdf::printMetaArgs(ostream& os) const
2302{
2304 for (int i=0 ; i<_pdfList.getSize() ; i++) {
2305 if (i>0) os << " * " ;
2306 RooArgSet* ncset = (RooArgSet*) niter.next() ;
2307 os << _pdfList.at(i)->GetName() ;
2308 if (ncset->getSize()>0) {
2309 if (string("nset")==ncset->GetName()) {
2310 os << *ncset ;
2311 } else {
2312 os << "|" ;
2313 RooFIter nciter = ncset->fwdIterator() ;
2314 RooAbsArg* arg ;
2315 Bool_t first(kTRUE) ;
2316 while((arg=(RooAbsArg*)nciter.next())) {
2317 if (!first) {
2318 os << "," ;
2319 } else {
2320 first = kFALSE ;
2321 }
2322 os << arg->GetName() ;
2323 }
2324 }
2325 }
2326 }
2327 os << " " ;
2328}
2329
2330
2331
2332////////////////////////////////////////////////////////////////////////////////
2333/// Implement support for node removal
2334
2335Bool_t RooProdPdf::redirectServersHook(const RooAbsCollection& /*newServerList*/, Bool_t /*mustReplaceAll*/, Bool_t nameChange, Bool_t /*isRecursive*/)
2336{
2337 if (nameChange && _pdfList.find("REMOVAL_DUMMY")) {
2338
2339 cxcoutD(LinkStateMgmt) << "RooProdPdf::redirectServersHook(" << GetName() << "): removing REMOVAL_DUMMY" << endl ;
2340
2341 // Remove node from _pdfList proxy and remove corresponding entry from normset list
2342 RooAbsArg* pdfDel = _pdfList.find("REMOVAL_DUMMY") ;
2343
2344 TObject* setDel = _pdfNSetList.At(_pdfList.index("REMOVAL_DUMMY")) ;
2345 _pdfList.remove(*pdfDel) ;
2346 _pdfNSetList.Remove(setDel) ;
2347
2348 // Clear caches
2349 _cacheMgr.reset() ;
2350 }
2351 return kFALSE ;
2352}
void Class()
Definition: Class.C:29
#define c(i)
Definition: RSha256.hxx:101
#define coutI(a)
Definition: RooMsgService.h:31
#define cxcoutD(a)
Definition: RooMsgService.h:79
#define coutW(a)
Definition: RooMsgService.h:33
#define dologD(a)
Definition: RooMsgService.h:63
#define coutE(a)
Definition: RooMsgService.h:34
#define TRACE_DESTROY
Definition: RooTrace.h:23
#define TRACE_CREATE
Definition: RooTrace.h:22
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 * R__STRTOK_R(char *str, const char *delim, char **saveptr)
Definition: Rtypes.h:484
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...
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
RooExpensiveObjectCache & expensiveObjectCache() const
Definition: RooAbsArg.cxx:2334
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
void printCompactTree(const char *indent="", const char *fileName=0, const char *namePat=0, RooAbsArg *client=0)
Print tree structure of expression tree on stdout, or to file if filename is specified.
Definition: RooAbsArg.cxx:1775
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
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
Definition: RooAbsArg.cxx:273
virtual CacheMode canNodeBeCached() const
Definition: RooAbsArg.h:317
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:75
friend class RooArgSet
Definition: RooAbsArg.h:471
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition: RooAbsArg.h:499
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
Definition: RooAbsArg.cxx:286
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
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2274
Bool_t getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
Definition: RooAbsArg.cxx:264
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Definition: RooAbsArg.h:88
RooArgSet * getComponents() const
Definition: RooAbsArg.cxx:680
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) 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.
RooAbsArg * first() const
Bool_t overlaps(const RooAbsCollection &otherColl) const
Check if this and other collection have common entries.
const char * GetName() const
Returns name of object.
void setName(const char *name)
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual void printValue(std::ostream &os) const
Print value of collection, i.e.
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 void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code.
Definition: RooAbsPdf.cxx:2069
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Interface function to create a generator context from a p.d.f.
Definition: RooAbsPdf.cxx:1666
TString _normRange
MC generator configuration specific for this object.
Definition: RooAbsPdf.h:353
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
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
friend class RooRealIntegral
Definition: RooAbsPdf.h:308
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in 'nset'.
Definition: RooAbsPdf.cxx:253
@ CanNotBeExtended
Definition: RooAbsPdf.h:223
const char * normRange() const
Definition: RooAbsPdf.h:257
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
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2057
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual Double_t getMax(const char *name=0) const
virtual Double_t getMin(const char *name=0) const
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
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:280
TString integralNameSuffix(const RooArgSet &iset, const RooArgSet *nset=0, const char *rangeName=0, Bool_t omitEmpty=kFALSE) const
Construct string with unique suffix name to give to integral object that encodes integrated observabl...
Definition: RooAbsReal.cxx:756
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:279
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
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:278
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition: RooAddition.h:26
const RooArgList & list1() const
Definition: RooAddition.h:41
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
Int_t index(const RooAbsArg *arg) const
Definition: RooArgList.h:76
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)
const RooNameSet * nameSet2ByIndex(Int_t index) const
T * getObjByIndex(Int_t index) const
Int_t lastIndex() const
const RooNameSet * nameSet1ByIndex(Int_t index) const
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:27
const RooArgSet * getSet(Int_t idx) const
Return RooArgSet stored in slot idx.
Definition: RooCmdArg.cxx:198
Int_t getInt(Int_t idx) const
Definition: RooCmdArg.h:80
virtual TObject * clone(const char *newname) const
Definition: RooConstVar.h:33
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
RooAbsArg * next()
RooAbsArg * getParameter(const char *name) const
Definition: RooFormulaVar.h:39
RooGenProdProj is an auxiliary class for RooProdPdf that calculates a general normalized projection o...
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:35
RooFIter fwdIterator() const
TObject * At(Int_t index) const
Return object stored in sequential position given by index.
void Delete(Option_t *o=0)
Remove all elements in collection and delete all elements NB: Collection does not own elements,...
virtual void Add(TObject *arg)
Definition: RooLinkedList.h:62
virtual Bool_t Remove(TObject *arg)
Remove object from collection.
Int_t IndexOf(const char *name) const
Return position of given object in list.
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()
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
RooNameSet is a utility class that stores the names the objects in a RooArget.
Definition: RooNameSet.h:24
RooArgSet * select(const RooArgSet &list) const
Construct a RooArgSet of objects in input 'list' whose names match to those in the internal name list...
Definition: RooNameSet.cxx:187
virtual ~CacheElem()
Destructor.
virtual void printCompactTreeHook(std::ostream &, const char *, Int_t, Int_t)
Hook function to print cache contents in tree printing of RooProdPdf.
std::unique_ptr< RooAbsReal > _rearrangedNum
Definition: RooProdPdf.h:141
RooArgList _ownedList
Definition: RooProdPdf.h:138
RooArgList _partList
Definition: RooProdPdf.h:135
std::unique_ptr< RooAbsReal > _rearrangedDen
Definition: RooProdPdf.h:142
virtual RooArgList containedArgs(Action)
Return RooAbsArg components contained in the cache.
RooLinkedList _normList
Definition: RooProdPdf.h:139
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
void rearrangeProduct(CacheElem &) const
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return the expected number of events associated with the extentable input p.d.f in the product.
void factorizeProduct(const RooArgSet &normSet, const RooArgSet &intSet, RooLinkedList &termList, RooLinkedList &normList, RooLinkedList &impDepList, RooLinkedList &crossDepList, RooLinkedList &intList) const
Factorize product in irreducible terms for given choice of integration/normalization.
Definition: RooProdPdf.cxx:545
void groupProductTerms(RooLinkedList &groupedTerms, RooArgSet &outerIntDeps, const RooLinkedList &terms, const RooLinkedList &norms, const RooLinkedList &imps, const RooLinkedList &ints, const RooLinkedList &cross) const
Group product into terms that can be calculated independently.
RooArgSet * _curNormSet
Registry of composite direct generator codes.
Definition: RooProdPdf.h:165
Int_t _extendedIndex
Definition: RooProdPdf.h:169
RooProdPdf()
Default constructor.
Definition: RooProdPdf.cxx:82
virtual Bool_t checkObservables(const RooArgSet *nset) const
Obsolete.
Double_t _cutOff
Definition: RooProdPdf.h:166
void printMetaArgs(std::ostream &os) const
Customized printing of arguments of a RooProdPdf to more intuitively reflect the contents of the prod...
Double_t evaluate() const
Calculate current value of object.
Definition: RooProdPdf.cxx:458
TNamed * _refRangeName
Definition: RooProdPdf.h:174
RooAICRegistry _genCode
Definition: RooProdPdf.h:163
RooListProxy _pdfList
Definition: RooProdPdf.h:167
Bool_t isBinnedDistribution(const RooArgSet &obs) const
If all components that depend on obs are binned that so is the product.
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Forward determination of safety of internal generator code to component p.d.f that would generate the...
Bool_t _useDefaultGen
Definition: RooProdPdf.h:172
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=0) const
Determine which part (if any) of given integral can be performed analytically.
RooAbsReal * specializeIntegral(RooAbsReal &orig, const char *targetRangeName) const
virtual ExtendMode extendMode() const
If this product contains exactly one extendable p.d.f return the extension abilities of that p....
RooObjCacheManager _cacheMgr
Definition: RooProdPdf.h:149
virtual void initGenerator(Int_t code)
Forward one-time initialization call to component generation initialization methods.
std::string makeRGPPName(const char *pfx, const RooArgSet &term, const RooArgSet &iset, const RooArgSet &nset, const char *isetRangeName) const
Make an appropriate automatic name for a RooGenProdProj object in getPartIntList()
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Return analytical integral defined by given scenario code.
virtual ~RooProdPdf()
Destructor.
Definition: RooProdPdf.cxx:436
virtual void getParametersHook(const RooArgSet *, RooArgSet *, Bool_t stripDisconnected) const
virtual Double_t getValV(const RooArgSet *set=0) const
Overload getVal() to intercept normalization set for use in evaluate()
Definition: RooProdPdf.cxx:447
RooAbsReal * makeCondPdfRatioCorr(RooAbsReal &term, const RooArgSet &termNset, const RooArgSet &termImpSet, const char *normRange, const char *refRange) const
For single normalization ranges.
RooArgSet * findPdfNSet(RooAbsPdf &pdf) const
Look up user specified normalization set for given input PDF component.
virtual void setCacheAndTrackHints(RooArgSet &)
Label OK'ed components of a RooProdPdf with cache-and-track, and label all RooProdPdf descendants wit...
virtual RooArgSet * getConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, Bool_t stripDisconnected) const
Return all parameter constraint p.d.f.s on parameters listed in constrainedParams The observables set...
void getPartIntList(const RooArgSet *nset, const RooArgSet *iset, pRooArgList &partList, pRooLinkedList &nsetList, Int_t &code, const char *isetRangeName=0) const
Return list of (partial) integrals of product terms for integration of p.d.f over observables iset wh...
Definition: RooProdPdf.cxx:730
RooLinkedList _pdfNSetList
Definition: RooProdPdf.h:168
friend class RooProdGenContext
Definition: RooProdPdf.h:158
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Query internal generation capabilities of component p.d.f.s and aggregate capabilities into master co...
virtual void selectNormalizationRange(const char *rangeName=0, Bool_t force=kFALSE)
Interface function used by test statistics to freeze choice of range for interpretation of conditiona...
RooArgSet * getConnectedParameters(const RooArgSet &observables) const
Return all parameter constraint p.d.f.s on parameters listed in constrainedParams The observables set...
Double_t calculate(const RooProdPdf::CacheElem &cache, Bool_t verbose=kFALSE) const
Calculate running product of pdfs terms, using the supplied normalization set in 'normSetList' for ea...
Definition: RooProdPdf.cxx:507
RooArgSet _defNormSet
Definition: RooProdPdf.h:177
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
RooAbsReal * specializeRatio(RooFormulaVar &input, const char *targetRangeName) const
Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Implement support for node removal.
const RooArgList & pdfList() const
Definition: RooProdPdf.h:69
void fixRefRange(const char *rangeName)
virtual Bool_t forceAnalyticalInt(const RooAbsArg &dep) const
Force RooRealIntegral to offer all observables for internal integration.
virtual TObject * clone(const char *newname) const
Definition: RooProdPdf.h:53
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Return generator context optimized for generating events from product p.d.f.s.
void initializeFromCmdArgList(const RooArgSet &fullPdfSet, const RooLinkedList &l)
Initialize RooProdPdf configuration from given list of RooCmdArg configuration arguments and set of '...
Definition: RooProdPdf.cxx:369
virtual void generateEvent(Int_t code)
Generate a single event with configuration specified by 'code' Defer internal generation to component...
std::vector< RooAbsReal * > processProductTerm(const RooArgSet *nset, const RooArgSet *iset, const char *isetRangeName, const RooArgSet *term, const RooArgSet &termNSet, const RooArgSet &termISet, Bool_t &isOwned, Bool_t forceWrap=kFALSE) const
Calculate integrals of factorized product terms over observables iset while normalized to observables...
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
RooProduct a RooAbsReal implementation that represent the product of a given set of other RooAbsReal ...
Definition: RooProduct.h:32
RooArgList components()
Definition: RooProduct.h:47
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooArgSet intVars() const
const RooAbsReal & integrand() const
const char * intRange()
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
const Int_t n
Definition: legend1.C:16
R__EXTERN PyObject * gIter
Definition: PyStrings.h:28
@ Minimization
Definition: RooGlobalFunc.h:57
@ Generation
Definition: RooGlobalFunc.h:57
@ Optimization
Definition: RooGlobalFunc.h:58
@ InputArguments
Definition: RooGlobalFunc.h:58
@ LinkStateMgmt
Definition: RooGlobalFunc.h:57
RooConstVar & RooConst(Double_t val)
Definition: first.py:1
STL namespace.
auto * l
Definition: textangle.C:4