Logo ROOT  
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(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{
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
94RooRealIntegral::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
222 // Initial fill of list of LValue leaf servers (put in intDepList)
223 RooArgSet exclLVServers("exclLVServers") ;
224 exclLVServers.add(intDepList) ;
225
226 // Obtain mutual exclusive dependence by iterative reduction
227 Bool_t converged(kFALSE) ;
228 while(!converged) {
229 converged=kTRUE ;
230
231 // Reduce exclLVServers to only those serving exclusively exclLVBranches
232 std::vector<RooAbsArg*> toBeRemoved;
233 for (auto server : exclLVServers) {
234 if (!servesExclusively(server,exclLVBranches,branchListVD)) {
235 toBeRemoved.push_back(server);
236 converged=kFALSE ;
237 }
238 }
239 exclLVServers.remove(toBeRemoved.begin(), toBeRemoved.end());
240
241 // Reduce exclLVBranches to only those depending exclusively on exclLVservers
242 // Attention: counting loop, since erasing from container
243 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
244 const RooAbsArg* branch = exclLVBranches[i];
245 RooArgSet* brDepList = branch->getObservables(&intDepList) ;
246 RooArgSet bsList(*brDepList,"bsList") ;
247 delete brDepList ;
248 bsList.remove(exclLVServers,kTRUE,kTRUE) ;
249 if (bsList.getSize()>0) {
250 exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
251 --i;
252 converged=kFALSE ;
253 }
254 }
255 }
256
257 // Eliminate exclLVBranches that do not depend on any LVServer
258 // Attention: Counting loop, since modifying container
259 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
260 const RooAbsArg* branch = exclLVBranches[i];
261 if (!branch->dependsOnValue(exclLVServers)) {
262 exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
263 --i;
264 }
265 }
266
267 // Replace exclusive lvalue branch servers with lvalue branches
268 // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
269 if (exclLVServers.getSize()>0 && !function.isBinnedDistribution(exclLVBranches)) {
270// cout << "activating LVservers " << exclLVServers << " for use in integration " << endl ;
271 intDepList.remove(exclLVServers) ;
272 intDepList.add(exclLVBranches) ;
273
274 //cout << "intDepList removing exclLVServers " << exclLVServers << endl ;
275 //cout << "intDepList adding exclLVBranches " << exclLVBranches << endl ;
276
277 }
278
279
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
281 // * C) Check for dependents that the PDF insists on integrating *
282 // analytically iself *
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
284
285 RooArgSet anIntOKDepList ;
286 for (auto arg : intDepList) {
287 if (function.forceAnalyticalInt(*arg)) {
288 anIntOKDepList.add(*arg) ;
289 }
290 }
291
292 if (anIntOKDepList.getSize()>0) {
293 oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << endl ;
294 }
295
296
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
298 // * D) Make list of servers that can be integrated analytically *
299 // Add all parameters/dependents as value/shape servers *
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
301
302 for (const auto arg : function.servers()) {
303
304 //cout << "considering server" << arg->GetName() << endl ;
305
306 // Dependent or parameter?
307 if (!arg->dependsOnValue(intDepList)) {
308
309 //cout << " server does not depend on observables, adding server as value server to integral" << endl ;
310
311 if (function.dependsOnValue(*arg)) {
312 addServer(*arg,kTRUE,kFALSE) ;
313 }
314
315 continue ;
316
317 } else {
318
319 // Add final dependents of arg as shape servers
320 RooArgSet argLeafServers ;
321 arg->leafNodeServerList(&argLeafServers,0,kFALSE) ;
322
323 //arg->printCompactTree() ;
324 //cout << "leaf nodes of server are " << argLeafServers << " depList = " << depList << endl ;
325
326 // Skip arg if it is neither value or shape server
327 if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
328 //cout << " server is neither value not shape server of function, ignoring" << endl ;
329 continue ;
330 }
331
332 for (const auto leaf : argLeafServers) {
333
334 //cout << " considering leafnode " << leaf->GetName() << " of server " << arg->GetName() << endl ;
335
336 if (depList.find(leaf->GetName()) && function.dependsOnValue(*leaf)) {
337
338 RooAbsRealLValue* leaflv = dynamic_cast<RooAbsRealLValue*>(leaf) ;
339 if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
340 oocxcoutD(&function,Integration) << function.GetName() << " : Observable " << leaf->GetName() << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf" << endl ;
341 if (leaflv->getBinning(rangeName).lowBoundFunc()) {
342 addServer(*leaflv->getBinning(rangeName).lowBoundFunc(),kTRUE,kFALSE) ;
343 }
344 if(leaflv->getBinning(rangeName).highBoundFunc()) {
345 addServer(*leaflv->getBinning(rangeName).highBoundFunc(),kTRUE,kFALSE) ;
346 }
347 } else {
348 oocxcoutD(&function,Integration) << function.GetName() << ": Adding observable " << leaf->GetName() << " of server "
349 << arg->GetName() << " as shape dependent" << endl ;
350 addServer(*leaf,kFALSE,kTRUE) ;
351 }
352 } else if (!depList.find(leaf->GetName())) {
353
354 if (function.dependsOnValue(*leaf)) {
355 oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as value dependent" << endl ;
356 addServer(*leaf,kTRUE,kFALSE) ;
357 } else {
358 oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as shape dependent" << endl ;
359 addServer(*leaf,kFALSE,kTRUE) ;
360 }
361 }
362 }
363 }
364
365 // If this dependent arg is self-normalized, stop here
366 //if (function.selfNormalized()) continue ;
367
368 Bool_t depOK(kFALSE) ;
369 // Check for integratable AbsRealLValue
370
371 //cout << "checking server " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
372
373 if (arg->isDerived()) {
374 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(arg) ;
375 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(arg) ;
376 //cout << "realArgLV = " << realArgLV << " intDepList = " << intDepList << endl ;
377 if ((realArgLV && intDepList.find(realArgLV->GetName()) && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
378
379 //cout << " arg " << arg->GetName() << " is derived LValue with valid jacobian" << endl ;
380
381 // Derived LValue with valid jacobian
382 depOK = kTRUE ;
383
384 // Now, check for overlaps
385 Bool_t overlapOK = kTRUE ;
386 for (const auto otherArg : function.servers()) {
387 // skip comparison with self
388 if (arg==otherArg) continue ;
389 if (otherArg->IsA()==RooConstVar::Class()) continue ;
390 if (arg->overlaps(*otherArg,kTRUE)) {
391 //cout << "arg " << arg->GetName() << " overlaps with " << otherArg->GetName() << endl ;
392 //overlapOK=kFALSE ;
393 }
394 }
395 // coverity[DEADCODE]
396 if (!overlapOK) depOK=kFALSE ;
397
398 //cout << "overlap check returns OK=" << (depOK?"T":"F") << endl ;
399
400 }
401 } else {
402 // Fundamental types are always OK
403 depOK = kTRUE ;
404 }
405
406 // Add server to list of dependents that are OK for analytical integration
407 if (depOK) {
408 anIntOKDepList.add(*arg,kTRUE) ;
409 oocxcoutI(&function,Integration) << function.GetName() << ": Observable " << arg->GetName() << " is suitable for analytical integration (if supported by p.d.f)" << endl ;
410 }
411 }
412 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
413 // * E) interact with function to make list of objects actually integrated analytically *
414 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
415
416 RooArgSet anIntDepList ;
417
418 RooArgSet *anaSet = new RooArgSet( _anaList, Form("UniqueCloneOf_%s",_anaList.GetName()));
421 _anaList.add(*anaSet);
422 delete anaSet;
423
424 // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaListx
425 if (_mode==0) {
427 }
428
429 if (_mode!=0) {
430 oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << endl ;
431 }
432
433
434 // WVE kludge: synchronize dset for use in analyticalIntegral
435 function.getVal(_funcNormSet) ;
436
437 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
438 // * F) Make list of numerical integration variables consisting of: *
439 // * - Category dependents of RealLValues in analytical integration *
440 // * - Expanded server lists of server that are not analytically integrated *
441 // * Make Jacobian list with analytically integrated RealLValues *
442 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
443
444 RooArgSet numIntDepList ;
445
446
447 // Loop over actually analytically integrated dependents
448 for (const auto arg : _anaList) {
449
450 // Process only derived RealLValues
451 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class()) && arg->isDerived() && !arg->isFundamental()) {
452
453 // Add to list of Jacobians to calculate
454 _jacList.add(*arg) ;
455
456 // Add category dependent of LValueReal used in integration
457 auto argDepList = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
458 for (const auto argDep : *argDepList) {
459 if (argDep->IsA()->InheritsFrom(RooAbsCategoryLValue::Class()) && intDepList.contains(*argDep)) {
460 numIntDepList.add(*argDep,kTRUE) ;
461 }
462 }
463 }
464 }
465
466
467 // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
468 if (_anaList.getSize()==0) {
469 if (exclLVServers.getSize()>0) {
470 //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << endl ;
471 intDepList.remove(exclLVBranches) ;
472 intDepList.add(exclLVServers) ;
473 }
474 }
475 //cout << "NUMINT intDepList = " << intDepList << endl ;
476
477 // Loop again over function servers to add remaining numeric integrations
478 for (const auto arg : function.servers()) {
479
480 //cout << "processing server for numeric integration " << arg->IsA()->GetName() << "::" << arg->GetName() << endl ;
481
482 // Process only servers that are not treated analytically
483 if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
484
485 // Process only derived RealLValues
486 if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
487 numIntDepList.add(*arg,kTRUE) ;
488 } else {
489
490 // WVE this will only get the observables, but not l-value transformations
491 // Expand server in final dependents
492 auto argDeps = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
493
494 if (argDeps->getSize()>0) {
495
496 // Add final dependents, that are not forcibly integrated analytically,
497 // to numerical integration list
498 for (const auto dep : *argDeps) {
499 if (!_anaList.find(dep->GetName())) {
500 numIntDepList.add(*dep,kTRUE) ;
501 }
502 }
503 }
504 }
505 }
506 }
507
508 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
509 // * G) Split numeric list in integration list and summation list *
510 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
511
512 // Split numeric integration list in summation and integration lists
513 for (const auto arg : numIntDepList) {
514 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
515 _intList.add(*arg) ;
516 } else if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
517 _sumList.add(*arg) ;
518 }
519 }
520
521 if (_anaList.getSize()>0) {
522 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << endl ;
523 }
524 if (_intList.getSize()>0) {
525 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << endl ;
526 }
527 if (_sumList.getSize()>0) {
528 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << endl ;
529 }
530
531
532 // Determine operating mode
533 if (numIntDepList.getSize()>0) {
534 // Numerical and optional Analytical integration
536 } else if (_anaList.getSize()>0) {
537 // Purely analytical integration
539 } else {
540 // No integration performed
542 }
543
544 // Determine auto-dirty status
546
547 // Create value caches for _intList and _sumList
550
551
552 if (_sumList.getSize()>0) {
553 RooSuperCategory *sumCat = new RooSuperCategory(Form("%s_sumCat",GetName()),"sumCat",_sumList) ;
554 _sumCatIter = sumCat->typeIterator() ;
555 _sumCat.addOwned(*sumCat) ;
556 }
557
559}
560
561
562
563////////////////////////////////////////////////////////////////////////////////
564/// Set appropriate cache operation mode for integral depending on cache operation
565/// mode of server objects
566
568{
569 // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
570 for (const auto server : _serverList) {
571 if (server->isValueServer(*this)) {
572 RooArgSet leafSet ;
573 server->leafNodeServerList(&leafSet) ;
574 for (const auto leaf : leafSet) {
575 if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
577 break ;
578 }
579 if (leaf->getAttribute("projectedDependent")) {
581 break ;
582 }
583 }
584 }
585 }
586}
587
588
589
590////////////////////////////////////////////////////////////////////////////////
591/// Utility function that returns true if 'object server' is a server
592/// to exactly one of the RooAbsArgs in 'exclLVBranches'
593
594Bool_t RooRealIntegral::servesExclusively(const RooAbsArg* server,const RooArgSet& exclLVBranches, const RooArgSet& allBranches) const
595{
596 // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
597
598 // Special case, no LV servers available
599 if (exclLVBranches.getSize()==0) return kFALSE ;
600
601 // If server has no clients and is not an LValue itself, return false
602 if (server->_clientList.empty() && exclLVBranches.find(server->GetName())) {
603 return kFALSE ;
604 }
605
606 // WVE must check for value relations only here!!!!
607
608 // Loop over all clients
609 Int_t numLVServ(0) ;
610 for (const auto client : server->valueClients()) {
611 // If client is not an LValue, recurse
612 if (!(exclLVBranches.find(client->GetName())==client)) {
613 if (allBranches.find(client->GetName())==client) {
614 if (!servesExclusively(client,exclLVBranches,allBranches)) {
615 // Client is a non-LValue that doesn't have an exclusive LValue server
616 return kFALSE ;
617 }
618 }
619 } else {
620 // Client is an LValue
621 numLVServ++ ;
622 }
623 }
624
625 return (numLVServ==1) ;
626}
627
628
629
630
631////////////////////////////////////////////////////////////////////////////////
632/// (Re)Initialize numerical integration engine if necessary. Return kTRUE if
633/// successful, or otherwise kFALSE.
634
636{
637 // if we already have an engine, check if it still works for the present limits.
638 if(0 != _numIntEngine) {
640 // otherwise, cleanup the old engine
641 delete _numIntEngine ;
642 _numIntEngine= 0;
643 if(0 != _numIntegrand) {
644 delete _numIntegrand;
645 _numIntegrand= 0;
646 }
647 }
648
649 // All done if there are no arguments to integrate numerically
650 if(0 == _intList.getSize()) return kTRUE;
651
652 // Bind the appropriate analytic integral (specified by _mode) of our RooRealVar object to
653 // those of its arguments that will be integrated out numerically.
654 if(_mode != 0) {
656 }
657 else {
659 }
660 if(0 == _numIntegrand || !_numIntegrand->isValid()) {
661 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << endl;
662 return kFALSE;
663 }
664
665 // Create appropriate numeric integrator using factory
668
669 if(0 == _numIntEngine || !_numIntEngine->isValid()) {
670 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << endl;
671 return kFALSE;
672 }
673
674 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
675 << _numIntEngine->IsA()->GetName() << " to calculate Int" << _intList << endl ;
676
677 if (_intList.getSize()>3) {
678 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 ;
679 }
680
682 return kTRUE;
683}
684
685
686
687////////////////////////////////////////////////////////////////////////////////
688/// Copy constructor
689
691 RooAbsReal(other,name),
692 _valid(other._valid),
693 _respectCompSelect(other._respectCompSelect),
694 _sumList("!sumList",this,other._sumList),
695 _intList("!intList",this,other._intList),
696 _anaList("!anaList",this,other._anaList),
697 _jacList("!jacList",this,other._jacList),
698 _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
699 _function("!func",this,other._function),
700 _iconfig(other._iconfig),
701 _sumCat("!sumCat",this,other._sumCat),
702 _sumCatIter(0),
703 _mode(other._mode),
704 _intOperMode(other._intOperMode),
705 _restartNumIntEngine(kFALSE),
706 _numIntEngine(0),
707 _numIntegrand(0),
708 _rangeName(other._rangeName),
709 _params(0),
710 _cacheNum(kFALSE)
711{
713
714 for (const auto arg : other._facList) {
715 RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
716 _facListOwned.addOwned(*argClone) ;
717 _facList.add(*argClone) ;
718 addServer(*argClone,kFALSE,kTRUE) ;
719 }
720
721 other._intList.snapshot(_saveInt) ;
722 other._sumList.snapshot(_saveSum) ;
723
725}
726
727
728
729////////////////////////////////////////////////////////////////////////////////
730
732 // Destructor
733{
734 if (_numIntEngine) delete _numIntEngine ;
735 if (_numIntegrand) delete _numIntegrand ;
736 if (_funcNormSet) delete _funcNormSet ;
737 if (_sumCatIter) delete _sumCatIter ;
738 if (_params) delete _params ;
739
741}
742
743
744
745
746
747////////////////////////////////////////////////////////////////////////////////
748
749RooAbsReal* RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const
750{
751 // Handle special case of no integration with default algorithm
752 if (iset.getSize()==0) {
753 return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
754 }
755
756 // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
757 RooArgSet isetAll(iset) ;
758 isetAll.add(_sumList) ;
759 isetAll.add(_intList) ;
760 isetAll.add(_anaList) ;
761 isetAll.add(_facList) ;
762
763 const RooArgSet* newNormSet(0) ;
764 RooArgSet* tmp(0) ;
765 if (nset && !_funcNormSet) {
766 newNormSet = nset ;
767 } else if (!nset && _funcNormSet) {
768 newNormSet = _funcNormSet ;
769 } else if (nset && _funcNormSet) {
770 tmp = new RooArgSet ;
771 tmp->add(*nset) ;
772 tmp->add(*_funcNormSet,kTRUE) ;
773 newNormSet = tmp ;
774 }
775 RooAbsReal* ret = _function.arg().createIntegral(isetAll,newNormSet,cfg,rangeName) ;
776
777 if (tmp) {
778 delete tmp ;
779 }
780
781 return ret ;
782}
783
784
785
786////////////////////////////////////////////////////////////////////////////////
787/// Return value of object. If the cache is clean, return the
788/// cached value, otherwise recalculate on the fly and refill
789/// the cache
790
792{
793// // fast-track clean-cache processing
794// if (_operMode==AClean) {
795// return _value ;
796// }
797
798 if (nset && nset!=_lastNSet) {
799 ((RooAbsReal*) this)->setProxyNormSet(nset) ;
800 _lastNSet = (RooArgSet*) nset ;
801 }
802
804 _value = traceEval(nset) ;
805 }
806
807 return _value ;
808}
809
810
811
812
813
814////////////////////////////////////////////////////////////////////////////////
815/// Perform the integration and return the result
816
818{
820
821 Double_t retVal(0) ;
822 switch (_intOperMode) {
823
824 case Hybrid:
825 {
826 // Cache numeric integrals in >1d expensive object cache
827 RooDouble* cacheVal(0) ;
830 }
831
832 if (cacheVal) {
833 retVal = *cacheVal ;
834 // cout << "using cached value of integral" << GetName() << endl ;
835 } else {
836
837
838 // Find any function dependents that are AClean
839 // and switch them temporarily to ADirty
840 Bool_t origState = inhibitDirty() ;
842
843 // try to initialize our numerical integration engine
844 if(!(_valid= initNumIntegrator())) {
845 coutE(Integration) << ClassName() << "::" << GetName()
846 << ":evaluate: cannot initialize numerical integrator" << endl;
847 return 0;
848 }
849
850 // Save current integral dependent values
853
854 // Evaluate sum/integral
855 retVal = sum() ;
856
857
858 // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
859 setDirtyInhibit(origState) ;
860
861 // Restore integral dependent values
864
865 // Cache numeric integrals in >1d expensive object cache
867 RooDouble* val = new RooDouble(retVal) ;
869 // cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << endl ;
870 }
871
872 }
873 break ;
874 }
875 case Analytic:
876 {
878 cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
879 << ")func = " << _function.arg().IsA()->GetName() << "::" << _function.arg().GetName()
880 << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << endl ;
881
882
883 break ;
884 }
885
886 case PassThrough:
887 {
888 //setDirtyInhibit(kTRUE) ;
889 retVal= _function.arg().getVal(_funcNormSet) ;
890 //setDirtyInhibit(kFALSE) ;
891 break ;
892 }
893 }
894
895
896 // Multiply answer with integration ranges of factorized variables
897 if (_facList.getSize()>0) {
898 for (const auto arg : _facList) {
899 // Multiply by fit range for 'real' dependents
900 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
901 RooAbsRealLValue* argLV = (RooAbsRealLValue*)arg ;
902 retVal *= (argLV->getMax() - argLV->getMin()) ;
903 }
904 // Multiply by number of states for category dependents
905 if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
907 retVal *= argLV->numTypes() ;
908 }
909 }
910 }
911
912
913 if (dologD(Tracing)) {
914 cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
915 switch(_mode) {
916 case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
917 case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
918 case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
919 }
920
921 ccxcoutD(Tracing) << "raw*fact = " << retVal << endl ;
922 }
923
924 return retVal ;
925}
926
927
928
929////////////////////////////////////////////////////////////////////////////////
930/// Return product of jacobian terms originating from analytical integration
931
933{
934 if (_jacList.getSize()==0) {
935 return 1 ;
936 }
937
938 Double_t jacProd(1) ;
939 for (const auto elm : _jacList) {
940 auto arg = static_cast<const RooAbsRealLValue*>(elm);
941 jacProd *= arg->jacobian() ;
942 }
943
944 // Take fabs() here: if jacobian is negative, min and max are swapped and analytical integral
945 // will be positive, so must multiply with positive jacobian.
946 return fabs(jacProd) ;
947}
948
949
950
951////////////////////////////////////////////////////////////////////////////////
952/// Perform summation of list of category dependents to be integrated
953
955{
956 if (_sumList.getSize()!=0) {
957 // Add integrals for all permutations of categories summed over
958 Double_t total(0) ;
959
960 _sumCatIter->Reset() ;
963 while((type=(RooCatType*)_sumCatIter->Next())) {
964 sumCat->setIndex(type->getVal()) ;
965 if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
967 }
968 }
969
970 return total ;
971
972 } else {
973 // Simply return integral
975 return ret ;
976 }
977}
978
979
980////////////////////////////////////////////////////////////////////////////////
981/// Perform hybrid numerical/analytical integration over all real-valued dependents
982
984{
985 if (!_numIntEngine) {
986 // Trivial case, fully analytical integration
987 return ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) ;
988 } else {
989 return _numIntEngine->calculate() ;
990 }
991}
992
993
994
995////////////////////////////////////////////////////////////////////////////////
996/// Intercept server redirects and reconfigure internal object accordingly
997
999 Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
1000{
1002
1004
1005 // Update contents value caches for _intList and _sumList
1010
1011 // Delete parameters cache if we have one
1012 if (_params) {
1013 delete _params ;
1014 _params = 0 ;
1015 }
1016
1017 return kFALSE ;
1018}
1019
1020
1021
1022////////////////////////////////////////////////////////////////////////////////
1023
1025{
1026 if (!_params) {
1027 _params = new RooArgSet("params") ;
1028
1029 RooArgSet params ;
1030 for (const auto server : _serverList) {
1031 if (server->isValueServer(*this)) _params->add(*server) ;
1032 }
1033 }
1034
1035 return *_params ;
1036}
1037
1038
1039
1040////////////////////////////////////////////////////////////////////////////////
1041/// Dummy
1042
1044{
1045 if (_operMode==ADirty) {
1046// cout << "RooRealIntegral::operModeHook(" << GetName() << " warning: mode set to ADirty" << endl ;
1047// if (TString(GetName()).Contains("FULL")) {
1048// cout << "blah" << endl ;
1049// }
1050 }
1051}
1052
1053
1054
1055////////////////////////////////////////////////////////////////////////////////
1056/// Check if current value is valid
1057
1058Bool_t RooRealIntegral::isValidReal(Double_t /*value*/, Bool_t /*printError*/) const
1059{
1060 return kTRUE ;
1061}
1062
1063////////////////////////////////////////////////////////////////////////////////
1064/// Check if component selection is allowed
1065
1067 return _respectCompSelect;
1068}
1069
1070////////////////////////////////////////////////////////////////////////////////
1071/// Set component selection to be allowed/forbidden
1072
1074 _respectCompSelect = allow;
1075}
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1079/// integration operation
1080
1081void RooRealIntegral::printMetaArgs(ostream& os) const
1082{
1083
1084 if (intVars().getSize()!=0) {
1085 os << "Int " ;
1086 }
1087 os << _function.arg().GetName() ;
1088 if (_funcNormSet) {
1089 os << "_Norm" ;
1090 os << *_funcNormSet ;
1091 os << " " ;
1092 }
1093
1094 // List internally integrated observables and factorizing observables as analytically integrated
1095 RooArgSet tmp(_anaList) ;
1096 tmp.add(_facList) ;
1097 if (tmp.getSize()>0) {
1098 os << "d[Ana]" ;
1099 os << tmp ;
1100 os << " " ;
1101 }
1102
1103 // List numerically integrated and summed observables as numerically integrated
1104 RooArgSet tmp2(_intList) ;
1105 tmp2.add(_sumList) ;
1106 if (tmp2.getSize()>0) {
1107 os << " d[Num]" ;
1108 os << tmp2 ;
1109 os << " " ;
1110 }
1111}
1112
1113
1114
1115////////////////////////////////////////////////////////////////////////////////
1116/// Print the state of this object to the specified output stream.
1117
1119{
1121 os << indent << "--- RooRealIntegral ---" << endl;
1122 os << indent << " Integrates ";
1124 TString deeper(indent);
1125 deeper.Append(" ");
1126 os << indent << " operating mode is "
1127 << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << endl ;
1128 os << indent << " Summed discrete args are " << _sumList << endl ;
1129 os << indent << " Numerically integrated args are " << _intList << endl;
1130 os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << endl ;
1131 os << indent << " Arguments included in Jacobian are " << _jacList << endl ;
1132 os << indent << " Factorized arguments are " << _facList << endl ;
1133 os << indent << " Function normalization set " ;
1134 if (_funcNormSet)
1135 _funcNormSet->Print("1") ;
1136 else
1137 os << "<none>" ;
1138
1139 os << endl ;
1140}
1141
1142
1143
1144////////////////////////////////////////////////////////////////////////////////
1145/// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1146
1148 _cacheAllNDim = ndim ;
1149}
1150
1151
1152////////////////////////////////////////////////////////////////////////////////
1153/// Return minimum dimensions of numeric integration for which values are cached.
1154
1156{
1157 return _cacheAllNDim ;
1158}
1159
1160
void Class()
Definition: Class.C:29
#define cxcoutI(a)
Definition: RooMsgService.h:86
#define cxcoutD(a)
Definition: RooMsgService.h:82
#define oocxcoutD(o, a)
Definition: RooMsgService.h:84
#define dologD(a)
Definition: RooMsgService.h:66
#define coutE(a)
Definition: RooMsgService.h:34
#define ccxcoutD(a)
Definition: RooMsgService.h:83
#define ccoutD(a)
Definition: RooMsgService.h:38
#define oocxcoutI(o, a)
Definition: RooMsgService.h:88
#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:365
static void indent(ostringstream &buf, int indent_level)
static unsigned int total
char name[80]
Definition: TGX11.cxx:109
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:71
RooExpensiveObjectCache & expensiveObjectCache() const
Definition: RooAbsArg.cxx:2176
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:240
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:734
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:83
static void setDirtyInhibit(Bool_t flag)
Control global dirty inhibit mode.
Definition: RooAbsArg.cxx:234
friend class RooArgSet
Definition: RooAbsArg.h:551
Bool_t inhibitDirty() const
Delete watch flag.
Definition: RooAbsArg.cxx:85
const RefCountList_t & valueClients() const
List of all value clients of this object. Value clients receive value updates.
Definition: RooAbsArg.h:153
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition: RooAbsArg.h:500
void addServer(RooAbsArg &server, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
Definition: RooAbsArg.cxx:358
Bool_t isValueOrShapeDirtyAndClear() const
Definition: RooAbsArg.h:424
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Change cache operation mode to given mode.
Definition: RooAbsArg.cxx:1726
OperMode _operMode
Mark batches as dirty (only meaningful for RooAbsReal).
Definition: RooAbsArg.h:635
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Definition: RooAbsArg.h:97
RefCountList_t _clientList
Definition: RooAbsArg.h:560
RefCountList_t _serverList
Definition: RooAbsArg.h:559
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...
Int_t numTypes(const char *=0) const
Return number of types defined (in range named rangeName if rangeName!=0)
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
Storage_t::size_type size() const
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.
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
Get maximum of currently defined range.
virtual const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const =0
Retrive binning configuration with given name or default binning.
virtual Bool_t isJacobianOK(const RooArgSet &depList) const
virtual Double_t jacobian() const
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
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:373
friend class RooRealBinding
Definition: RooAbsReal.h:442
RooArgSet * _lastNSet
Definition: RooAbsReal.h:555
Double_t traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
Definition: RooAbsReal.cxx:341
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:402
Double_t _value
Definition: RooAbsReal.h:446
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Definition: RooAbsReal.cxx:477
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
static Bool_t _globalSelectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition: RooAbsReal.h:536
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:562
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:309
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:134
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,...
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::addOwned()
RooNameReg is a registry for const char* names.
Definition: RooNameReg.h:24
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition: RooNameReg.cxx:103
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) const
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
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 * _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
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
bool empty() const
Check if empty.
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 can join 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...
const T & arg() const
Return reference to object held in proxy.
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
Basic string class.
Definition: TString.h:131
TString & Append(const char *cs)
Definition: TString.h:559
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:69
@ InputArguments
Definition: RooGlobalFunc.h:68
@ Integration
Definition: RooGlobalFunc.h:67