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