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