Logo ROOT   6.16/01
Reference Guide
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
22RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
23The class performs none of the actual integration, but only manages the logic
24of what variables can be integrated analytically, accounts for eventual jacobian
25terms and defines what numerical integrations needs to be done to complement the
26analytical integral.
27The actual analytical integrations (if any) are done in the PDF themselves, the numerical
28integration 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"
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
54using namespace std;
55
57;
58
59
61
62
63////////////////////////////////////////////////////////////////////////////////
64
66 _valid(kFALSE),
67 _respectCompSelect(kFALSE),
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{
83}
84
85
86
87////////////////////////////////////////////////////////////////////////////////
88/// Construct integral of 'function' over observables in 'depList'
89/// in range 'rangeName' with normalization observables 'funcNormSet'
90/// (for p.d.f.s). In the integral is performed to the maximum extent
91/// possible the internal (analytical) integrals advertised by function.
92/// The other integrations are performed numerically. The optional
93/// config object prescribes how these numeric integrations are configured.
94///
95
96RooRealIntegral::RooRealIntegral(const char *name, const char *title,
97 const RooAbsReal& function, const RooArgSet& depList,
98 const RooArgSet* funcNormSet, const RooNumIntConfig* config,
99 const char* rangeName) :
100 RooAbsReal(name,title),
101 _valid(kTRUE),
102 _respectCompSelect(kFALSE),
103 _sumList("!sumList","Categories to be summed numerically",this,kFALSE,kFALSE),
104 _intList("!intList","Variables to be integrated numerically",this,kFALSE,kFALSE),
105 _anaList("!anaList","Variables to be integrated analytically",this,kFALSE,kFALSE),
106 _jacList("!jacList","Jacobian product term",this,kFALSE,kFALSE),
107 _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
108 _facListIter(_facList.createIterator()),
109 _jacListIter(_jacList.createIterator()),
110 _function("!func","Function to be integrated",this,
111 const_cast<RooAbsReal&>(function),kFALSE,kFALSE),
112 _iconfig((RooNumIntConfig*)config),
113 _sumCat("!sumCat","SuperCategory for summation",this,kFALSE,kFALSE),
114 _sumCatIter(0),
115 _mode(0),
116 _intOperMode(Hybrid),
117 _restartNumIntEngine(kFALSE),
118 _numIntEngine(0),
119 _numIntegrand(0),
120 _rangeName((TNamed*)RooNameReg::ptr(rangeName)),
121 _params(0),
122 _cacheNum(kFALSE)
123{
124 // A) Check that all dependents are lvalues
125 //
126 // B) Check if list of dependents can be re-expressed in
127 // lvalues that are higher in the expression tree
128 //
129 // C) Check for dependents that the PDF insists on integrating
130 // analytically iself
131 //
132 // D) Make list of servers that can be integrated analytically
133 // Add all parameters/dependents as value/shape servers
134 //
135 // E) Interact with function to make list of objects actually integrated analytically
136 //
137 // F) Make list of numerical integration variables consisting of:
138 // - Category dependents of RealLValues in analytical integration
139 // - Leaf nodes server lists of function server that are not analytically integrated
140 // - Make Jacobian list for analytically integrated RealLValues
141 //
142 // G) Split numeric list in integration list and summation list
143 //
144
145 oocxcoutI(&function,Integration) << "RooRealIntegral::ctor(" << GetName() << ") Constructing integral of function "
146 << function.GetName() << " over observables" << depList << " with normalization "
147 << (funcNormSet?*funcNormSet:RooArgSet()) << " with range identifier "
148 << (rangeName?rangeName:"<none>") << endl ;
149
150
151 // Choose same expensive object cache as integrand
152 setExpensiveObjectCache(function.expensiveObjectCache()) ;
153// cout << "RRI::ctor(" << GetName() << ") setting expensive object cache to " << &expensiveObjectCache() << " as taken from " << function.GetName() << endl ;
154
155 // Use objects integrator configuration if none is specified
156 if (!_iconfig) _iconfig = (RooNumIntConfig*) function.getIntegratorConfig() ;
157
158 // Save private copy of funcNormSet, if supplied, excluding factorizing terms
159 if (funcNormSet) {
160 _funcNormSet = new RooArgSet ;
161 TIterator* iter = funcNormSet->createIterator() ;
162 RooAbsArg* nArg ;
163 while ((nArg=(RooAbsArg*)iter->Next())) {
164 if (function.dependsOn(*nArg)) {
165 _funcNormSet->addClone(*nArg) ;
166 }
167 }
168 delete iter ;
169 } else {
170 _funcNormSet = 0 ;
171 }
172
173 //_funcNormSet = funcNormSet ? (RooArgSet*)funcNormSet->snapshot(kFALSE) : 0 ;
174
175 // Make internal copy of dependent list
176 RooArgSet intDepList(depList) ;
177
178 RooAbsArg *arg ;
179
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
181 // * A) Check that all dependents are lvalues and filter out any
182 // dependents that the PDF doesn't explicitly depend on
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
184
185 TIterator* depIter = intDepList.createIterator() ;
186 while((arg=(RooAbsArg*)depIter->Next())) {
187 if(!arg->isLValue()) {
188 coutE(InputArguments) << ClassName() << "::" << GetName() << ": cannot integrate non-lvalue ";
189 arg->Print("1");
190 _valid= kFALSE;
191 }
192 if (!function.dependsOn(*arg)) {
193 RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
194 _facListOwned.addOwned(*argClone) ;
195 _facList.add(*argClone) ;
196 addServer(*argClone,kFALSE,kTRUE) ;
197 }
198 }
199
200 if (_facList.getSize()>0) {
201 oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables are " << _facList << endl ;
202 }
203
204
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
206 // * B) Check if list of dependents can be re-expressed in *
207 // * lvalues that are higher in the expression tree *
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
209
210
211 // Initial fill of list of LValue branches
212 RooArgSet exclLVBranches("exclLVBranches") ;
213 RooArgSet branchList,branchListVD ;
214 function.branchNodeServerList(&branchList) ;
215
216 TIterator* bIter = branchList.createIterator() ;
217 RooAbsArg* branch ;
218 while((branch=(RooAbsArg*)bIter->Next())) {
219 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(branch) ;
220 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(branch) ;
221 if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
222 exclLVBranches.add(*branch) ;
223// cout << "exclv branch = " << endl ;
224// branch->printCompactTree() ;
225 }
226 if (dependsOnValue(*branch)) {
227 branchListVD.add(*branch) ;
228 } else {
229// cout << "value of self does not depend on branch " << branch->GetName() << endl ;
230 }
231 }
232 delete bIter ;
233 exclLVBranches.remove(depList,kTRUE,kTRUE) ;
234// cout << "exclLVBranches = " << exclLVBranches << endl ;
235
236 // Initial fill of list of LValue leaf servers (put in intDepList)
237 RooArgSet exclLVServers("exclLVServers") ;
238 exclLVServers.add(intDepList) ;
239
240// cout << "begin exclLVServers = " << exclLVServers << endl ;
241
242 // Obtain mutual exclusive dependence by iterative reduction
243 TIterator *sIter = exclLVServers.createIterator() ;
244 bIter = exclLVBranches.createIterator() ;
245 RooAbsArg *server ;
246 Bool_t converged(kFALSE) ;
247 while(!converged) {
248 converged=kTRUE ;
249
250 // Reduce exclLVServers to only those serving exclusively exclLVBranches
251 sIter->Reset() ;
252 while ((server=(RooAbsArg*)sIter->Next())) {
253 if (!servesExclusively(server,exclLVBranches,branchListVD)) {
254 exclLVServers.remove(*server) ;
255// cout << "removing " << server->GetName() << " from exclLVServers because servesExclusively(" << server->GetName() << "," << exclLVBranches << ") faile" << endl ;
256 converged=kFALSE ;
257 }
258 }
259
260 // Reduce exclLVBranches to only those depending exclusisvely on exclLVservers
261 bIter->Reset() ;
262 while((branch=(RooAbsArg*)bIter->Next())) {
263 RooArgSet* brDepList = branch->getObservables(&intDepList) ;
264 RooArgSet bsList(*brDepList,"bsList") ;
265 delete brDepList ;
266 bsList.remove(exclLVServers,kTRUE,kTRUE) ;
267 if (bsList.getSize()>0) {
268 exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
269// cout << "removing " << branch->GetName() << " from exclLVBranches" << endl ;
270 converged=kFALSE ;
271 }
272 }
273 }
274
275 // Eliminate exclLVBranches that do not depend on any LVServer
276 bIter->Reset() ;
277 while((branch=(RooAbsArg*)bIter->Next())) {
278 if (!branch->dependsOnValue(exclLVServers)) {
279 //cout << "LV branch " << branch->GetName() << " does not depend on any LVServer (" << exclLVServers << ") and will be removed" << endl ;
280 exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
281 }
282 }
283
284 delete sIter ;
285 delete bIter ;
286
287// cout << "end exclLVServers = " << exclLVServers << endl ;
288
289 // Replace exclusive lvalue branch servers with lvalue branches
290 // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
291 if (exclLVServers.getSize()>0 && !function.isBinnedDistribution(exclLVBranches)) {
292// cout << "activating LVservers " << exclLVServers << " for use in integration " << endl ;
293 intDepList.remove(exclLVServers) ;
294 intDepList.add(exclLVBranches) ;
295
296 //cout << "intDepList removing exclLVServers " << exclLVServers << endl ;
297 //cout << "intDepList adding exclLVBranches " << exclLVBranches << endl ;
298
299 }
300
301
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
303 // * C) Check for dependents that the PDF insists on integrating *
304 // analytically iself *
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
306
307 RooArgSet anIntOKDepList ;
308 depIter->Reset() ;
309 while((arg=(RooAbsArg*)depIter->Next())) {
310 if (function.forceAnalyticalInt(*arg)) {
311 anIntOKDepList.add(*arg) ;
312 }
313 }
314 delete depIter ;
315
316 if (anIntOKDepList.getSize()>0) {
317 oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << endl ;
318 }
319
320
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
322 // * D) Make list of servers that can be integrated analytically *
323 // Add all parameters/dependents as value/shape servers *
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
325
326 sIter = function.serverIterator() ;
327 while((arg=(RooAbsArg*)sIter->Next())) {
328
329 //cout << "considering server" << arg->GetName() << endl ;
330
331 // Dependent or parameter?
332 if (!arg->dependsOnValue(intDepList)) {
333
334 //cout << " server does not depend on observables, adding server as value server to integral" << endl ;
335
336 if (function.dependsOnValue(*arg)) {
337 addServer(*arg,kTRUE,kFALSE) ;
338 }
339
340 continue ;
341
342 } else {
343
344 // Add final dependents of arg as shape servers
345 RooArgSet argLeafServers ;
346 arg->leafNodeServerList(&argLeafServers,0,kFALSE) ;
347
348 //arg->printCompactTree() ;
349 //cout << "leaf nodes of server are " << argLeafServers << " depList = " << depList << endl ;
350
351 // Skip arg if it is neither value or shape server
352 if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
353 //cout << " server is neither value not shape server of function, ignoring" << endl ;
354 continue ;
355 }
356
357 TIterator* lIter = argLeafServers.createIterator() ;
358 RooAbsArg* leaf ;
359 while((leaf=(RooAbsArg*)lIter->Next())) {
360
361 //cout << " considering leafnode " << leaf->GetName() << " of server " << arg->GetName() << endl ;
362
363 if (depList.find(leaf->GetName()) && function.dependsOnValue(*leaf)) {
364
365 RooAbsRealLValue* leaflv = dynamic_cast<RooAbsRealLValue*>(leaf) ;
366 if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
367 oocxcoutD(&function,Integration) << function.GetName() << " : Observable " << leaf->GetName() << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf" << endl ;
368 if (leaflv->getBinning(rangeName).lowBoundFunc()) {
369 addServer(*leaflv->getBinning(rangeName).lowBoundFunc(),kTRUE,kFALSE) ;
370 }
371 if(leaflv->getBinning(rangeName).highBoundFunc()) {
372 addServer(*leaflv->getBinning(rangeName).highBoundFunc(),kTRUE,kFALSE) ;
373 }
374 } else {
375 oocxcoutD(&function,Integration) << function.GetName() << ": Adding observable " << leaf->GetName() << " of server "
376 << arg->GetName() << " as shape dependent" << endl ;
377 addServer(*leaf,kFALSE,kTRUE) ;
378 }
379 } else if (!depList.find(leaf->GetName())) {
380
381 if (function.dependsOnValue(*leaf)) {
382 oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as value dependent" << endl ;
383 addServer(*leaf,kTRUE,kFALSE) ;
384 } else {
385 oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as shape dependent" << endl ;
386 addServer(*leaf,kFALSE,kTRUE) ;
387 }
388 }
389 }
390 delete lIter ;
391 }
392
393 // If this dependent arg is self-normalized, stop here
394 //if (function.selfNormalized()) continue ;
395
396 Bool_t depOK(kFALSE) ;
397 // Check for integratable AbsRealLValue
398
399 //cout << "checking server " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
400
401 if (arg->isDerived()) {
402 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(arg) ;
403 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(arg) ;
404 //cout << "realArgLV = " << realArgLV << " intDepList = " << intDepList << endl ;
405 if ((realArgLV && intDepList.find(realArgLV->GetName()) && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
406
407 //cout << " arg " << arg->GetName() << " is derived LValue with valid jacobian" << endl ;
408
409 // Derived LValue with valid jacobian
410 depOK = kTRUE ;
411
412 // Now, check for overlaps
413 Bool_t overlapOK = kTRUE ;
414 RooAbsArg *otherArg ;
415 TIterator* sIter2 = function.serverIterator() ;
416 while((otherArg=(RooAbsArg*)sIter2->Next())) {
417 // skip comparison with self
418 if (arg==otherArg) continue ;
419 if (otherArg->IsA()==RooConstVar::Class()) continue ;
420 if (arg->overlaps(*otherArg,kTRUE)) {
421 //cout << "arg " << arg->GetName() << " overlaps with " << otherArg->GetName() << endl ;
422 //overlapOK=kFALSE ;
423 }
424 }
425 // coverity[DEADCODE]
426 if (!overlapOK) depOK=kFALSE ;
427
428 //cout << "overlap check returns OK=" << (depOK?"T":"F") << endl ;
429
430 delete sIter2 ;
431 }
432 } else {
433 // Fundamental types are always OK
434 depOK = kTRUE ;
435 }
436
437 // Add server to list of dependents that are OK for analytical integration
438 if (depOK) {
439 anIntOKDepList.add(*arg,kTRUE) ;
440 oocxcoutI(&function,Integration) << function.GetName() << ": Observable " << arg->GetName() << " is suitable for analytical integration (if supported by p.d.f)" << endl ;
441 }
442 }
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
444 // * E) interact with function to make list of objects actually integrated analytically *
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
446
447 RooArgSet anIntDepList ;
448
449 RooArgSet *anaSet = new RooArgSet( _anaList, Form("UniqueCloneOf_%s",_anaList.GetName()));
452 _anaList.add(*anaSet);
453 delete anaSet;
454
455 // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaListx
456 if (_mode==0) {
458 }
459
460 if (_mode!=0) {
461 oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << endl ;
462 }
463
464
465 // WVE kludge: synchronize dset for use in analyticalIntegral
466 function.getVal(_funcNormSet) ;
467
468 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
469 // * F) Make list of numerical integration variables consisting of: *
470 // * - Category dependents of RealLValues in analytical integration *
471 // * - Expanded server lists of server that are not analytically integrated *
472 // * Make Jacobian list with analytically integrated RealLValues *
473 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
474
475 RooArgSet numIntDepList ;
476
477
478 // Loop over actually analytically integrated dependents
479 TIterator* aiIter = _anaList.createIterator() ;
480 while ((arg=(RooAbsArg*)aiIter->Next())) {
481
482 // Process only derived RealLValues
483 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class()) && arg->isDerived() && !arg->isFundamental()) {
484
485 // Add to list of Jacobians to calculate
486 _jacList.add(*arg) ;
487
488 // Add category dependent of LValueReal used in integration
489 RooAbsArg *argDep ;
490 RooArgSet *argDepList = arg->getObservables(&intDepList) ;
491 TIterator *adIter = argDepList->createIterator() ;
492 while ((argDep=(RooAbsArg*)adIter->Next())) {
493 if (argDep->IsA()->InheritsFrom(RooAbsCategoryLValue::Class()) && intDepList.contains(*argDep)) {
494 numIntDepList.add(*argDep,kTRUE) ;
495 }
496 }
497 delete adIter ;
498 delete argDepList ;
499 }
500 }
501 delete aiIter ;
502
503
504 // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
505 if (_anaList.getSize()==0) {
506 if (exclLVServers.getSize()>0) {
507 //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << endl ;
508 intDepList.remove(exclLVBranches) ;
509 intDepList.add(exclLVServers) ;
510 }
511 }
512 //cout << "NUMINT intDepList = " << intDepList << endl ;
513
514 // Loop again over function servers to add remaining numeric integrations
515 sIter->Reset() ;
516 while((arg=(RooAbsArg*)sIter->Next())) {
517
518 //cout << "processing server for numeric integration " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
519
520 // Process only servers that are not treated analytically
521 if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
522
523 // Process only derived RealLValues
524 if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
525 numIntDepList.add(*arg,kTRUE) ;
526 } else {
527
528 // WVE this will only get the observables, but not l-value transformations
529 // Expand server in final dependents
530 RooArgSet *argDeps = arg->getObservables(&intDepList) ;
531
532 if (argDeps->getSize()>0) {
533
534 // Add final dependents, that are not forcibly integrated analytically,
535 // to numerical integration list
536 TIterator* iter = argDeps->createIterator() ;
537 RooAbsArg* dep ;
538 while((dep=(RooAbsArg*)iter->Next())) {
539 if (!_anaList.find(dep->GetName())) {
540 numIntDepList.add(*dep,kTRUE) ;
541 }
542 }
543 delete iter ;
544
545 }
546 delete argDeps ;
547 }
548
549 }
550 }
551 delete sIter ;
552
553 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
554 // * G) Split numeric list in integration list and summation list *
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
556
557 // Split numeric integration list in summation and integration lists
558 TIterator* numIter=numIntDepList.createIterator() ;
559 while ((arg=(RooAbsArg*)numIter->Next())) {
560
561 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
562 _intList.add(*arg) ;
563 } else if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
564 _sumList.add(*arg) ;
565 }
566 }
567 delete numIter ;
568
569 if (_anaList.getSize()>0) {
570 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << endl ;
571 }
572 if (_intList.getSize()>0) {
573 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << endl ;
574 }
575 if (_sumList.getSize()>0) {
576 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << endl ;
577 }
578
579
580 // Determine operating mode
581 if (numIntDepList.getSize()>0) {
582 // Numerical and optional Analytical integration
584 } else if (_anaList.getSize()>0) {
585 // Purely analytical integration
587 } else {
588 // No integration performed
590 }
591
592 // Determine auto-dirty status
594
595 // Create value caches for _intList and _sumList
598
599
600 if (_sumList.getSize()>0) {
601 RooSuperCategory *sumCat = new RooSuperCategory(Form("%s_sumCat",GetName()),"sumCat",_sumList) ;
602 _sumCatIter = sumCat->typeIterator() ;
603 _sumCat.addOwned(*sumCat) ;
604 }
605
607}
608
609
610
611////////////////////////////////////////////////////////////////////////////////
612/// Set appropriate cache operation mode for integral depending on cache operation
613/// mode of server objects
614
616{
617 // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
618 TIterator* siter = serverIterator() ;
619 RooAbsArg* server ;
620 while((server=(RooAbsArg*)siter->Next())){
621 if (server->isValueServer(*this)) {
622 RooArgSet leafSet ;
623 server->leafNodeServerList(&leafSet) ;
624 TIterator* liter = leafSet.createIterator() ;
625 RooAbsArg* leaf ;
626 while((leaf=(RooAbsArg*)liter->Next())) {
627 if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
629 break ;
630 }
631 if (leaf->getAttribute("projectedDependent")) {
633 break ;
634 }
635 }
636 delete liter ;
637 }
638 }
639 delete siter ;
640}
641
642
643
644////////////////////////////////////////////////////////////////////////////////
645/// Utility function that returns true if 'object server' is a server
646/// to exactly one of the RooAbsArgs in 'exclLVBranches'
647
648Bool_t RooRealIntegral::servesExclusively(const RooAbsArg* server,const RooArgSet& exclLVBranches, const RooArgSet& allBranches) const
649{
650 // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
651
652 // Special case, no LV servers available
653 if (exclLVBranches.getSize()==0) return kFALSE ;
654
655 // If server has no clients and is not an LValue itself, return false
656 if (server->_clientList.GetSize()==0 && exclLVBranches.find(server->GetName())) {
657 return kFALSE ;
658 }
659
660 // WVE must check for value relations only here!!!!
661
662
663// cout << "servesExclusively: does " << server->GetName() << " serve only one of " << exclLVBranches << endl ;
664
665 // Loop over all clients
666 Int_t numLVServ(0) ;
667 RooAbsArg* client ;
668 TIterator* cIter = server->valueClientIterator() ;
669 while((client=(RooAbsArg*)cIter->Next())) {
670// cout << "now checking value client " << client->GetName() << " of server " << server->GetName() << endl ;
671 // If client is not an LValue, recurse
672 if (!(exclLVBranches.find(client->GetName())==client)) {
673// cout << " client " << client->GetName() << "is not an lvalue" << endl ;
674 if (allBranches.find(client->GetName())==client) {
675// cout << " ... recursing call" << endl ;
676 if (!servesExclusively(client,exclLVBranches,allBranches)) {
677 // Client is a non-LValue that doesn't have an exclusive LValue server
678 delete cIter ;
679// cout << "client " << client->GetName() << " is a non-lvalue that doesn't have an exclusive lvalue server" << endl ;
680 return kFALSE ;
681 }
682 }
683 } else {
684 // Client is an LValue
685// cout << "client " << client->GetName() << " of server " << server->GetName() << " is an LValue " << endl ;
686 numLVServ++ ;
687 }
688 }
689
690 delete cIter ;
691// cout << "numLVserv = " << numLVServ << endl ;
692 return (numLVServ==1) ;
693}
694
695
696
697
698////////////////////////////////////////////////////////////////////////////////
699/// (Re)Initialize numerical integration engine if necessary. Return kTRUE if
700/// successful, or otherwise kFALSE.
701
703{
704 // if we already have an engine, check if it still works for the present limits.
705 if(0 != _numIntEngine) {
707 // otherwise, cleanup the old engine
708 delete _numIntEngine ;
709 _numIntEngine= 0;
710 if(0 != _numIntegrand) {
711 delete _numIntegrand;
712 _numIntegrand= 0;
713 }
714 }
715
716 // All done if there are no arguments to integrate numerically
717 if(0 == _intList.getSize()) return kTRUE;
718
719 // Bind the appropriate analytic integral (specified by _mode) of our RooRealVar object to
720 // those of its arguments that will be integrated out numerically.
721 if(_mode != 0) {
723 }
724 else {
726 }
727 if(0 == _numIntegrand || !_numIntegrand->isValid()) {
728 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << endl;
729 return kFALSE;
730 }
731
732 // Create appropriate numeric integrator using factory
735
736 if(0 == _numIntEngine || !_numIntEngine->isValid()) {
737 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << endl;
738 return kFALSE;
739 }
740
741 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
742 << _numIntEngine->IsA()->GetName() << " to calculate Int" << _intList << endl ;
743
744 if (_intList.getSize()>3) {
745 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 ;
746 }
747
749 return kTRUE;
750}
751
752
753
754////////////////////////////////////////////////////////////////////////////////
755/// Copy constructor
756
758 RooAbsReal(other,name),
759 _valid(other._valid),
760 _respectCompSelect(other._respectCompSelect),
761 _sumList("!sumList",this,other._sumList),
762 _intList("!intList",this,other._intList),
763 _anaList("!anaList",this,other._anaList),
764 _jacList("!jacList",this,other._jacList),
765 _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
766 _facListIter(_facList.createIterator()),
767 _jacListIter(_jacList.createIterator()),
768 _function("!func",this,other._function),
769 _iconfig(other._iconfig),
770 _sumCat("!sumCat",this,other._sumCat),
771 _sumCatIter(0),
772 _mode(other._mode),
773 _intOperMode(other._intOperMode),
774 _restartNumIntEngine(kFALSE),
775 _numIntEngine(0),
776 _numIntegrand(0),
777 _rangeName(other._rangeName),
778 _params(0),
779 _cacheNum(kFALSE)
780{
782
783 other._facListIter->Reset() ;
784 RooAbsArg* arg ;
785 while((arg=(RooAbsArg*)other._facListIter->Next())) {
786 RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
787 _facListOwned.addOwned(*argClone) ;
788 _facList.add(*argClone) ;
789 addServer(*argClone,kFALSE,kTRUE) ;
790 }
791
792 other._intList.snapshot(_saveInt) ;
793 other._sumList.snapshot(_saveSum) ;
794
796}
797
798
799
800////////////////////////////////////////////////////////////////////////////////
801
803 // Destructor
804{
805 if (_numIntEngine) delete _numIntEngine ;
806 if (_numIntegrand) delete _numIntegrand ;
807 if (_funcNormSet) delete _funcNormSet ;
808 delete _facListIter ;
809 delete _jacListIter ;
810 if (_sumCatIter) delete _sumCatIter ;
811 if (_params) delete _params ;
812
814}
815
816
817
818
819
820////////////////////////////////////////////////////////////////////////////////
821
822RooAbsReal* RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const
823{
824 // Handle special case of no integration with default algorithm
825 if (iset.getSize()==0) {
826 return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
827 }
828
829 // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
830 RooArgSet isetAll(iset) ;
831 isetAll.add(_sumList) ;
832 isetAll.add(_intList) ;
833 isetAll.add(_anaList) ;
834 isetAll.add(_facList) ;
835
836 const RooArgSet* newNormSet(0) ;
837 RooArgSet* tmp(0) ;
838 if (nset && !_funcNormSet) {
839 newNormSet = nset ;
840 } else if (!nset && _funcNormSet) {
841 newNormSet = _funcNormSet ;
842 } else if (nset && _funcNormSet) {
843 tmp = new RooArgSet ;
844 tmp->add(*nset) ;
845 tmp->add(*_funcNormSet,kTRUE) ;
846 newNormSet = tmp ;
847 }
848 RooAbsReal* ret = _function.arg().createIntegral(isetAll,newNormSet,cfg,rangeName) ;
849
850 if (tmp) {
851 delete tmp ;
852 }
853
854 return ret ;
855}
856
857
858
859////////////////////////////////////////////////////////////////////////////////
860/// Return value of object. If the cache is clean, return the
861/// cached value, otherwise recalculate on the fly and refill
862/// the cache
863
865{
866// // fast-track clean-cache processing
867// if (_operMode==AClean) {
868// return _value ;
869// }
870
871 if (nset && nset!=_lastNSet) {
872 ((RooAbsReal*) this)->setProxyNormSet(nset) ;
873 _lastNSet = (RooArgSet*) nset ;
874 }
875
877 _value = traceEval(nset) ;
878 }
879
880 return _value ;
881}
882
883
884
885
886
887////////////////////////////////////////////////////////////////////////////////
888/// Perform the integration and return the result
889
891{
892
896 }
897
898 Double_t retVal(0) ;
899 switch (_intOperMode) {
900
901 case Hybrid:
902 {
903 // Cache numeric integrals in >1d expensive object cache
904 RooDouble* cacheVal(0) ;
907 }
908
909 if (cacheVal) {
910 retVal = *cacheVal ;
911 // cout << "using cached value of integral" << GetName() << endl ;
912 } else {
913
914
915 // Find any function dependents that are AClean
916 // and switch them temporarily to ADirty
917 Bool_t origState = inhibitDirty() ;
919
920 // try to initialize our numerical integration engine
921 if(!(_valid= initNumIntegrator())) {
922 coutE(Integration) << ClassName() << "::" << GetName()
923 << ":evaluate: cannot initialize numerical integrator" << endl;
925 return 0;
926 }
927
928 // Save current integral dependent values
931
932 // Evaluate sum/integral
933 retVal = sum() ;
934
935
936 // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
937 setDirtyInhibit(origState) ;
938
939 // Restore integral dependent values
942
943 // Cache numeric integrals in >1d expensive object cache
945 RooDouble* val = new RooDouble(retVal) ;
947 // cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << endl ;
948 }
949
950 }
951 break ;
952 }
953 case Analytic:
954 {
956 cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
957 << ")func = " << _function.arg().IsA()->GetName() << "::" << _function.arg().GetName()
958 << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << endl ;
959
960
961 break ;
962 }
963
964 case PassThrough:
965 {
966 //setDirtyInhibit(kTRUE) ;
967 retVal= _function.arg().getVal(_funcNormSet) ;
968 //setDirtyInhibit(kFALSE) ;
969 break ;
970 }
971 }
972
973
974 // Multiply answer with integration ranges of factorized variables
975 if (_facList.getSize()>0) {
976 RooAbsArg *arg ;
978 while((arg=(RooAbsArg*)_facListIter->Next())) {
979 // Multiply by fit range for 'real' dependents
980 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
981 RooAbsRealLValue* argLV = (RooAbsRealLValue*)arg ;
982 retVal *= (argLV->getMax() - argLV->getMin()) ;
983 }
984 // Multiply by number of states for category dependents
985 if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
987 retVal *= argLV->numTypes() ;
988 }
989 }
990 }
991
992
993 if (dologD(Tracing)) {
994 cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
995 switch(_mode) {
996 case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
997 case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
998 case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
999 }
1000
1001 ccxcoutD(Tracing) << "raw*fact = " << retVal << endl ;
1002 }
1003
1005
1006 return retVal ;
1007}
1008
1009
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Return product of jacobian terms originating from analytical integration
1013
1015{
1016 if (_jacList.getSize()==0) {
1017 return 1 ;
1018 }
1019
1020 Double_t jacProd(1) ;
1021 _jacListIter->Reset() ;
1022 RooAbsRealLValue* arg ;
1023 while ((arg=(RooAbsRealLValue*)_jacListIter->Next())) {
1024 jacProd *= arg->jacobian() ;
1025 }
1026
1027 // Take fabs() here: if jacobian is negative, min and max are swapped and analytical integral
1028 // will be positive, so must multiply with positive jacobian.
1029 return fabs(jacProd) ;
1030}
1031
1032
1033
1034////////////////////////////////////////////////////////////////////////////////
1035/// Perform summation of list of category dependents to be integrated
1036
1038{
1039 if (_sumList.getSize()!=0) {
1040 // Add integrals for all permutations of categories summed over
1041 Double_t total(0) ;
1042
1043 _sumCatIter->Reset() ;
1044 RooCatType* type ;
1046 while((type=(RooCatType*)_sumCatIter->Next())) {
1047 sumCat->setIndex(type->getVal()) ;
1048 if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
1050 }
1051 }
1052
1053 return total ;
1054
1055 } else {
1056 // Simply return integral
1057 Double_t ret = integrate() / jacobianProduct() ;
1058 return ret ;
1059 }
1060}
1061
1062
1063////////////////////////////////////////////////////////////////////////////////
1064/// Perform hybrid numerical/analytical integration over all real-valued dependents
1065
1067{
1068 if (!_numIntEngine) {
1069 // Trivial case, fully analytical integration
1070 return ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) ;
1071 } else {
1072 return _numIntEngine->calculate() ;
1073 }
1074}
1075
1076
1077
1078////////////////////////////////////////////////////////////////////////////////
1079/// Intercept server redirects and reconfigure internal object accordingly
1080
1082 Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
1083{
1085
1087
1088 // Update contents value caches for _intList and _sumList
1093
1094 // Delete parameters cache if we have one
1095 if (_params) {
1096 delete _params ;
1097 _params = 0 ;
1098 }
1099
1100 return kFALSE ;
1101}
1102
1103
1104
1105////////////////////////////////////////////////////////////////////////////////
1106
1108{
1109 if (!_params) {
1110 _params = new RooArgSet("params") ;
1111
1112 TIterator* siter = serverIterator() ;
1113 RooArgSet params ;
1114 RooAbsArg* server ;
1115 while((server = (RooAbsArg*)siter->Next())) {
1116 if (server->isValueServer(*this)) _params->add(*server) ;
1117 }
1118 delete siter ;
1119 }
1120
1121 return *_params ;
1122}
1123
1124
1125
1126////////////////////////////////////////////////////////////////////////////////
1127/// Dummy
1128
1130{
1131 if (_operMode==ADirty) {
1132// cout << "RooRealIntegral::operModeHook(" << GetName() << " warning: mode set to ADirty" << endl ;
1133// if (TString(GetName()).Contains("FULL")) {
1134// cout << "blah" << endl ;
1135// }
1136 }
1137}
1138
1139
1140
1141////////////////////////////////////////////////////////////////////////////////
1142/// Check if current value is valid
1143
1144Bool_t RooRealIntegral::isValidReal(Double_t /*value*/, Bool_t /*printError*/) const
1145{
1146 return kTRUE ;
1147}
1148
1149////////////////////////////////////////////////////////////////////////////////
1150/// Check if component selection is allowed
1151
1153 return _respectCompSelect;
1154}
1155
1156////////////////////////////////////////////////////////////////////////////////
1157/// Set component selection to be allowed/forbidden
1158
1160 _respectCompSelect = allow;
1161}
1162
1163////////////////////////////////////////////////////////////////////////////////
1164/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1165/// integration operation
1166
1167void RooRealIntegral::printMetaArgs(ostream& os) const
1168{
1169
1170 if (intVars().getSize()!=0) {
1171 os << "Int " ;
1172 }
1173 os << _function.arg().GetName() ;
1174 if (_funcNormSet) {
1175 os << "_Norm" ;
1176 os << *_funcNormSet ;
1177 os << " " ;
1178 }
1179
1180 // List internally integrated observables and factorizing observables as analytically integrated
1181 RooArgSet tmp(_anaList) ;
1182 tmp.add(_facList) ;
1183 if (tmp.getSize()>0) {
1184 os << "d[Ana]" ;
1185 os << tmp ;
1186 os << " " ;
1187 }
1188
1189 // List numerically integrated and summed observables as numerically integrated
1190 RooArgSet tmp2(_intList) ;
1191 tmp2.add(_sumList) ;
1192 if (tmp2.getSize()>0) {
1193 os << " d[Num]" ;
1194 os << tmp2 ;
1195 os << " " ;
1196 }
1197}
1198
1199
1200
1201////////////////////////////////////////////////////////////////////////////////
1202/// Print the state of this object to the specified output stream.
1203
1204void RooRealIntegral::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
1205{
1206 RooAbsReal::printMultiline(os,contents,verbose,indent) ;
1207 os << indent << "--- RooRealIntegral ---" << endl;
1208 os << indent << " Integrates ";
1210 TString deeper(indent);
1211 deeper.Append(" ");
1212 os << indent << " operating mode is "
1213 << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << endl ;
1214 os << indent << " Summed discrete args are " << _sumList << endl ;
1215 os << indent << " Numerically integrated args are " << _intList << endl;
1216 os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << endl ;
1217 os << indent << " Arguments included in Jacobian are " << _jacList << endl ;
1218 os << indent << " Factorized arguments are " << _facList << endl ;
1219 os << indent << " Function normalization set " ;
1220 if (_funcNormSet)
1221 _funcNormSet->Print("1") ;
1222 else
1223 os << "<none>" ;
1224
1225 os << endl ;
1226}
1227
1228
1229
1230////////////////////////////////////////////////////////////////////////////////
1231/// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1232
1234 _cacheAllNDim = ndim ;
1235}
1236
1237
1238////////////////////////////////////////////////////////////////////////////////
1239/// Return minimum dimensions of numeric integration for which values are cached.
1240
1242{
1243 return _cacheAllNDim ;
1244}
1245
1246
void Class()
Definition: Class.C:29
#define cxcoutI(a)
Definition: RooMsgService.h:83
#define cxcoutD(a)
Definition: RooMsgService.h:79
#define oocxcoutD(o, a)
Definition: RooMsgService.h:81
#define dologD(a)
Definition: RooMsgService.h:63
#define coutE(a)
Definition: RooMsgService.h:34
#define ccxcoutD(a)
Definition: RooMsgService.h:80
#define ccoutD(a)
Definition: RooMsgService.h:37
#define oocxcoutI(o, a)
Definition: RooMsgService.h:85
#define TRACE_DESTROY
Definition: RooTrace.h:23
#define TRACE_CREATE
Definition: RooTrace.h:22
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
static unsigned int total
int type
Definition: TGX11.cxx:120
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
RooExpensiveObjectCache & expensiveObjectCache() const
Definition: RooAbsArg.cxx:2334
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
virtual Bool_t isLValue() const
Definition: RooAbsArg.h:170
void addServer(RooAbsArg &server, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
Definition: RooAbsArg.cxx:335
Bool_t isValueServer(const RooAbsArg &arg) const
Definition: RooAbsArg.h:134
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
Definition: RooAbsArg.cxx:736
Bool_t isShapeServer(const RooAbsArg &arg) const
Definition: RooAbsArg.h:142
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
Definition: RooAbsArg.cxx:466
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:75
static void setDirtyInhibit(Bool_t flag)
Control global dirty inhibit mode.
Definition: RooAbsArg.cxx:214
friend class RooArgSet
Definition: RooAbsArg.h:471
Bool_t inhibitDirty() const
Delete watch flag.
Definition: RooAbsArg.cxx:84
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition: RooAbsArg.h:499
virtual Bool_t isFundamental() const
Definition: RooAbsArg.h:157
Bool_t isValueOrShapeDirtyAndClear() const
Definition: RooAbsArg.h:370
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
Bool_t overlaps(const RooAbsArg &testArg, Bool_t valueOnly=kFALSE) const
Test if any of the nodes of tree are shared with that of the given tree.
Definition: RooAbsArg.cxx:795
virtual Bool_t isDerived() const
Definition: RooAbsArg.h:81
TIterator * valueClientIterator() const
Definition: RooAbsArg.h:104
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Change cache operation mode to given mode.
Definition: RooAbsArg.cxx:1746
OperMode _operMode
Definition: RooAbsArg.h:575
Bool_t getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
Definition: RooAbsArg.cxx:264
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Definition: RooAbsArg.h:88
TIterator * serverIterator() const
Definition: RooAbsArg.h:112
OperMode operMode() const
Definition: RooAbsArg.h:399
RooRefCountList _clientList
Definition: RooAbsArg.h:476
virtual Bool_t isParameterized() const
Definition: RooAbsBinning.h:79
virtual RooAbsReal * highBoundFunc() const
Definition: RooAbsBinning.h:87
virtual RooAbsReal * lowBoundFunc() const
Definition: RooAbsBinning.h:83
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
TIterator * typeIterator() const
Return iterator over all defined states.
Int_t numTypes(const char *=0) const
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooAbsArg * first() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
const char * GetName() const
Returns name of object.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Bool_t isValid() const
Definition: RooAbsFunc.h:33
Double_t calculate(const Double_t *yvec=0)
Calculate integral value with given array of parameter values.
virtual Bool_t checkLimits() const
Bool_t isValid() const
Abstract base class for objects that are lvalues, i.e.
Definition: RooAbsLValue.h:26
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual Double_t getMax(const char *name=0) const
virtual const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const =0
virtual Bool_t isJacobianOK(const RooArgSet &depList) const
virtual Double_t jacobian() const
virtual Double_t getMin(const char *name=0) const
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
virtual Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
Definition: RooAbsReal.cxx:319
friend class RooRealBinding
Definition: RooAbsReal.h:387
RooArgSet * _lastNSet
Definition: RooAbsReal.h:484
Double_t traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
Definition: RooAbsReal.cxx:287
virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition: RooAbsReal.cxx:348
Double_t _value
Definition: RooAbsReal.h:391
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Definition: RooAbsReal.cxx:423
Double_t getVal(const RooArgSet *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
static Bool_t _globalSelectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition: RooAbsReal.h:482
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:502
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:278
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t addOwned(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling addOwned() for each element in the source...
Definition: RooArgSet.h:92
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
virtual void addClone(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling addOwned() for each element in the source...
Definition: RooArgSet.h:96
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state.
Definition: RooCatType.h:22
RooDouble is a minimal implementation of a TObject holding a Double_t value.
Definition: RooDouble.h:22
Bool_t registerObject(const char *ownerName, const char *objectName, TObject &cacheObject, TIterator *paramIter)
Register object associated with given name and given associated parameters with given values in cache...
const TObject * retrieveObject(const char *name, TClass *tclass, const RooArgSet &params)
Retrieve object from cache that was registered under given name with given parameters,...
Int_t GetSize() const
Definition: RooLinkedList.h:60
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::addOwned()
RooNameReg is a registry for 'const char*' name.
Definition: RooNameReg.h:23
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition: RooNameReg.cxx:135
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
RooAbsIntegrator * createIntegrator(RooAbsFunc &func, const RooNumIntConfig &config, Int_t ndim=0, Bool_t isBinned=kFALSE)
Construct a numeric integrator instance that operates on function 'func' and is configured with 'conf...
static RooNumIntFactory & instance()
Static method returning reference to singleton instance of factory.
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
Lightweight RooAbsFunc interface adaptor that binds an analytic integral of a RooAbsReal object (spec...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooNumIntConfig * _iconfig
do not persist
virtual void operModeHook()
Cache all integrals with given numeric dimension.
RooArgSet * _params
Bool_t getAllowComponentSelection() const
Check if component selection is allowed.
Bool_t servesExclusively(const RooAbsArg *server, const RooArgSet &exclLVBranches, const RooArgSet &allBranches) const
Utility function that returns true if 'object server' is a server to exactly one of the RooAbsArgs in...
RooRealProxy _function
Iterator over lvalue observables with Jacobian.
RooArgSet intVars() const
void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print the state of this object to the specified output stream.
RooSetProxy _intList
TIterator * _facListIter
TIterator * _sumCatIter
do not persist
RooSetProxy _facList
RooAbsFunc * _numIntegrand
do not persist
static void setCacheAllNumeric(Int_t ndim)
Global switch to cache all integral values that integrate at least ndim dimensions numerically.
IntOperMode _intOperMode
RooArgSet * _funcNormSet
RooAbsIntegrator * _numIntEngine
do not persist
virtual ~RooRealIntegral()
void printMetaArgs(std::ostream &os) const
Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the...
const RooArgSet & parameters() const
virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const
Check if current value is valid.
TNamed * _rangeName
do not persist
TIterator * _jacListIter
Iterator over factorizing observables.
virtual RooAbsReal * createIntegral(const RooArgSet &iset, const RooArgSet *nset=0, const RooNumIntConfig *cfg=0, const char *rangeName=0) const
Create an object that represents the integral of the function over one or more observables listed in ...
RooSetProxy _jacList
virtual Double_t sum() const
Perform summation of list of category dependents to be integrated.
RooSetProxy _anaList
Bool_t _restartNumIntEngine
RooSetProxy _sumList
Bool_t initNumIntegrator() const
(Re)Initialize numerical integration engine if necessary.
virtual Bool_t redirectServersHook(const RooAbsCollection &newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
Intercept server redirects and reconfigure internal object accordingly.
virtual Double_t getValV(const RooArgSet *set=0) const
Return value of object.
virtual Double_t jacobianProduct() const
Return product of jacobian terms originating from analytical integration.
RooListProxy _sumCat
Bool_t _cacheNum
cache for set of parameters
static Int_t getCacheAllNumeric()
Return minimum dimensions of numeric integration for which values are cached.
RooArgSet _saveSum
do not persist
static Int_t _cacheAllNDim
void setAllowComponentSelection(Bool_t allow)
Set component selection to be allowed/forbidden.
Double_t evaluate() const
Perform the integration and return the result.
virtual Double_t integrate() const
Perform hybrid numerical/analytical integration over all real-valued dependents.
void autoSelectDirtyMode()
Set appropriate cache operation mode for integral depending on cache operation mode of server objects...
Bool_t _respectCompSelect
RooArgSet _facListOwned
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual void removeAll()
Remove all argument inset using remove(const RooAbsArg&).
RooSuperCategory consolidates several RooAbsCategoryLValue objects into a single category.
virtual Bool_t inRange(const char *rangeName) const
Return true of all of the input category states are in the given range.
virtual Bool_t setIndex(Int_t index, Bool_t printError=kTRUE)
Set the value of the super category by specifying the state index code by setting the states of the c...
Iterator abstract base class.
Definition: TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
@ NumIntegration
Definition: RooGlobalFunc.h:59
@ InputArguments
Definition: RooGlobalFunc.h:58
@ Integration
Definition: RooGlobalFunc.h:57
static constexpr double liter
STL namespace.