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 "RooMsgService.h"
34#include "RooArgSet.h"
35#include "RooAbsRealLValue.h"
37#include "RooRealBinding.h"
38#include "RooInvTransform.h"
39#include "RooSuperCategory.h"
40#include "RooNumIntFactory.h"
41#include "RooNumIntConfig.h"
42#include "RooNameReg.h"
44#include "RooConstVar.h"
45#include "RooDouble.h"
46#include "RooTrace.h"
47#include "RooHelpers.h"
48
49#include "TClass.h"
50
51#include <iostream>
52#include <memory>
53
55
56namespace {
57
58struct ServerToAdd {
59 ServerToAdd(RooAbsArg *theArg, bool isShape) : arg{theArg}, isShapeServer{isShape} {}
60 RooAbsArg *arg = nullptr;
61 bool isShapeServer = false;
62};
63
64void addObservableToServers(RooAbsReal const &function, RooAbsArg &leaf, std::vector<ServerToAdd> &serversToAdd,
65 const char *rangeName)
66{
67 auto leaflv = dynamic_cast<RooAbsRealLValue *>(&leaf);
68 if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
69 oocxcoutD(&function, Integration)
70 << function.GetName() << " : Observable " << leaf.GetName()
71 << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf"
72 << std::endl;
73 if (leaflv->getBinning(rangeName).lowBoundFunc()) {
74 serversToAdd.emplace_back(leaflv->getBinning(rangeName).lowBoundFunc(), false);
75 }
76 if (leaflv->getBinning(rangeName).highBoundFunc()) {
77 serversToAdd.emplace_back(leaflv->getBinning(rangeName).highBoundFunc(), false);
78 }
79 } else {
80 oocxcoutD(&function, Integration) << function.GetName() << ": Adding observable " << leaf.GetName()
81 << " as shape dependent" << std::endl;
82 serversToAdd.emplace_back(&leaf, true);
83 }
84}
85
86void addParameterToServers(RooAbsReal const &function, RooAbsArg &leaf, std::vector<ServerToAdd> &serversToAdd,
87 bool isShapeServer)
88{
89 if (!isShapeServer) {
90 oocxcoutD(&function, Integration) << function.GetName() << ": Adding parameter " << leaf.GetName()
91 << " as value dependent" << std::endl;
92 } else {
93 oocxcoutD(&function, Integration) << function.GetName() << ": Adding parameter " << leaf.GetName()
94 << " as shape dependent" << std::endl;
95 }
96 serversToAdd.emplace_back(&leaf, isShapeServer);
97}
98
99enum class MarkedState { Dependent, Independent, AlreadyAdded };
100
101/// Mark all args that recursively are value clients of "dep".
102void unmarkDepValueClients(RooAbsArg const &dep, RooArgSet const &args, std::vector<MarkedState> &marked)
103{
104 assert(args.size() == marked.size());
105 auto index = args.index(dep);
106 if (index >= 0) {
107 marked[index] = MarkedState::Dependent;
108 for (RooAbsArg *client : dep.valueClients()) {
109 unmarkDepValueClients(*client, args, marked);
110 }
111 }
112}
113
114std::vector<ServerToAdd>
115getValueAndShapeServers(RooAbsReal const &function, RooArgSet const &depList, const char *rangeName)
116{
117 std::vector<ServerToAdd> serversToAdd;
118
119 // Get the full computation graph and sort it topologically
120 RooArgList allArgsList;
121 function.treeNodeServerList(&allArgsList, nullptr, true, true, /*valueOnly=*/false, false);
122 RooArgSet allArgs{allArgsList};
123 allArgs.sortTopologically();
124
125 // Figure out what are all the value servers only
126 RooArgList allValueArgsList;
127 function.treeNodeServerList(&allValueArgsList, nullptr, true, true, /*valueOnly=*/true, false);
128 RooArgSet allValueArgs{allValueArgsList};
129
130 // All "marked" args will be added as value servers to the integral
131 std::vector<MarkedState> marked(allArgs.size(), MarkedState::Independent);
132 marked.back() = MarkedState::Dependent; // We don't want to consider the function itself
133
134 // Mark all args that are (indirect) value servers of the integration
135 // variable or the integration variable itself. If something was marked,
136 // it means the integration variable was in the compute graph and we will
137 // add it to the server list.
138 for (RooAbsArg *dep : depList) {
139 if (RooAbsArg * depInArgs = allArgs.find(dep->GetName())) {
140 unmarkDepValueClients(*depInArgs, allArgs, marked);
141 addObservableToServers(function, *depInArgs, serversToAdd, rangeName);
142 }
143 }
144
145 // We are adding all independent direct servers of the args depending on the
146 // integration variables
147 for (std::size_t i = 0; i < allArgs.size(); ++i) {
148 if(marked[i] == MarkedState::Dependent) {
149 for(RooAbsArg* server : allArgs[i]->servers()) {
150 int index = allArgs.index(server->GetName());
151 if(index >= 0 && marked[index] == MarkedState::Independent) {
152 addParameterToServers(function, *server, serversToAdd, !allValueArgs.find(*server));
153 marked[index] = MarkedState::AlreadyAdded;
154 }
155 }
156 }
157 }
158
159 return serversToAdd;
160}
161
162void fillAnIntOKDepList(RooAbsReal const &function, RooArgSet const &intDepList, RooArgSet &anIntOKDepList)
163{
164 for (const auto arg : function.servers()) {
165
166 // Dependent or parameter?
167 if (!arg->dependsOnValue(intDepList)) {
168 continue;
169 } else if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
170 // Skip arg if it is neither value or shape server
171 continue;
172 }
173
174 bool depOK(false);
175 // Check for integratable AbsRealLValue
176
177 if (arg->isDerived()) {
178 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue *>(arg);
179 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue *>(arg);
180 if ((realArgLV && intDepList.find(realArgLV->GetName()) && (realArgLV->isJacobianOK(intDepList) != 0)) ||
181 catArgLV) {
182
183 // Derived LValue with valid jacobian
184 depOK = true;
185
186 // Now, check for overlaps
187 bool overlapOK = true;
188 for (const auto otherArg : function.servers()) {
189 // skip comparison with self
190 if (arg == otherArg)
191 continue;
192 if (otherArg->IsA() == RooConstVar::Class())
193 continue;
194 if (arg->overlaps(*otherArg, true)) {
195 }
196 }
197 // coverity[DEADCODE]
198 if (!overlapOK)
199 depOK = false;
200 }
201 } else {
202 // Fundamental types are always OK
203 depOK = true;
204 }
205
206 // Add server to list of dependents that are OK for analytical integration
207 if (depOK) {
208 anIntOKDepList.add(*arg, true);
209 oocxcoutI(&function, Integration) << function.GetName() << ": Observable " << arg->GetName()
210 << " is suitable for analytical integration (if supported by p.d.f)"
211 << std::endl;
212 }
213 }
214}
215
216} // namespace
217
218
219using namespace std;
220
221
223
224
225////////////////////////////////////////////////////////////////////////////////
226
228{
230}
231
232
233
234////////////////////////////////////////////////////////////////////////////////
235/// Construct integral of 'function' over observables in 'depList'
236/// in range 'rangeName' with normalization observables 'funcNormSet'
237/// (for p.d.f.s). In the integral is performed to the maximum extent
238/// possible the internal (analytical) integrals advertised by function.
239/// The other integrations are performed numerically. The optional
240/// config object prescribes how these numeric integrations are configured.
241///
242/// \Note If pdf component selection was globally overridden to always include
243/// all components (either with RooAbsReal::globalSelectComp(bool) or a
244/// RooAbsReal::GlobalSelectComponentRAII), then any created integral will
245/// ignore component selections during its lifetime. This is especially useful
246/// when creating normalization or projection integrals.
247RooRealIntegral::RooRealIntegral(const char *name, const char *title,
248 const RooAbsReal& function, const RooArgSet& depList,
249 const RooArgSet* funcNormSet, const RooNumIntConfig* config,
250 const char* rangeName) :
251 RooAbsReal(name,title),
252 _valid(true),
253 _respectCompSelect{!_globalSelectComp},
254 _sumList("!sumList","Categories to be summed numerically",this,false,false),
255 _intList("!intList","Variables to be integrated numerically",this,false,false),
256 _anaList("!anaList","Variables to be integrated analytically",this,false,false),
257 _jacList("!jacList","Jacobian product term",this,false,false),
258 _facList("!facList","Variables independent of function",this,false,true),
259 _function("!func","Function to be integrated",this,false,false),
260 _iconfig((RooNumIntConfig*)config),
261 _sumCat("!sumCat","SuperCategory for summation",this,false,false),
262 _rangeName((TNamed*)RooNameReg::ptr(rangeName))
263{
264 // A) Check that all dependents are lvalues
265 //
266 // B) Check if list of dependents can be re-expressed in
267 // lvalues that are higher in the expression tree
268 //
269 // C) Check for dependents that the PDF insists on integrating
270 // analytically itself
271 //
272 // D) Make list of servers that can be integrated analytically
273 // Add all parameters/dependents as value/shape servers
274 //
275 // E) Interact with function to make list of objects actually integrated analytically
276 //
277 // F) Make list of numerical integration variables consisting of:
278 // - Category dependents of RealLValues in analytical integration
279 // - Leaf nodes server lists of function server that are not analytically integrated
280 // - Make Jacobian list for analytically integrated RealLValues
281 //
282 // G) Split numeric list in integration list and summation list
283 //
284
285 oocxcoutI(&function,Integration) << "RooRealIntegral::ctor(" << GetName() << ") Constructing integral of function "
286 << function.GetName() << " over observables" << depList << " with normalization "
287 << (funcNormSet?*funcNormSet:RooArgSet()) << " with range identifier "
288 << (rangeName?rangeName:"<none>") << std::endl ;
289
290
291 // Choose same expensive object cache as integrand
292 setExpensiveObjectCache(function.expensiveObjectCache()) ;
293// cout << "RRI::ctor(" << GetName() << ") setting expensive object cache to " << &expensiveObjectCache() << " as taken from " << function.GetName() << std::endl ;
294
295 // Use objects integrator configuration if none is specified
296 if (!_iconfig) _iconfig = (RooNumIntConfig*) function.getIntegratorConfig() ;
297
298 // Save private copy of funcNormSet, if supplied, excluding factorizing terms
299 if (funcNormSet) {
300 _funcNormSet = std::make_unique<RooArgSet>();
301 for (const auto nArg : *funcNormSet) {
302 if (function.dependsOn(*nArg)) {
303 _funcNormSet->addClone(*nArg) ;
304 }
305 }
306 }
307
308 //_funcNormSet = funcNormSet ? (RooArgSet*)funcNormSet->snapshot(false) : 0 ;
309
310 // Make internal copy of dependent list
311 RooArgSet intDepList(depList) ;
312
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
314 // * A) Check that all dependents are lvalues and filter out any
315 // dependents that the PDF doesn't explicitly depend on
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
317
318 for (auto arg : intDepList) {
319 if(!arg->isLValue()) {
320 coutE(InputArguments) << ClassName() << "::" << GetName() << ": cannot integrate non-lvalue ";
321 arg->Print("1");
322 _valid= false;
323 }
324 if (!function.dependsOn(*arg)) {
325 std::unique_ptr<RooAbsArg> argClone{static_cast<RooAbsArg*>(arg->Clone())};
326 _facList.add(*argClone) ;
327 addServer(*argClone,false,true) ;
328 _facListOwned.addOwned(std::move(argClone));
329 }
330 }
331
332 if (!_facList.empty()) {
333 oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables are " << _facList << std::endl ;
334 }
335
336
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
338 // * B) Check if list of dependents can be re-expressed in *
339 // * lvalues that are higher in the expression tree *
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
341
342
343 // Initial fill of list of LValue branches
344 RooArgSet exclLVBranches("exclLVBranches") ;
345 RooArgSet branchList,branchListVD ;
346 function.branchNodeServerList(&branchList) ;
347
348 for (auto branch: branchList) {
349 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(branch) ;
350 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(branch) ;
351 if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
352 exclLVBranches.add(*branch) ;
353// cout << "exclv branch = " << std::endl ;
354// branch->printCompactTree() ;
355 }
356 if (dependsOnValue(*branch)) {
357 branchListVD.add(*branch) ;
358 } else {
359// cout << "value of self does not depend on branch " << branch->GetName() << std::endl ;
360 }
361 }
362 exclLVBranches.remove(depList,true,true) ;
363
364 // Initial fill of list of LValue leaf servers (put in intDepList, but the
365 // instances that are in the actual computation graph of the function)
366 RooArgSet exclLVServers("exclLVServers") ;
367 function.getObservables(&intDepList, exclLVServers);
368
369 // Obtain mutual exclusive dependence by iterative reduction
370 bool converged(false) ;
371 while(!converged) {
372 converged=true ;
373
374 // Reduce exclLVServers to only those serving exclusively exclLVBranches
375 std::vector<RooAbsArg*> toBeRemoved;
376 for (auto server : exclLVServers) {
377 if (!servesExclusively(server,exclLVBranches,branchListVD)) {
378 toBeRemoved.push_back(server);
379 converged=false ;
380 }
381 }
382 exclLVServers.remove(toBeRemoved.begin(), toBeRemoved.end());
383
384 // Reduce exclLVBranches to only those depending exclusively on exclLVservers
385 // Attention: counting loop, since erasing from container
386 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
387 const RooAbsArg* branch = exclLVBranches[i];
388 RooArgSet brDepList;
389 branch->getObservables(&intDepList, brDepList);
390 RooArgSet bsList(brDepList,"bsList") ;
391 bsList.remove(exclLVServers,true,true) ;
392 if (!bsList.empty()) {
393 exclLVBranches.remove(*branch,true,true) ;
394 --i;
395 converged=false ;
396 }
397 }
398 }
399
400 // Eliminate exclLVBranches that do not depend on any LVServer
401 // Attention: Counting loop, since modifying container
402 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
403 const RooAbsArg* branch = exclLVBranches[i];
404 if (!branch->dependsOnValue(exclLVServers)) {
405 exclLVBranches.remove(*branch,true,true) ;
406 --i;
407 }
408 }
409
410 // Replace exclusive lvalue branch servers with lvalue branches
411 // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
412 if (!exclLVServers.empty() && !function.isBinnedDistribution(exclLVBranches)) {
413// cout << "activating LVservers " << exclLVServers << " for use in integration " << std::endl ;
414 intDepList.remove(exclLVServers) ;
415 intDepList.add(exclLVBranches) ;
416
417 //cout << "intDepList removing exclLVServers " << exclLVServers << std::endl ;
418 //cout << "intDepList adding exclLVBranches " << exclLVBranches << std::endl ;
419
420 }
421
422
423 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
424 // * C) Check for dependents that the PDF insists on integrating *
425 // analytically itself *
426 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
427
428 RooArgSet anIntOKDepList ;
429 for (auto arg : intDepList) {
430 if (function.forceAnalyticalInt(*arg)) {
431 anIntOKDepList.add(*arg) ;
432 }
433 }
434
435 if (!anIntOKDepList.empty()) {
436 oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << std::endl ;
437 }
438
439
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
441 // * D) Make list of servers that can be integrated analytically *
442 // Add all parameters/dependents as value/shape servers *
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
444
445 auto serversToAdd = getValueAndShapeServers(function, depList, rangeName);
446 fillAnIntOKDepList(function, intDepList, anIntOKDepList);
447 // We will not add the servers just now, because it makes only sense to add
448 // them once we have made sure that this integral is not operating in
449 // pass-through mode. It will be done at the end of this constructor.
450
451 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
452 // * E) interact with function to make list of objects actually integrated analytically *
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
454
455 _mode = function.getAnalyticalIntegralWN(anIntOKDepList,_anaList,_funcNormSet.get(),RooNameReg::str(_rangeName)) ;
456
457 // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaList
458 if (_mode==0) {
460 }
461
462 if (_mode!=0) {
463 oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << std::endl ;
464 }
465
466 // WVE kludge: synchronize dset for use in analyticalIntegral
467 // LM : I think this is needed only if _funcNormSet is not an empty set
468 if (_funcNormSet && !_funcNormSet->empty()) {
469 function.getVal(_funcNormSet.get()) ;
470 }
471
472 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
473 // * F) Make list of numerical integration variables consisting of: *
474 // * - Category dependents of RealLValues in analytical integration *
475 // * - Expanded server lists of server that are not analytically integrated *
476 // * Make Jacobian list with analytically integrated RealLValues *
477 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
478
479 RooArgSet numIntDepList ;
480
481
482 // Loop over actually analytically integrated dependents
483 for (const auto arg : _anaList) {
484
485 // Process only derived RealLValues
486 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class()) && arg->isDerived() && !arg->isFundamental()) {
487
488 // Add to list of Jacobians to calculate
489 _jacList.add(*arg) ;
490
491 // Add category dependent of LValueReal used in integration
492 auto argDepList = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
493 for (const auto argDep : *argDepList) {
494 if (argDep->IsA()->InheritsFrom(RooAbsCategoryLValue::Class()) && intDepList.contains(*argDep)) {
495 numIntDepList.add(*argDep,true) ;
496 }
497 }
498 }
499 }
500
501
502 // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
503 if (_anaList.empty()) {
504 if (!exclLVServers.empty()) {
505 //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << std::endl ;
506 intDepList.remove(exclLVBranches) ;
507 intDepList.add(exclLVServers) ;
508 }
509 }
510 //cout << "NUMINT intDepList = " << intDepList << std::endl ;
511
512 // Loop again over function servers to add remaining numeric integrations
513 for (const auto arg : function.servers()) {
514
515 //cout << "processing server for numeric integration " << arg->ClassName() << "::" << arg->GetName() << std::endl ;
516
517 // Process only servers that are not treated analytically
518 if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
519
520 // Process only derived RealLValues
521 if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
522 numIntDepList.add(*arg,true) ;
523 } else {
524
525 // WVE this will only get the observables, but not l-value transformations
526 // Expand server in final dependents
527 auto argDeps = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
528
529 if (!argDeps->empty()) {
530
531 // Add final dependents, that are not forcibly integrated analytically,
532 // to numerical integration list
533 for (const auto dep : *argDeps) {
534 if (!_anaList.find(dep->GetName())) {
535 numIntDepList.add(*dep,true) ;
536 }
537 }
538 }
539 }
540 }
541 }
542
543 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
544 // * G) Split numeric list in integration list and summation list *
545 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
546
547 // Split numeric integration list in summation and integration lists
548 for (const auto arg : numIntDepList) {
549 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
550 _intList.add(*arg) ;
551 } else if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
552 _sumList.add(*arg) ;
553 }
554 }
555
556 if (!_anaList.empty()) {
557 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << std::endl ;
558 }
559 if (!_intList.empty()) {
560 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << std::endl ;
561 }
562 if (!_sumList.empty()) {
563 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << std::endl ;
564 }
565
566
567 // Determine operating mode
568 if (!numIntDepList.empty()) {
569 // Numerical and optional Analytical integration
571 } else if (!_anaList.empty()) {
572 // Purely analytical integration
574 } else {
575 // No integration performed, where the function is a direct value server
577 _function._valueServer = true;
578 }
579 // We are only setting the function proxy now that it's clear if it's a value
580 // server or not.
581 _function.setArg(const_cast<RooAbsReal&>(function));
582
583 // Determine auto-dirty status
585
586 // Create value caches for _intList and _sumList
589
590
591 if (!_sumList.empty()) {
592 _sumCat.addOwned(std::make_unique<RooSuperCategory>(Form("%s_sumCat",GetName()),"sumCat",_sumList));
593 }
594
595 // Only if we are not in pass-through mode we need to add the shape and value
596 // servers separately.
598 for(auto const& toAdd : serversToAdd) {
599 addServer(*toAdd.arg, !toAdd.isShapeServer, toAdd.isShapeServer);
600 }
601 }
602
604}
605
606
607
608////////////////////////////////////////////////////////////////////////////////
609/// Set appropriate cache operation mode for integral depending on cache operation
610/// mode of server objects
611
613{
614 // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
615 for (const auto server : _serverList) {
616 if (server->isValueServer(*this)) {
617 RooArgSet leafSet ;
618 server->leafNodeServerList(&leafSet) ;
619 for (const auto leaf : leafSet) {
620 if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
622 break ;
623 }
624 if (leaf->getAttribute("projectedDependent")) {
626 break ;
627 }
628 }
629 }
630 }
631}
632
633
634
635////////////////////////////////////////////////////////////////////////////////
636/// Utility function that returns true if 'object server' is a server
637/// to exactly one of the RooAbsArgs in 'exclLVBranches'
638
639bool RooRealIntegral::servesExclusively(const RooAbsArg* server,const RooArgSet& exclLVBranches, const RooArgSet& allBranches) const
640{
641 // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
642
643 // Special case, no LV servers available
644 if (exclLVBranches.empty()) return false ;
645
646 // If server has no clients and is not an LValue itself, return false
647 if (server->_clientList.empty() && exclLVBranches.find(server->GetName())) {
648 return false ;
649 }
650
651 // WVE must check for value relations only here!!!!
652
653 // Loop over all clients
654 Int_t numLVServ(0) ;
655 for (const auto client : server->valueClients()) {
656 // If client is not an LValue, recurse
657 if (!(exclLVBranches.find(client->GetName())==client)) {
658 if (allBranches.find(client->GetName())==client) {
659 if (!servesExclusively(client,exclLVBranches,allBranches)) {
660 // Client is a non-LValue that doesn't have an exclusive LValue server
661 return false ;
662 }
663 }
664 } else {
665 // Client is an LValue
666 numLVServ++ ;
667 }
668 }
669
670 return (numLVServ==1) ;
671}
672
673
674
675
676////////////////////////////////////////////////////////////////////////////////
677/// (Re)Initialize numerical integration engine if necessary. Return true if
678/// successful, or otherwise false.
679
681{
682 // if we already have an engine, check if it still works for the present limits.
683 if(_numIntEngine) {
684 if(_numIntEngine->isValid() && _numIntEngine->checkLimits() && !_restartNumIntEngine ) return true;
685 // otherwise, cleanup the old engine
686 _numIntEngine.reset();
687 _numIntegrand.reset();
688 }
689
690 // All done if there are no arguments to integrate numerically
691 if(_intList.empty()) return true;
692
693 // Bind the appropriate analytic integral of our RooRealVar object to
694 // those of its arguments that will be integrated out numerically.
695 if(_mode != 0) {
696 std::unique_ptr<RooAbsReal> analyticalPart{_function->createIntegral(_anaList, *funcNormSet(), RooNameReg::str(_rangeName))};
697 _numIntegrand = std::make_unique<RooRealBinding>(*analyticalPart,_intList,nullptr,false,_rangeName);
698 const_cast<RooRealIntegral*>(this)->addOwnedComponents(std::move(analyticalPart));
699 }
700 else {
701 _numIntegrand = std::make_unique<RooRealBinding>(*_function,_intList,funcNormSet(),false,_rangeName);
702 }
703 if(nullptr == _numIntegrand || !_numIntegrand->isValid()) {
704 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << std::endl;
705 return false;
706 }
707
708 // Create appropriate numeric integrator using factory
709 bool isBinned = _function->isBinnedDistribution(_intList) ;
711
712 if(_numIntEngine == nullptr || !_numIntEngine->isValid()) {
713 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << std::endl;
714 return false;
715 }
716
717 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
718 << _numIntEngine->ClassName() << " to calculate Int" << _intList << std::endl ;
719
720 if (_intList.size()>3) {
721 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") evaluation requires " << _intList.size() << "-D numeric integration step. Evaluation may be slow, sufficient numeric precision for fitting & minimization is not guaranteed" << std::endl ;
722 }
723
724 _restartNumIntEngine = false ;
725 return true;
726}
727
728
729
730////////////////////////////////////////////////////////////////////////////////
731/// Copy constructor
732
734 RooAbsReal(other,name),
735 _valid(other._valid),
736 _respectCompSelect(other._respectCompSelect),
737 _sumList("!sumList",this,other._sumList),
738 _intList("!intList",this,other._intList),
739 _anaList("!anaList",this,other._anaList),
740 _jacList("!jacList",this,other._jacList),
741 _facList("!facList","Variables independent of function",this,false,true),
742 _function("!func",this,other._function),
743 _iconfig(other._iconfig),
744 _sumCat("!sumCat",this,other._sumCat),
745 _mode(other._mode),
746 _intOperMode(other._intOperMode),
747 _restartNumIntEngine(false),
748 _rangeName(other._rangeName),
749 _cacheNum(false)
750{
751 if(other._funcNormSet) {
752 _funcNormSet = std::make_unique<RooArgSet>();
753 other._funcNormSet->snapshot(*_funcNormSet, false);
754 }
755
756 for (const auto arg : other._facList) {
757 std::unique_ptr<RooAbsArg> argClone{static_cast<RooAbsArg*>(arg->Clone())};
758 _facList.add(*argClone) ;
759 addServer(*argClone,false,true) ;
760 _facListOwned.addOwned(std::move(argClone));
761 }
762
763 other._intList.snapshot(_saveInt) ;
764 other._sumList.snapshot(_saveSum) ;
765
767}
768
769
770////////////////////////////////////////////////////////////////////////////////
771
773{
775}
776
777
778////////////////////////////////////////////////////////////////////////////////
779
780RooFit::OwningPtr<RooAbsReal> RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const
781{
782 // Handle special case of no integration with default algorithm
783 if (iset.empty()) {
784 return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
785 }
786
787 // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
788 RooArgSet isetAll(iset) ;
789 isetAll.add(_sumList) ;
790 isetAll.add(_intList) ;
791 isetAll.add(_anaList) ;
792 isetAll.add(_facList) ;
793
794 const RooArgSet* newNormSet(nullptr) ;
795 std::unique_ptr<RooArgSet> tmp;
796 if (nset && !_funcNormSet) {
797 newNormSet = nset ;
798 } else if (!nset && _funcNormSet) {
799 newNormSet = _funcNormSet.get();
800 } else if (nset && _funcNormSet) {
801 tmp = std::make_unique<RooArgSet>();
802 tmp->add(*nset) ;
803 tmp->add(*_funcNormSet,true) ;
804 newNormSet = tmp.get();
805 }
806 return _function->createIntegral(isetAll,newNormSet,cfg,rangeName);
807}
808
809
810
811////////////////////////////////////////////////////////////////////////////////
812/// Return value of object. If the cache is clean, return the
813/// cached value, otherwise recalculate on the fly and refill
814/// the cache
815
816double RooRealIntegral::getValV(const RooArgSet* nset) const
817{
818// // fast-track clean-cache processing
819// if (_operMode==AClean) {
820// return _value ;
821// }
822
823 if (nset && nset->uniqueId().value() != _lastNormSetId) {
824 const_cast<RooRealIntegral*>(this)->setProxyNormSet(nset);
825 _lastNormSetId = nset->uniqueId().value();
826 }
827
829 _value = traceEval(nset) ;
830 }
831
832 return _value ;
833}
834
835
836
837
838
839////////////////////////////////////////////////////////////////////////////////
840/// Perform the integration and return the result
841
843{
845
846 double retVal(0) ;
847 switch (_intOperMode) {
848
849 case Hybrid:
850 {
851 // Cache numeric integrals in >1d expensive object cache
852 RooDouble* cacheVal(nullptr) ;
855 }
856
857 if (cacheVal) {
858 retVal = *cacheVal ;
859 // cout << "using cached value of integral" << GetName() << std::endl ;
860 } else {
861
862
863 // Find any function dependents that are AClean
864 // and switch them temporarily to ADirty
865 bool origState = inhibitDirty() ;
866 setDirtyInhibit(true) ;
867
868 // try to initialize our numerical integration engine
869 if(!(_valid= initNumIntegrator())) {
870 coutE(Integration) << ClassName() << "::" << GetName()
871 << ":evaluate: cannot initialize numerical integrator" << std::endl;
872 return 0;
873 }
874
875 // Save current integral dependent values
878
879 // Evaluate sum/integral
880 retVal = sum() ;
881
882
883 // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
884 setDirtyInhibit(origState) ;
885
886 // Restore integral dependent values
889
890 // Cache numeric integrals in >1d expensive object cache
892 RooDouble* val = new RooDouble(retVal) ;
894 // cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << std::endl ;
895 }
896
897 }
898 break ;
899 }
900 case Analytic:
901 {
903 cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
904 << ")func = " << _function->ClassName() << "::" << _function->GetName()
905 << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << std::endl ;
906
907
908 break ;
909 }
910
911 case PassThrough:
912 {
913 // In pass through mode, the RooRealIntegral should have registered the
914 // function as a value server, because we directly depend on its value.
915 assert(_function.isValueServer());
916 // There should be no other servers besides the actual function and the
917 // factorized observables that the function doesn't depend on but are
918 // integrated over later.
919 assert(servers().size() == _facList.size() + 1);
920
921 //setDirtyInhibit(true) ;
922 retVal= _function->getVal(funcNormSet());
923 //setDirtyInhibit(false) ;
924 break ;
925 }
926 }
927
928
929 // Multiply answer with integration ranges of factorized variables
930 if (!_facList.empty()) {
931 for (const auto arg : _facList) {
932 // Multiply by fit range for 'real' dependents
933 if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
934 RooAbsRealLValue* argLV = (RooAbsRealLValue*)arg ;
935 retVal *= (argLV->getMax(intRange()) - argLV->getMin(intRange())) ;
936 }
937 // Multiply by number of states for category dependents
938 if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
940 retVal *= argLV->numTypes() ;
941 }
942 }
943 }
944
945
946 if (dologD(Tracing)) {
947 cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
948 switch(_intOperMode) {
949 case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
950 case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
951 case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
952 }
953
954 ccxcoutD(Tracing) << "raw*fact = " << retVal << std::endl ;
955 }
956
957 return retVal ;
958}
959
960
961
962////////////////////////////////////////////////////////////////////////////////
963/// Return product of jacobian terms originating from analytical integration
964
966{
967 if (_jacList.empty()) {
968 return 1 ;
969 }
970
971 double jacProd(1) ;
972 for (const auto elm : _jacList) {
973 auto arg = static_cast<const RooAbsRealLValue*>(elm);
974 jacProd *= arg->jacobian() ;
975 }
976
977 // Take std::abs() here: if jacobian is negative, min and max are swapped and analytical integral
978 // will be positive, so must multiply with positive jacobian.
979 return std::abs(jacProd) ;
980}
981
982
983
984////////////////////////////////////////////////////////////////////////////////
985/// Perform summation of list of category dependents to be integrated
986
988{
989 if (!_sumList.empty()) {
990 // Add integrals for all permutations of categories summed over
991 double total(0) ;
992
994 for (const auto& nameIdx : *sumCat) {
995 sumCat->setIndex(nameIdx);
996 if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
998 }
999 }
1000
1001 return total ;
1002
1003 } else {
1004 // Simply return integral
1005 double ret = integrate() / jacobianProduct() ;
1006 return ret ;
1007 }
1008}
1009
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Perform hybrid numerical/analytical integration over all real-valued dependents
1013
1015{
1016 if (!_numIntEngine) {
1017 // Trivial case, fully analytical integration
1019 } else {
1020 return _numIntEngine->calculate() ;
1021 }
1022}
1023
1024
1025
1026////////////////////////////////////////////////////////////////////////////////
1027/// Intercept server redirects and reconfigure internal object accordingly
1028
1030 bool mustReplaceAll, bool nameChange, bool isRecursive)
1031{
1032 _restartNumIntEngine = true ;
1033
1035
1036 // Update contents value caches for _intList and _sumList
1041
1042 // Delete parameters cache if we have one
1043 _params.reset();
1044
1045 return RooAbsReal::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
1046}
1047
1048
1049
1050////////////////////////////////////////////////////////////////////////////////
1051
1053{
1054 if (!_params) {
1055 _params = std::make_unique<RooArgSet>("params") ;
1056
1057 RooArgSet params ;
1058 for (const auto server : _serverList) {
1059 if (server->isValueServer(*this)) _params->add(*server) ;
1060 }
1061 }
1062
1063 return *_params ;
1064}
1065
1066
1067////////////////////////////////////////////////////////////////////////////////
1068/// Check if current value is valid
1069
1070bool RooRealIntegral::isValidReal(double /*value*/, bool /*printError*/) const
1071{
1072 return true ;
1073}
1074
1075////////////////////////////////////////////////////////////////////////////////
1076/// Check if component selection is allowed
1077
1079 return _respectCompSelect;
1080}
1081
1082////////////////////////////////////////////////////////////////////////////////
1083/// Set component selection to be allowed/forbidden
1084
1086 _respectCompSelect = allow;
1087}
1088
1090{
1091 if (!_sumList.empty() || !_intList.empty()) {
1092 std::stringstream errorMsg;
1093 errorMsg << "Only analytical integrals are supported for AD for class" << _function.GetName();
1094 coutE(Minimization) << errorMsg.str() << std::endl;
1095 throw std::runtime_error(errorMsg.str().c_str());
1096 }
1097
1099}
1100
1101////////////////////////////////////////////////////////////////////////////////
1102/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1103/// integration operation
1104
1105void RooRealIntegral::printMetaArgs(ostream& os) const
1106{
1107
1108 if (!intVars().empty()) {
1109 os << "Int " ;
1110 }
1111 os << _function->GetName() ;
1112 if (_funcNormSet) {
1113 os << "_Norm" << *_funcNormSet << " " ;
1114 }
1115
1116 // List internally integrated observables and factorizing observables as analytically integrated
1117 RooArgSet tmp(_anaList) ;
1118 tmp.add(_facList) ;
1119 if (!tmp.empty()) {
1120 os << "d[Ana]" << tmp << " ";
1121 }
1122
1123 // List numerically integrated and summed observables as numerically integrated
1124 RooArgSet tmp2(_intList) ;
1125 tmp2.add(_sumList) ;
1126 if (!tmp2.empty()) {
1127 os << " d[Num]" << tmp2 << " ";
1128 }
1129}
1130
1131
1132
1133////////////////////////////////////////////////////////////////////////////////
1134/// Print the state of this object to the specified output stream.
1135
1136void RooRealIntegral::printMultiline(ostream& os, Int_t contents, bool verbose, TString indent) const
1137{
1138 RooAbsReal::printMultiline(os,contents,verbose,indent) ;
1139 os << indent << "--- RooRealIntegral ---" << std::endl;
1140 os << indent << " Integrates ";
1142 TString deeper(indent);
1143 deeper.Append(" ");
1144 os << indent << " operating mode is "
1145 << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << std::endl ;
1146 os << indent << " Summed discrete args are " << _sumList << std::endl ;
1147 os << indent << " Numerically integrated args are " << _intList << std::endl;
1148 os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << std::endl ;
1149 os << indent << " Arguments included in Jacobian are " << _jacList << std::endl ;
1150 os << indent << " Factorized arguments are " << _facList << std::endl ;
1151 os << indent << " Function normalization set " ;
1152 if (_funcNormSet)
1153 _funcNormSet->Print("1") ;
1154 else
1155 os << "<none>" ;
1156
1157 os << std::endl ;
1158}
1159
1160
1161
1162////////////////////////////////////////////////////////////////////////////////
1163/// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1164
1166 _cacheAllNDim = ndim ;
1167}
1168
1169
1170////////////////////////////////////////////////////////////////////////////////
1171/// Return minimum dimensions of numeric integration for which values are cached.
1172
1174{
1175 return _cacheAllNDim ;
1176}
1177
1178
1179std::unique_ptr<RooAbsArg>
1181{
1183}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#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
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
static unsigned int total
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2468
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:80
RooExpensiveObjectCache & expensiveObjectCache() const
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
void setOperMode(OperMode mode, bool recurseADirty=true)
Set the operation mode of this node.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
const RefCountList_t & valueClients() const
List of all value clients of this object. Value clients receive value updates.
Definition RooAbsArg.h:199
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition RooAbsArg.h:504
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
static void setDirtyInhibit(bool flag)
Control global dirty inhibit mode.
virtual std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:209
bool dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr) const
Check whether this object depends on values from an element in the serverList.
Definition RooAbsArg.h:109
void addServer(RooAbsArg &server, bool valueProp=true, bool shapeProp=false, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
bool isValueOrShapeDirtyAndClear() const
Definition RooAbsArg.h:456
bool inhibitDirty() const
Delete watch flag.
void setProxyNormSet(const RooArgSet *nset)
Forward a change in the cached normalization argset to all the registered proxies.
RefCountList_t _clientList
Definition RooAbsArg.h:636
RefCountList_t _serverList
Definition RooAbsArg.h:635
bool isValueServer(const RooAbsArg &arg) const
Check if this is serving values to arg.
Definition RooAbsArg.h:226
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:485
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
static TClass * Class()
Int_t numTypes(const char *=nullptr) const
Return number of types defined (in range named rangeName if rangeName!=nullptr)
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
RooFit::UniqueId< RooAbsCollection > const & uniqueId() const
Returns a unique ID that is different for every instantiated RooAbsCollection.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
Storage_t const & get() const
Const access to the underlying stl container.
void sortTopologically()
Sort collection topologically: the servers of any RooAbsArg will be before that RooAbsArg in the coll...
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
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...
Storage_t::size_type size() const
RooAbsArg * first() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
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 bool isJacobianOK(const RooArgSet &depList) const
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
static TClass * Class()
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:59
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Structure printing.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Function that is called at the end of redirectServers().
double _value
Cache for current value of object.
Definition RooAbsReal.h:543
double traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
virtual std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const
This function defines the analytical integral translation for the class.
RooFit::UniqueId< RooArgSet >::Value_t _lastNormSetId
Component selection flag for RooAbsPdf::plotCompOn.
Definition RooAbsReal.h:550
virtual double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
virtual bool isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:353
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
bool _valueServer
If true contents is value server of owner.
Definition RooArgProxy.h:81
bool isValueServer() const
Returns true of contents is value server of owner.
Definition RooArgProxy.h:61
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
bool addOwned(RooAbsArg &var, bool silent=false) override
Overloaded RooCollection_t::addOwned() method insert object into owning set and registers object as s...
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
static TClass * Class()
RooDouble is a minimal implementation of a TObject holding a double value.
Definition RooDouble.h:22
static TClass * Class()
bool registerObject(const char *ownerName, const char *objectName, TObject &cacheObject, const RooArgSet &params)
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,...
A class to maintain the context for squashing of RooFit models into code.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
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.
Definition RooNameReg.h:37
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
std::unique_ptr< RooAbsIntegrator > createIntegrator(RooAbsFunc &func, const RooNumIntConfig &config, Int_t ndim=0, bool isBinned=false) 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,...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooNumIntConfig * _iconfig
bool initNumIntegrator() const
(Re)Initialize numerical integration engine if necessary.
RooArgSet const * funcNormSet() const
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooArgSet *nset=nullptr, const RooNumIntConfig *cfg=nullptr, const char *rangeName=nullptr) const override
Create an object that represents the integral of the function over one or more observables listed in ...
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
void setAllowComponentSelection(bool allow)
Set component selection to be allowed/forbidden.
RooRealProxy _function
Function being integration.
RooArgSet intVars() const
RooSetProxy _intList
Set of continuous observables over which is integrated numerically.
virtual double sum() const
Perform summation of list of category dependents to be integrated.
RooSetProxy _facList
Set of observables on which function does not depends, which are integrated nevertheless.
std::unique_ptr< RooArgSet > _params
! cache for set of parameters
static void setCacheAllNumeric(Int_t ndim)
Global switch to cache all integral values that integrate at least ndim dimensions numerically.
IntOperMode _intOperMode
integration operation mode
bool _cacheNum
Cache integral if numeric.
double evaluate() const override
Perform the integration and return the result.
const RooArgSet & parameters() const
std::unique_ptr< RooAbsFunc > _numIntegrand
! do not persist
RooSetProxy _jacList
Set of lvalue observables over which is analytically integration that have a non-unit Jacobian.
bool isValidReal(double value, bool printError=false) const override
Check if current value is valid.
RooArgSet _saveInt
! do not persist
double getValV(const RooArgSet *set=nullptr) const override
Return value of object.
RooSetProxy _anaList
Set of observables over which is integrated/summed analytically.
bool _restartNumIntEngine
! do not persist
bool 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...
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override
Intercept server redirects and reconfigure internal object accordingly.
RooSetProxy _sumList
Set of discrete observable over which is summed numerically.
~RooRealIntegral() override
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the...
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print the state of this object to the specified output stream.
std::unique_ptr< RooAbsIntegrator > _numIntEngine
! do not persist
virtual double integrate() const
Perform hybrid numerical/analytical integration over all real-valued dependents.
RooListProxy _sumCat
! do not persist
virtual double jacobianProduct() const
Return product of jacobian terms originating from analytical integration.
static Int_t getCacheAllNumeric()
Return minimum dimensions of numeric integration for which values are cached.
RooArgSet _saveSum
! do not persist
static Int_t _cacheAllNDim
! Cache all integrals with given numeric dimension
std::unique_ptr< RooArgSet > _funcNormSet
Optional normalization set passed to function.
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
void autoSelectDirtyMode()
Set appropriate cache operation mode for integral depending on cache operation mode of server objects...
RooArgSet _facListOwned
Owned components in _facList.
const char * intRange() const
bool getAllowComponentSelection() const
Check if component selection is allowed.
The RooSuperCategory can join several RooAbsCategoryLValue objects into a single category.
bool setIndex(value_type index, bool printError=true) override
Set the value of the super category to the specified index.
bool 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
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:576
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition RExports.h:167
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:43
constexpr Value_t value() const
Return numerical value of ID.
Definition UniqueId.h:59