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