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