Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooRealIntegral.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 /**
18 \file RooRealIntegral.cxx
19 \class RooRealIntegral
20 \ingroup Roofitcore
21 
22 RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
23 The class performs none of the actual integration, but only manages the logic
24 of what variables can be integrated analytically, accounts for eventual jacobian
25 terms and defines what numerical integrations needs to be done to complement the
26 analytical integral.
27 The actual analytical integrations (if any) are done in the PDF themselves, the numerical
28 integration is performed in the various implementations of the RooAbsIntegrator base class.
29 **/
30 
31 #include "RooFit.h"
32 
33 #include "TClass.h"
34 #include "RooMsgService.h"
35 #include "Riostream.h"
36 #include "TObjString.h"
37 #include "TH1.h"
38 #include "RooRealIntegral.h"
39 #include "RooArgSet.h"
40 #include "RooAbsRealLValue.h"
41 #include "RooAbsCategoryLValue.h"
42 #include "RooRealBinding.h"
43 #include "RooRealAnalytic.h"
44 #include "RooInvTransform.h"
45 #include "RooSuperCategory.h"
46 #include "RooNumIntFactory.h"
47 #include "RooNumIntConfig.h"
48 #include "RooNameReg.h"
50 #include "RooConstVar.h"
51 #include "RooDouble.h"
52 #include "RooTrace.h"
53 
54 using namespace std;
55 
56 ClassImp(RooRealIntegral);
57 ;
58 
59 
60 Int_t RooRealIntegral::_cacheAllNDim(2) ;
61 
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 
65 RooRealIntegral::RooRealIntegral() :
66  _valid(kFALSE),
67  _respectCompSelect(true),
68  _funcNormSet(0),
69  _iconfig(0),
70  _sumCatIter(0),
71  _mode(0),
72  _intOperMode(Hybrid),
73  _restartNumIntEngine(kFALSE),
74  _numIntEngine(0),
75  _numIntegrand(0),
76  _rangeName(0),
77  _params(0),
78  _cacheNum(kFALSE)
79 {
80  TRACE_CREATE
81 }
82 
83 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Construct integral of 'function' over observables in 'depList'
87 /// in range 'rangeName' with normalization observables 'funcNormSet'
88 /// (for p.d.f.s). In the integral is performed to the maximum extent
89 /// possible the internal (analytical) integrals advertised by function.
90 /// The other integrations are performed numerically. The optional
91 /// config object prescribes how these numeric integrations are configured.
92 ///
93 
94 RooRealIntegral::RooRealIntegral(const char *name, const char *title,
95  const RooAbsReal& function, const RooArgSet& depList,
96  const RooArgSet* funcNormSet, const RooNumIntConfig* config,
97  const char* rangeName) :
98  RooAbsReal(name,title),
99  _valid(kTRUE),
100  _respectCompSelect(true),
101  _sumList("!sumList","Categories to be summed numerically",this,kFALSE,kFALSE),
102  _intList("!intList","Variables to be integrated numerically",this,kFALSE,kFALSE),
103  _anaList("!anaList","Variables to be integrated analytically",this,kFALSE,kFALSE),
104  _jacList("!jacList","Jacobian product term",this,kFALSE,kFALSE),
105  _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
106  _function("!func","Function to be integrated",this,
107  const_cast<RooAbsReal&>(function),kFALSE,kFALSE),
108  _iconfig((RooNumIntConfig*)config),
109  _sumCat("!sumCat","SuperCategory for summation",this,kFALSE,kFALSE),
110  _sumCatIter(0),
111  _mode(0),
112  _intOperMode(Hybrid),
113  _restartNumIntEngine(kFALSE),
114  _numIntEngine(0),
115  _numIntegrand(0),
116  _rangeName((TNamed*)RooNameReg::ptr(rangeName)),
117  _params(0),
118  _cacheNum(kFALSE)
119 {
120  // A) Check that all dependents are lvalues
121  //
122  // B) Check if list of dependents can be re-expressed in
123  // lvalues that are higher in the expression tree
124  //
125  // C) Check for dependents that the PDF insists on integrating
126  // analytically iself
127  //
128  // D) Make list of servers that can be integrated analytically
129  // Add all parameters/dependents as value/shape servers
130  //
131  // E) Interact with function to make list of objects actually integrated analytically
132  //
133  // F) Make list of numerical integration variables consisting of:
134  // - Category dependents of RealLValues in analytical integration
135  // - Leaf nodes server lists of function server that are not analytically integrated
136  // - Make Jacobian list for analytically integrated RealLValues
137  //
138  // G) Split numeric list in integration list and summation list
139  //
140 
141  oocxcoutI(&function,Integration) << "RooRealIntegral::ctor(" << GetName() << ") Constructing integral of function "
142  << function.GetName() << " over observables" << depList << " with normalization "
143  << (funcNormSet?*funcNormSet:RooArgSet()) << " with range identifier "
144  << (rangeName?rangeName:"<none>") << endl ;
145 
146 
147  // Choose same expensive object cache as integrand
148  setExpensiveObjectCache(function.expensiveObjectCache()) ;
149 // cout << "RRI::ctor(" << GetName() << ") setting expensive object cache to " << &expensiveObjectCache() << " as taken from " << function.GetName() << endl ;
150 
151  // Use objects integrator configuration if none is specified
152  if (!_iconfig) _iconfig = (RooNumIntConfig*) function.getIntegratorConfig() ;
153 
154  // Save private copy of funcNormSet, if supplied, excluding factorizing terms
155  if (funcNormSet) {
156  _funcNormSet = new RooArgSet ;
157  for (const auto nArg : *funcNormSet) {
158  if (function.dependsOn(*nArg)) {
159  _funcNormSet->addClone(*nArg) ;
160  }
161  }
162  } else {
163  _funcNormSet = 0 ;
164  }
165 
166  //_funcNormSet = funcNormSet ? (RooArgSet*)funcNormSet->snapshot(kFALSE) : 0 ;
167 
168  // Make internal copy of dependent list
169  RooArgSet intDepList(depList) ;
170 
171  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
172  // * A) Check that all dependents are lvalues and filter out any
173  // dependents that the PDF doesn't explicitly depend on
174  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
175 
176  for (auto arg : intDepList) {
177  if(!arg->isLValue()) {
178  coutE(InputArguments) << ClassName() << "::" << GetName() << ": cannot integrate non-lvalue ";
179  arg->Print("1");
180  _valid= kFALSE;
181  }
182  if (!function.dependsOn(*arg)) {
183  RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
184  _facListOwned.addOwned(*argClone) ;
185  _facList.add(*argClone) ;
186  addServer(*argClone,kFALSE,kTRUE) ;
187  }
188  }
189 
190  if (_facList.getSize()>0) {
191  oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables are " << _facList << endl ;
192  }
193 
194 
195  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
196  // * B) Check if list of dependents can be re-expressed in *
197  // * lvalues that are higher in the expression tree *
198  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
199 
200 
201  // Initial fill of list of LValue branches
202  RooArgSet exclLVBranches("exclLVBranches") ;
203  RooArgSet branchList,branchListVD ;
204  function.branchNodeServerList(&branchList) ;
205 
206  for (auto branch: branchList) {
207  RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(branch) ;
208  RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(branch) ;
209  if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
210  exclLVBranches.add(*branch) ;
211 // cout << "exclv branch = " << endl ;
212 // branch->printCompactTree() ;
213  }
214  if (dependsOnValue(*branch)) {
215  branchListVD.add(*branch) ;
216  } else {
217 // cout << "value of self does not depend on branch " << branch->GetName() << endl ;
218  }
219  }
220  exclLVBranches.remove(depList,kTRUE,kTRUE) ;
221 // cout << "exclLVBranches = " << exclLVBranches << endl ;
222 
223  // Initial fill of list of LValue leaf servers (put in intDepList)
224  RooArgSet exclLVServers("exclLVServers") ;
225  exclLVServers.add(intDepList) ;
226 
227 // cout << "begin exclLVServers = " << exclLVServers << endl ;
228 
229  // Obtain mutual exclusive dependence by iterative reduction
230  Bool_t converged(kFALSE) ;
231  while(!converged) {
232  converged=kTRUE ;
233 
234  // Reduce exclLVServers to only those serving exclusively exclLVBranches
235  std::vector<RooAbsArg*> toBeRemoved;
236  for (auto server : exclLVServers) {
237  if (!servesExclusively(server,exclLVBranches,branchListVD)) {
238  toBeRemoved.push_back(server);
239  // cout << "removing " << server->GetName() << " from exclLVServers because servesExclusively(" << server->GetName() << "," << exclLVBranches << ") faile" << endl ;
240  converged=kFALSE ;
241  }
242  exclLVServers.remove(toBeRemoved.begin(), toBeRemoved.end());
243  }
244 
245  // Reduce exclLVBranches to only those depending exclusisvely on exclLVservers
246  for (auto branch : exclLVBranches) {
247  RooArgSet* brDepList = branch->getObservables(&intDepList) ;
248  RooArgSet bsList(*brDepList,"bsList") ;
249  delete brDepList ;
250  bsList.remove(exclLVServers,kTRUE,kTRUE) ;
251  if (bsList.getSize()>0) {
252  exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
253 // cout << "removing " << branch->GetName() << " from exclLVBranches" << endl ;
254  converged=kFALSE ;
255  }
256  }
257  }
258 
259  // Eliminate exclLVBranches that do not depend on any LVServer
260  for (auto branch : exclLVBranches) {
261  if (!branch->dependsOnValue(exclLVServers)) {
262  //cout << "LV branch " << branch->GetName() << " does not depend on any LVServer (" << exclLVServers << ") and will be removed" << endl ;
263  exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
264  }
265  }
266 
267 // cout << "end exclLVServers = " << exclLVServers << endl ;
268 
269  // Replace exclusive lvalue branch servers with lvalue branches
270  // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
271  if (exclLVServers.getSize()>0 && !function.isBinnedDistribution(exclLVBranches)) {
272 // cout << "activating LVservers " << exclLVServers << " for use in integration " << endl ;
273  intDepList.remove(exclLVServers) ;
274  intDepList.add(exclLVBranches) ;
275 
276  //cout << "intDepList removing exclLVServers " << exclLVServers << endl ;
277  //cout << "intDepList adding exclLVBranches " << exclLVBranches << endl ;
278 
279  }
280 
281 
282  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
283  // * C) Check for dependents that the PDF insists on integrating *
284  // analytically iself *
285  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
286 
287  RooArgSet anIntOKDepList ;
288  for (auto arg : intDepList) {
289  if (function.forceAnalyticalInt(*arg)) {
290  anIntOKDepList.add(*arg) ;
291  }
292  }
293 
294  if (anIntOKDepList.getSize()>0) {
295  oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << endl ;
296  }
297 
298 
299  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
300  // * D) Make list of servers that can be integrated analytically *
301  // Add all parameters/dependents as value/shape servers *
302  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
303 
304  for (const auto arg : function.servers()) {
305 
306  //cout << "considering server" << arg->GetName() << endl ;
307 
308  // Dependent or parameter?
309  if (!arg->dependsOnValue(intDepList)) {
310 
311  //cout << " server does not depend on observables, adding server as value server to integral" << endl ;
312 
313  if (function.dependsOnValue(*arg)) {
314  addServer(*arg,kTRUE,kFALSE) ;
315  }
316 
317  continue ;
318 
319  } else {
320 
321  // Add final dependents of arg as shape servers
322  RooArgSet argLeafServers ;
323  arg->leafNodeServerList(&argLeafServers,0,kFALSE) ;
324 
325  //arg->printCompactTree() ;
326  //cout << "leaf nodes of server are " << argLeafServers << " depList = " << depList << endl ;
327 
328  // Skip arg if it is neither value or shape server
329  if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
330  //cout << " server is neither value not shape server of function, ignoring" << endl ;
331  continue ;
332  }
333 
334  for (const auto leaf : argLeafServers) {
335 
336  //cout << " considering leafnode " << leaf->GetName() << " of server " << arg->GetName() << endl ;
337 
338  if (depList.find(leaf->GetName()) && function.dependsOnValue(*leaf)) {
339 
340  RooAbsRealLValue* leaflv = dynamic_cast<RooAbsRealLValue*>(leaf) ;
341  if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
342  oocxcoutD(&function,Integration) << function.GetName() << " : Observable " << leaf->GetName() << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf" << endl ;
343  if (leaflv->getBinning(rangeName).lowBoundFunc()) {
344  addServer(*leaflv->getBinning(rangeName).lowBoundFunc(),kTRUE,kFALSE) ;
345  }
346  if(leaflv->getBinning(rangeName).highBoundFunc()) {
347  addServer(*leaflv->getBinning(rangeName).highBoundFunc(),kTRUE,kFALSE) ;
348  }
349  } else {
350  oocxcoutD(&function,Integration) << function.GetName() << ": Adding observable " << leaf->GetName() << " of server "
351  << arg->GetName() << " as shape dependent" << endl ;
352  addServer(*leaf,kFALSE,kTRUE) ;
353  }
354  } else if (!depList.find(leaf->GetName())) {
355 
356  if (function.dependsOnValue(*leaf)) {
357  oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as value dependent" << endl ;
358  addServer(*leaf,kTRUE,kFALSE) ;
359  } else {
360  oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as shape dependent" << endl ;
361  addServer(*leaf,kFALSE,kTRUE) ;
362  }
363  }
364  }
365  }
366 
367  // If this dependent arg is self-normalized, stop here
368  //if (function.selfNormalized()) continue ;
369 
370  Bool_t depOK(kFALSE) ;
371  // Check for integratable AbsRealLValue
372 
373  //cout << "checking server " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
374 
375  if (arg->isDerived()) {
376  RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(arg) ;
377  RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(arg) ;
378  //cout << "realArgLV = " << realArgLV << " intDepList = " << intDepList << endl ;
379  if ((realArgLV && intDepList.find(realArgLV->GetName()) && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
380 
381  //cout << " arg " << arg->GetName() << " is derived LValue with valid jacobian" << endl ;
382 
383  // Derived LValue with valid jacobian
384  depOK = kTRUE ;
385 
386  // Now, check for overlaps
387  Bool_t overlapOK = kTRUE ;
388  for (const auto otherArg : function.servers()) {
389  // skip comparison with self
390  if (arg==otherArg) continue ;
391  if (otherArg->IsA()==RooConstVar::Class()) continue ;
392  if (arg->overlaps(*otherArg,kTRUE)) {
393  //cout << "arg " << arg->GetName() << " overlaps with " << otherArg->GetName() << endl ;
394  //overlapOK=kFALSE ;
395  }
396  }
397  // coverity[DEADCODE]
398  if (!overlapOK) depOK=kFALSE ;
399 
400  //cout << "overlap check returns OK=" << (depOK?"T":"F") << endl ;
401 
402  }
403  } else {
404  // Fundamental types are always OK
405  depOK = kTRUE ;
406  }
407 
408  // Add server to list of dependents that are OK for analytical integration
409  if (depOK) {
410  anIntOKDepList.add(*arg,kTRUE) ;
411  oocxcoutI(&function,Integration) << function.GetName() << ": Observable " << arg->GetName() << " is suitable for analytical integration (if supported by p.d.f)" << endl ;
412  }
413  }
414  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
415  // * E) interact with function to make list of objects actually integrated analytically *
416  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
417 
418  RooArgSet anIntDepList ;
419 
420  RooArgSet *anaSet = new RooArgSet( _anaList, Form("UniqueCloneOf_%s",_anaList.GetName()));
421  _mode = ((RooAbsReal&)_function.arg()).getAnalyticalIntegralWN(anIntOKDepList,*anaSet,_funcNormSet,RooNameReg::str(_rangeName)) ;
422  _anaList.removeAll() ;
423  _anaList.add(*anaSet);
424  delete anaSet;
425 
426  // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaListx
427  if (_mode==0) {
428  _anaList.removeAll() ;
429  }
430 
431  if (_mode!=0) {
432  oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << endl ;
433  }
434 
435 
436  // WVE kludge: synchronize dset for use in analyticalIntegral
437  function.getVal(_funcNormSet) ;
438 
439  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
440  // * F) Make list of numerical integration variables consisting of: *
441  // * - Category dependents of RealLValues in analytical integration *
442  // * - Expanded server lists of server that are not analytically integrated *
443  // * Make Jacobian list with analytically integrated RealLValues *
444  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
445 
446  RooArgSet numIntDepList ;
447 
448 
449  // Loop over actually analytically integrated dependents
450  for (const auto arg : _anaList) {
451 
452  // Process only derived RealLValues
453  if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class()) && arg->isDerived() && !arg->isFundamental()) {
454 
455  // Add to list of Jacobians to calculate
456  _jacList.add(*arg) ;
457 
458  // Add category dependent of LValueReal used in integration
459  auto argDepList = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
460  for (const auto argDep : *argDepList) {
461  if (argDep->IsA()->InheritsFrom(RooAbsCategoryLValue::Class()) && intDepList.contains(*argDep)) {
462  numIntDepList.add(*argDep,kTRUE) ;
463  }
464  }
465  }
466  }
467 
468 
469  // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
470  if (_anaList.getSize()==0) {
471  if (exclLVServers.getSize()>0) {
472  //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << endl ;
473  intDepList.remove(exclLVBranches) ;
474  intDepList.add(exclLVServers) ;
475  }
476  }
477  //cout << "NUMINT intDepList = " << intDepList << endl ;
478 
479  // Loop again over function servers to add remaining numeric integrations
480  for (const auto arg : function.servers()) {
481 
482  //cout << "processing server for numeric integration " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
483 
484  // Process only servers that are not treated analytically
485  if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
486 
487  // Process only derived RealLValues
488  if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
489  numIntDepList.add(*arg,kTRUE) ;
490  } else {
491 
492  // WVE this will only get the observables, but not l-value transformations
493  // Expand server in final dependents
494  auto argDeps = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
495 
496  if (argDeps->getSize()>0) {
497 
498  // Add final dependents, that are not forcibly integrated analytically,
499  // to numerical integration list
500  for (const auto dep : *argDeps) {
501  if (!_anaList.find(dep->GetName())) {
502  numIntDepList.add(*dep,kTRUE) ;
503  }
504  }
505  }
506  }
507  }
508  }
509 
510  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
511  // * G) Split numeric list in integration list and summation list *
512  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
513 
514  // Split numeric integration list in summation and integration lists
515  for (const auto arg : numIntDepList) {
516  if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
517  _intList.add(*arg) ;
518  } else if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
519  _sumList.add(*arg) ;
520  }
521  }
522 
523  if (_anaList.getSize()>0) {
524  oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << endl ;
525  }
526  if (_intList.getSize()>0) {
527  oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << endl ;
528  }
529  if (_sumList.getSize()>0) {
530  oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << endl ;
531  }
532 
533 
534  // Determine operating mode
535  if (numIntDepList.getSize()>0) {
536  // Numerical and optional Analytical integration
537  _intOperMode = Hybrid ;
538  } else if (_anaList.getSize()>0) {
539  // Purely analytical integration
540  _intOperMode = Analytic ;
541  } else {
542  // No integration performed
543  _intOperMode = PassThrough ;
544  }
545 
546  // Determine auto-dirty status
547  autoSelectDirtyMode() ;
548 
549  // Create value caches for _intList and _sumList
550  _intList.snapshot(_saveInt) ;
551  _sumList.snapshot(_saveSum) ;
552 
553 
554  if (_sumList.getSize()>0) {
555  RooSuperCategory *sumCat = new RooSuperCategory(Form("%s_sumCat",GetName()),"sumCat",_sumList) ;
556  _sumCatIter = sumCat->typeIterator() ;
557  _sumCat.addOwned(*sumCat) ;
558  }
559 
560  TRACE_CREATE
561 }
562 
563 
564 
565 ////////////////////////////////////////////////////////////////////////////////
566 /// Set appropriate cache operation mode for integral depending on cache operation
567 /// mode of server objects
568 
569 void RooRealIntegral::autoSelectDirtyMode()
570 {
571  // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
572  for (const auto server : _serverList) {
573  if (server->isValueServer(*this)) {
574  RooArgSet leafSet ;
575  server->leafNodeServerList(&leafSet) ;
576  for (const auto leaf : leafSet) {
577  if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
578  setOperMode(ADirty) ;
579  break ;
580  }
581  if (leaf->getAttribute("projectedDependent")) {
582  setOperMode(ADirty) ;
583  break ;
584  }
585  }
586  }
587  }
588 }
589 
590 
591 
592 ////////////////////////////////////////////////////////////////////////////////
593 /// Utility function that returns true if 'object server' is a server
594 /// to exactly one of the RooAbsArgs in 'exclLVBranches'
595 
596 Bool_t RooRealIntegral::servesExclusively(const RooAbsArg* server,const RooArgSet& exclLVBranches, const RooArgSet& allBranches) const
597 {
598  // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
599 
600  // Special case, no LV servers available
601  if (exclLVBranches.getSize()==0) return kFALSE ;
602 
603  // If server has no clients and is not an LValue itself, return false
604  if (server->_clientList.empty() && exclLVBranches.find(server->GetName())) {
605  return kFALSE ;
606  }
607 
608  // WVE must check for value relations only here!!!!
609 
610  // Loop over all clients
611  Int_t numLVServ(0) ;
612  for (const auto client : server->valueClients()) {
613  // If client is not an LValue, recurse
614  if (!(exclLVBranches.find(client->GetName())==client)) {
615  if (allBranches.find(client->GetName())==client) {
616  if (!servesExclusively(client,exclLVBranches,allBranches)) {
617  // Client is a non-LValue that doesn't have an exclusive LValue server
618  return kFALSE ;
619  }
620  }
621  } else {
622  // Client is an LValue
623  numLVServ++ ;
624  }
625  }
626 
627  return (numLVServ==1) ;
628 }
629 
630 
631 
632 
633 ////////////////////////////////////////////////////////////////////////////////
634 /// (Re)Initialize numerical integration engine if necessary. Return kTRUE if
635 /// successful, or otherwise kFALSE.
636 
637 Bool_t RooRealIntegral::initNumIntegrator() const
638 {
639  // if we already have an engine, check if it still works for the present limits.
640  if(0 != _numIntEngine) {
641  if(_numIntEngine->isValid() && _numIntEngine->checkLimits() && !_restartNumIntEngine ) return kTRUE;
642  // otherwise, cleanup the old engine
643  delete _numIntEngine ;
644  _numIntEngine= 0;
645  if(0 != _numIntegrand) {
646  delete _numIntegrand;
647  _numIntegrand= 0;
648  }
649  }
650 
651  // All done if there are no arguments to integrate numerically
652  if(0 == _intList.getSize()) return kTRUE;
653 
654  // Bind the appropriate analytic integral (specified by _mode) of our RooRealVar object to
655  // those of its arguments that will be integrated out numerically.
656  if(_mode != 0) {
657  _numIntegrand= new RooRealAnalytic(_function.arg(),_intList,_mode,_funcNormSet,_rangeName);
658  }
659  else {
660  _numIntegrand= new RooRealBinding(_function.arg(),_intList,_funcNormSet,kFALSE,_rangeName);
661  }
662  if(0 == _numIntegrand || !_numIntegrand->isValid()) {
663  coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << endl;
664  return kFALSE;
665  }
666 
667  // Create appropriate numeric integrator using factory
668  Bool_t isBinned = _function.arg().isBinnedDistribution(_intList) ;
669  _numIntEngine = RooNumIntFactory::instance().createIntegrator(*_numIntegrand,*_iconfig,0,isBinned) ;
670 
671  if(0 == _numIntEngine || !_numIntEngine->isValid()) {
672  coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << endl;
673  return kFALSE;
674  }
675 
676  cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
677  << _numIntEngine->IsA()->GetName() << " to calculate Int" << _intList << endl ;
678 
679  if (_intList.getSize()>3) {
680  cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") evaluation requires " << _intList.getSize() << "-D numeric integration step. Evaluation may be slow, sufficient numeric precision for fitting & minimization is not guaranteed" << endl ;
681  }
682 
683  _restartNumIntEngine = kFALSE ;
684  return kTRUE;
685 }
686 
687 
688 
689 ////////////////////////////////////////////////////////////////////////////////
690 /// Copy constructor
691 
692 RooRealIntegral::RooRealIntegral(const RooRealIntegral& other, const char* name) :
693  RooAbsReal(other,name),
694  _valid(other._valid),
695  _respectCompSelect(other._respectCompSelect),
696  _sumList("!sumList",this,other._sumList),
697  _intList("!intList",this,other._intList),
698  _anaList("!anaList",this,other._anaList),
699  _jacList("!jacList",this,other._jacList),
700  _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
701  _function("!func",this,other._function),
702  _iconfig(other._iconfig),
703  _sumCat("!sumCat",this,other._sumCat),
704  _sumCatIter(0),
705  _mode(other._mode),
706  _intOperMode(other._intOperMode),
707  _restartNumIntEngine(kFALSE),
708  _numIntEngine(0),
709  _numIntegrand(0),
710  _rangeName(other._rangeName),
711  _params(0),
712  _cacheNum(kFALSE)
713 {
714  _funcNormSet = other._funcNormSet ? (RooArgSet*)other._funcNormSet->snapshot(kFALSE) : 0 ;
715 
716  for (const auto arg : other._facList) {
717  RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
718  _facListOwned.addOwned(*argClone) ;
719  _facList.add(*argClone) ;
720  addServer(*argClone,kFALSE,kTRUE) ;
721  }
722 
723  other._intList.snapshot(_saveInt) ;
724  other._sumList.snapshot(_saveSum) ;
725 
726  TRACE_CREATE
727 }
728 
729 
730 
731 ////////////////////////////////////////////////////////////////////////////////
732 
733 RooRealIntegral::~RooRealIntegral()
734  // Destructor
735 {
736  if (_numIntEngine) delete _numIntEngine ;
737  if (_numIntegrand) delete _numIntegrand ;
738  if (_funcNormSet) delete _funcNormSet ;
739  if (_sumCatIter) delete _sumCatIter ;
740  if (_params) delete _params ;
741 
742  TRACE_DESTROY
743 }
744 
745 
746 
747 
748 
749 ////////////////////////////////////////////////////////////////////////////////
750 
751 RooAbsReal* RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const
752 {
753  // Handle special case of no integration with default algorithm
754  if (iset.getSize()==0) {
755  return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
756  }
757 
758  // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
759  RooArgSet isetAll(iset) ;
760  isetAll.add(_sumList) ;
761  isetAll.add(_intList) ;
762  isetAll.add(_anaList) ;
763  isetAll.add(_facList) ;
764 
765  const RooArgSet* newNormSet(0) ;
766  RooArgSet* tmp(0) ;
767  if (nset && !_funcNormSet) {
768  newNormSet = nset ;
769  } else if (!nset && _funcNormSet) {
770  newNormSet = _funcNormSet ;
771  } else if (nset && _funcNormSet) {
772  tmp = new RooArgSet ;
773  tmp->add(*nset) ;
774  tmp->add(*_funcNormSet,kTRUE) ;
775  newNormSet = tmp ;
776  }
777  RooAbsReal* ret = _function.arg().createIntegral(isetAll,newNormSet,cfg,rangeName) ;
778 
779  if (tmp) {
780  delete tmp ;
781  }
782 
783  return ret ;
784 }
785 
786 
787 
788 ////////////////////////////////////////////////////////////////////////////////
789 /// Return value of object. If the cache is clean, return the
790 /// cached value, otherwise recalculate on the fly and refill
791 /// the cache
792 
793 Double_t RooRealIntegral::getValV(const RooArgSet* nset) const
794 {
795 // // fast-track clean-cache processing
796 // if (_operMode==AClean) {
797 // return _value ;
798 // }
799 
800  if (nset && nset!=_lastNSet) {
801  ((RooAbsReal*) this)->setProxyNormSet(nset) ;
802  _lastNSet = (RooArgSet*) nset ;
803  }
804 
805  if (isValueOrShapeDirtyAndClear()) {
806  _value = traceEval(nset) ;
807  }
808 
809  return _value ;
810 }
811 
812 
813 
814 
815 
816 ////////////////////////////////////////////////////////////////////////////////
817 /// Perform the integration and return the result
818 
819 Double_t RooRealIntegral::evaluate() const
820 {
821  GlobalSelectComponentRAII selCompRAII(_globalSelectComp || !_respectCompSelect);
822 
823  Double_t retVal(0) ;
824  switch (_intOperMode) {
825 
826  case Hybrid:
827  {
828  // Cache numeric integrals in >1d expensive object cache
829  RooDouble* cacheVal(0) ;
830  if ((_cacheNum && _intList.getSize()>0) || _intList.getSize()>=_cacheAllNDim) {
831  cacheVal = (RooDouble*) expensiveObjectCache().retrieveObject(GetName(),RooDouble::Class(),parameters()) ;
832  }
833 
834  if (cacheVal) {
835  retVal = *cacheVal ;
836  // cout << "using cached value of integral" << GetName() << endl ;
837  } else {
838 
839 
840  // Find any function dependents that are AClean
841  // and switch them temporarily to ADirty
842  Bool_t origState = inhibitDirty() ;
843  setDirtyInhibit(kTRUE) ;
844 
845  // try to initialize our numerical integration engine
846  if(!(_valid= initNumIntegrator())) {
847  coutE(Integration) << ClassName() << "::" << GetName()
848  << ":evaluate: cannot initialize numerical integrator" << endl;
849  return 0;
850  }
851 
852  // Save current integral dependent values
853  _saveInt = _intList ;
854  _saveSum = _sumList ;
855 
856  // Evaluate sum/integral
857  retVal = sum() ;
858 
859 
860  // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
861  setDirtyInhibit(origState) ;
862 
863  // Restore integral dependent values
864  _intList=_saveInt ;
865  _sumList=_saveSum ;
866 
867  // Cache numeric integrals in >1d expensive object cache
868  if ((_cacheNum && _intList.getSize()>0) || _intList.getSize()>=_cacheAllNDim) {
869  RooDouble* val = new RooDouble(retVal) ;
870  expensiveObjectCache().registerObject(_function.arg().GetName(),GetName(),*val,parameters()) ;
871  // cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << endl ;
872  }
873 
874  }
875  break ;
876  }
877  case Analytic:
878  {
879  retVal = ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) / jacobianProduct() ;
880  cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
881  << ")func = " << _function.arg().IsA()->GetName() << "::" << _function.arg().GetName()
882  << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << endl ;
883 
884 
885  break ;
886  }
887 
888  case PassThrough:
889  {
890  //setDirtyInhibit(kTRUE) ;
891  retVal= _function.arg().getVal(_funcNormSet) ;
892  //setDirtyInhibit(kFALSE) ;
893  break ;
894  }
895  }
896 
897 
898  // Multiply answer with integration ranges of factorized variables
899  if (_facList.getSize()>0) {
900  for (const auto arg : _facList) {
901  // Multiply by fit range for 'real' dependents
902  if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
903  RooAbsRealLValue* argLV = (RooAbsRealLValue*)arg ;
904  retVal *= (argLV->getMax() - argLV->getMin()) ;
905  }
906  // Multiply by number of states for category dependents
907  if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
908  RooAbsCategoryLValue* argLV = (RooAbsCategoryLValue*)arg ;
909  retVal *= argLV->numTypes() ;
910  }
911  }
912  }
913 
914 
915  if (dologD(Tracing)) {
916  cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
917  switch(_mode) {
918  case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
919  case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
920  case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
921  }
922 
923  ccxcoutD(Tracing) << "raw*fact = " << retVal << endl ;
924  }
925 
926  return retVal ;
927 }
928 
929 
930 
931 ////////////////////////////////////////////////////////////////////////////////
932 /// Return product of jacobian terms originating from analytical integration
933 
934 Double_t RooRealIntegral::jacobianProduct() const
935 {
936  if (_jacList.getSize()==0) {
937  return 1 ;
938  }
939 
940  Double_t jacProd(1) ;
941  for (const auto elm : _jacList) {
942  auto arg = static_cast<const RooAbsRealLValue*>(elm);
943  jacProd *= arg->jacobian() ;
944  }
945 
946  // Take fabs() here: if jacobian is negative, min and max are swapped and analytical integral
947  // will be positive, so must multiply with positive jacobian.
948  return fabs(jacProd) ;
949 }
950 
951 
952 
953 ////////////////////////////////////////////////////////////////////////////////
954 /// Perform summation of list of category dependents to be integrated
955 
956 Double_t RooRealIntegral::sum() const
957 {
958  if (_sumList.getSize()!=0) {
959  // Add integrals for all permutations of categories summed over
960  Double_t total(0) ;
961 
962  _sumCatIter->Reset() ;
963  RooCatType* type ;
964  RooSuperCategory* sumCat = (RooSuperCategory*) _sumCat.first() ;
965  while((type=(RooCatType*)_sumCatIter->Next())) {
966  sumCat->setIndex(type->getVal()) ;
967  if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
968  total += integrate() / jacobianProduct() ;
969  }
970  }
971 
972  return total ;
973 
974  } else {
975  // Simply return integral
976  Double_t ret = integrate() / jacobianProduct() ;
977  return ret ;
978  }
979 }
980 
981 
982 ////////////////////////////////////////////////////////////////////////////////
983 /// Perform hybrid numerical/analytical integration over all real-valued dependents
984 
985 Double_t RooRealIntegral::integrate() const
986 {
987  if (!_numIntEngine) {
988  // Trivial case, fully analytical integration
989  return ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) ;
990  } else {
991  return _numIntEngine->calculate() ;
992  }
993 }
994 
995 
996 
997 ////////////////////////////////////////////////////////////////////////////////
998 /// Intercept server redirects and reconfigure internal object accordingly
999 
1000 Bool_t RooRealIntegral::redirectServersHook(const RooAbsCollection& /*newServerList*/,
1001  Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
1002 {
1003  _restartNumIntEngine = kTRUE ;
1004 
1005  autoSelectDirtyMode() ;
1006 
1007  // Update contents value caches for _intList and _sumList
1008  _saveInt.removeAll() ;
1009  _saveSum.removeAll() ;
1010  _intList.snapshot(_saveInt) ;
1011  _sumList.snapshot(_saveSum) ;
1012 
1013  // Delete parameters cache if we have one
1014  if (_params) {
1015  delete _params ;
1016  _params = 0 ;
1017  }
1018 
1019  return kFALSE ;
1020 }
1021 
1022 
1023 
1024 ////////////////////////////////////////////////////////////////////////////////
1025 
1026 const RooArgSet& RooRealIntegral::parameters() const
1027 {
1028  if (!_params) {
1029  _params = new RooArgSet("params") ;
1030 
1031  RooArgSet params ;
1032  for (const auto server : _serverList) {
1033  if (server->isValueServer(*this)) _params->add(*server) ;
1034  }
1035  }
1036 
1037  return *_params ;
1038 }
1039 
1040 
1041 
1042 ////////////////////////////////////////////////////////////////////////////////
1043 /// Dummy
1044 
1045 void RooRealIntegral::operModeHook()
1046 {
1047  if (_operMode==ADirty) {
1048 // cout << "RooRealIntegral::operModeHook(" << GetName() << " warning: mode set to ADirty" << endl ;
1049 // if (TString(GetName()).Contains("FULL")) {
1050 // cout << "blah" << endl ;
1051 // }
1052  }
1053 }
1054 
1055 
1056 
1057 ////////////////////////////////////////////////////////////////////////////////
1058 /// Check if current value is valid
1059 
1060 Bool_t RooRealIntegral::isValidReal(Double_t /*value*/, Bool_t /*printError*/) const
1061 {
1062  return kTRUE ;
1063 }
1064 
1065 ////////////////////////////////////////////////////////////////////////////////
1066 /// Check if component selection is allowed
1067 
1068 Bool_t RooRealIntegral::getAllowComponentSelection() const {
1069  return _respectCompSelect;
1070 }
1071 
1072 ////////////////////////////////////////////////////////////////////////////////
1073 /// Set component selection to be allowed/forbidden
1074 
1075 void RooRealIntegral::setAllowComponentSelection(Bool_t allow){
1076  _respectCompSelect = allow;
1077 }
1078 
1079 ////////////////////////////////////////////////////////////////////////////////
1080 /// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1081 /// integration operation
1082 
1083 void RooRealIntegral::printMetaArgs(ostream& os) const
1084 {
1085 
1086  if (intVars().getSize()!=0) {
1087  os << "Int " ;
1088  }
1089  os << _function.arg().GetName() ;
1090  if (_funcNormSet) {
1091  os << "_Norm" ;
1092  os << *_funcNormSet ;
1093  os << " " ;
1094  }
1095 
1096  // List internally integrated observables and factorizing observables as analytically integrated
1097  RooArgSet tmp(_anaList) ;
1098  tmp.add(_facList) ;
1099  if (tmp.getSize()>0) {
1100  os << "d[Ana]" ;
1101  os << tmp ;
1102  os << " " ;
1103  }
1104 
1105  // List numerically integrated and summed observables as numerically integrated
1106  RooArgSet tmp2(_intList) ;
1107  tmp2.add(_sumList) ;
1108  if (tmp2.getSize()>0) {
1109  os << " d[Num]" ;
1110  os << tmp2 ;
1111  os << " " ;
1112  }
1113 }
1114 
1115 
1116 
1117 ////////////////////////////////////////////////////////////////////////////////
1118 /// Print the state of this object to the specified output stream.
1119 
1120 void RooRealIntegral::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
1121 {
1122  RooAbsReal::printMultiline(os,contents,verbose,indent) ;
1123  os << indent << "--- RooRealIntegral ---" << endl;
1124  os << indent << " Integrates ";
1125  _function.arg().printStream(os,kName|kArgs,kSingleLine,indent);
1126  TString deeper(indent);
1127  deeper.Append(" ");
1128  os << indent << " operating mode is "
1129  << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << endl ;
1130  os << indent << " Summed discrete args are " << _sumList << endl ;
1131  os << indent << " Numerically integrated args are " << _intList << endl;
1132  os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << endl ;
1133  os << indent << " Arguments included in Jacobian are " << _jacList << endl ;
1134  os << indent << " Factorized arguments are " << _facList << endl ;
1135  os << indent << " Function normalization set " ;
1136  if (_funcNormSet)
1137  _funcNormSet->Print("1") ;
1138  else
1139  os << "<none>" ;
1140 
1141  os << endl ;
1142 }
1143 
1144 
1145 
1146 ////////////////////////////////////////////////////////////////////////////////
1147 /// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1148 
1149 void RooRealIntegral::setCacheAllNumeric(Int_t ndim) {
1150  _cacheAllNDim = ndim ;
1151 }
1152 
1153 
1154 ////////////////////////////////////////////////////////////////////////////////
1155 /// Return minimum dimensions of numeric integration for which values are cached.
1156 
1157 Int_t RooRealIntegral::getCacheAllNumeric()
1158 {
1159  return _cacheAllNDim ;
1160 }
1161 
1162