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
22Performs 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
34#include <RooAbsRealLValue.h>
35#include <RooConstVar.h>
36#include <RooDouble.h>
38#include <RooInvTransform.h>
39#include <RooMsgService.h>
40#include <RooNameReg.h>
41#include <RooNumIntConfig.h>
42#include <RooNumIntFactory.h>
43#include <RooRealBinding.h>
44#include <RooSuperCategory.h>
45#include <RooTrace.h>
46
47#include <iostream>
48#include <memory>
49
50
51namespace {
52
53/// Utility function that returns true if 'object server' is a server
54/// to exactly one of the RooAbsArgs in 'exclLVBranches'
56{
57 // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
58
59 // Special case, no LV servers available
60 if (exclLVBranches.empty())
61 return false;
62
63 // If server has no clients and is not an LValue itself, return false
64 if (!server->hasClients() && exclLVBranches.find(server->GetName())) {
65 return false;
66 }
67
68 // WVE must check for value relations only here!!!!
69
70 // Loop over all clients
72 for (const auto client : server->valueClients()) {
73 // If client is not an LValue, recurse
74 if (!(exclLVBranches.find(client->GetName()) == client)) {
75 if (allBranches.find(client->GetName()) == client) {
77 // Client is a non-LValue that doesn't have an exclusive LValue server
78 return false;
79 }
80 }
81 } else {
82 // Client is an LValue
83 numLVServ++;
84 }
85 }
86
87 return (numLVServ == 1);
88}
89
90struct ServerToAdd {
91 ServerToAdd(RooAbsArg *theArg, bool isShape) : arg{theArg}, isShapeServer{isShape} {}
92 RooAbsArg *arg = nullptr;
93 bool isShapeServer = false;
94};
95
96void addObservableToServers(RooAbsReal const &function, RooAbsArg &leaf, std::vector<ServerToAdd> &serversToAdd,
97 const char *rangeName)
98{
99 auto leaflv = dynamic_cast<RooAbsRealLValue *>(&leaf);
100 if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
101 oocxcoutD(&function, Integration)
102 << function.GetName() << " : Observable " << leaf.GetName()
103 << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf"
104 << std::endl;
105 if (leaflv->getBinning(rangeName).lowBoundFunc()) {
106 serversToAdd.emplace_back(leaflv->getBinning(rangeName).lowBoundFunc(), false);
107 }
108 if (leaflv->getBinning(rangeName).highBoundFunc()) {
109 serversToAdd.emplace_back(leaflv->getBinning(rangeName).highBoundFunc(), false);
110 }
111 } else {
112 oocxcoutD(&function, Integration) << function.GetName() << ": Adding observable " << leaf.GetName()
113 << " as shape dependent" << std::endl;
114 serversToAdd.emplace_back(&leaf, true);
115 }
116}
117
118void addParameterToServers(RooAbsReal const &function, RooAbsArg &leaf, std::vector<ServerToAdd> &serversToAdd,
119 bool isShapeServer)
120{
121 if (!isShapeServer) {
122 oocxcoutD(&function, Integration) << function.GetName() << ": Adding parameter " << leaf.GetName()
123 << " as value dependent" << std::endl;
124 } else {
125 oocxcoutD(&function, Integration) << function.GetName() << ": Adding parameter " << leaf.GetName()
126 << " as shape dependent" << std::endl;
127 }
128 serversToAdd.emplace_back(&leaf, isShapeServer);
129}
130
131enum class MarkedState { Dependent, Independent, AlreadyAdded };
132
133/// Mark all args that recursively are value clients of "dep".
134void unmarkDepValueClients(RooAbsArg const &dep, RooArgSet const &args, std::vector<MarkedState> &marked)
135{
136 assert(args.size() == marked.size());
137 auto index = args.index(dep);
138 if (index >= 0) {
139 marked[index] = MarkedState::Dependent;
140 for (RooAbsArg *client : dep.valueClients()) {
141 unmarkDepValueClients(*client, args, marked);
142 }
143 }
144}
145
146std::vector<ServerToAdd>
147getValueAndShapeServers(RooAbsReal const &function, RooArgSet const &depList, const char *rangeName)
148{
149 std::vector<ServerToAdd> serversToAdd;
150
151 // Get the full computation graph and sort it topologically
153 function.treeNodeServerList(&allArgsList, nullptr, true, true, /*valueOnly=*/false, false);
155 allArgs.sortTopologically();
156
157 // Figure out what are all the value servers only
159 function.treeNodeServerList(&allValueArgsList, nullptr, true, true, /*valueOnly=*/true, false);
161
162 // All "marked" args will be added as value servers to the integral
163 std::vector<MarkedState> marked(allArgs.size(), MarkedState::Independent);
164 marked.back() = MarkedState::Dependent; // We don't want to consider the function itself
165
166 // Mark all args that are (indirect) value servers of the integration
167 // variable or the integration variable itself. If something was marked,
168 // it means the integration variable was in the compute graph and we will
169 // add it to the server list.
170 for (RooAbsArg *dep : depList) {
171 if (RooAbsArg *depInArgs = allArgs.find(dep->GetName())) {
174 }
175 }
176
177 // We are adding all independent direct servers of the args depending on the
178 // integration variables
179 for (std::size_t i = 0; i < allArgs.size(); ++i) {
180 if (marked[i] == MarkedState::Dependent) {
181 for (RooAbsArg *server : allArgs[i]->servers()) {
182 int index = allArgs.index(server->GetName());
183 if (index >= 0 && marked[index] == MarkedState::Independent) {
185 marked[index] = MarkedState::AlreadyAdded;
186 }
187 }
188 }
189 }
190
191 return serversToAdd;
192}
193
195 const RooArgSet &allBranches)
196{
197 // If any of the branches in the computation graph of the function depend on
198 // the integrated variable, we can't do analytical integration. The only
199 // case where this would work is if the branch is an l-value with known
200 // Jacobian, but this case is already handled in step B) in the constructor
201 // by reexpressing the original integration variables in terms of
202 // higher-order l-values if possible.
204 for (RooAbsArg *intDep : intDeps) {
205 bool depOK = true;
206 for (RooAbsArg *branch : allBranches) {
207 // It's ok if the branch is the integration variable itself
208 if (intDep->namePtr() != branch->namePtr() && branch->dependsOnValue(*intDep)) {
209 depOK = false;
210 }
211 if (!depOK) break;
212 }
213 if (depOK) {
215 }
216 }
217
218 for (const auto arg : function.servers()) {
219
220 // Dependent or parameter?
221 if (!arg->dependsOnValue(filteredIntDeps)) {
222 continue;
223 } else if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
224 // Skip arg if it is neither value or shape server
225 continue;
226 }
227
228 bool depOK(false);
229 // Check for integratable AbsRealLValue
230
231 if (arg->isDerived()) {
232 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue *>(arg);
233 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue *>(arg);
234 if ((realArgLV && filteredIntDeps.find(realArgLV->GetName()) &&
235 (realArgLV->isJacobianOK(filteredIntDeps) != 0)) ||
236 catArgLV) {
237
238 // Derived LValue with valid jacobian
239 depOK = true;
240
241 // Now, check for overlaps
242 bool overlapOK = true;
243 for (const auto otherArg : function.servers()) {
244 // skip comparison with self
245 if (arg == otherArg)
246 continue;
247 if (dynamic_cast<RooConstVar const *>(otherArg))
248 continue;
249 if (arg->overlaps(*otherArg, true)) {
250 }
251 }
252 // coverity[DEADCODE]
253 if (!overlapOK)
254 depOK = false;
255 }
256 } else {
257 // Fundamental types are always OK
258 depOK = true;
259 }
260
261 // Add server to list of dependents that are OK for analytical integration
262 if (depOK) {
263 anIntOKDepList.add(*arg, true);
264 oocxcoutI(&function, Integration)
265 << function.GetName() << ": Observable " << arg->GetName()
266 << " is suitable for analytical integration (if supported by p.d.f)" << std::endl;
267 }
268 }
269}
270
271} // namespace
272
274
275////////////////////////////////////////////////////////////////////////////////
276
281
282////////////////////////////////////////////////////////////////////////////////
283/// Construct integral of 'function' over observables in 'depList'
284/// in range 'rangeName' with normalization observables 'funcNormSet'
285/// (for p.d.f.s). In the integral is performed to the maximum extent
286/// possible the internal (analytical) integrals advertised by function.
287/// The other integrations are performed numerically. The optional
288/// config object prescribes how these numeric integrations are configured.
289///
290/// \Note If pdf component selection was globally overridden to always include
291/// all components (either with RooAbsReal::globalSelectComp(bool) or a
292/// RooAbsReal::GlobalSelectComponentRAII), then any created integral will
293/// ignore component selections during its lifetime. This is especially useful
294/// when creating normalization or projection integrals.
295RooRealIntegral::RooRealIntegral(const char *name, const char *title,
296 const RooAbsReal& function, const RooArgSet& depList,
297 const RooArgSet* funcNormSet, const RooNumIntConfig* config,
298 const char* rangeName) :
299 RooAbsReal(name,title),
300 _valid(true),
301 _respectCompSelect{!_globalSelectComp},
302 _sumList("!sumList","Categories to be summed numerically",this,false,false),
303 _intList("!intList","Variables to be integrated numerically",this,false,false),
304 _anaList("!anaList","Variables to be integrated analytically",this,false,false),
305 _jacList("!jacList","Jacobian product term",this,false,false),
306 _facList("!facList","Variables independent of function",this,false,true),
307 _function("!func","Function to be integrated",this,false,false),
308 _iconfig(const_cast<RooNumIntConfig*>(config)),
309 _sumCat("!sumCat","SuperCategory for summation",this,false,false),
310 _rangeName(const_cast<TNamed*>(RooNameReg::ptr(rangeName)))
311{
312 // A) Check that all dependents are lvalues
313 //
314 // B) Check if list of dependents can be re-expressed in
315 // lvalues that are higher in the expression tree
316 //
317 // C) Check for dependents that the PDF insists on integrating
318 // analytically itself
319 //
320 // D) Make list of servers that can be integrated analytically
321 // Add all parameters/dependents as value/shape servers
322 //
323 // E) Interact with function to make list of objects actually integrated analytically
324 //
325 // F) Make list of numerical integration variables consisting of:
326 // - Category dependents of RealLValues in analytical integration
327 // - Leaf nodes server lists of function server that are not analytically integrated
328 // - Make Jacobian list for analytically integrated RealLValues
329 //
330 // G) Split numeric list in integration list and summation list
331 //
332
333 oocxcoutI(&function,Integration) << "RooRealIntegral::ctor(" << GetName() << ") Constructing integral of function "
334 << function.GetName() << " over observables" << depList << " with normalization "
335 << (funcNormSet?*funcNormSet:RooArgSet()) << " with range identifier "
336 << (rangeName?rangeName:"<none>") << std::endl ;
337
338
339 // Choose same expensive object cache as integrand
340 setExpensiveObjectCache(function.expensiveObjectCache()) ;
341// std::cout << "RRI::ctor(" << GetName() << ") setting expensive object cache to " << &expensiveObjectCache() << " as taken from " << function.GetName() << std::endl ;
342
343 // Use objects integrator configuration if none is specified
344 if (!_iconfig) _iconfig = const_cast<RooNumIntConfig*>(function.getIntegratorConfig());
345
346 // Save private copy of funcNormSet, if supplied, excluding factorizing terms
347 if (funcNormSet) {
348 _funcNormSet = std::make_unique<RooArgSet>();
349 for (const auto nArg : *funcNormSet) {
350 if (function.dependsOn(*nArg)) {
351 _funcNormSet->addClone(*nArg) ;
352 }
353 }
354 }
355
356 //_funcNormSet = funcNormSet ? (RooArgSet*)funcNormSet->snapshot(false) : 0 ;
357
358 // Make internal copy of dependent list
360
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
362 // * A) Check that all dependents are lvalues and filter out any
363 // dependents that the PDF doesn't explicitly depend on
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
365
366 for (auto arg : intDepList) {
367 if(!arg->isLValue()) {
368 coutE(InputArguments) << ClassName() << "::" << GetName() << ": cannot integrate non-lvalue ";
369 arg->Print("1");
370 _valid= false;
371 }
372 if (!function.dependsOn(*arg)) {
373 std::unique_ptr<RooAbsArg> argClone{static_cast<RooAbsArg*>(arg->Clone())};
375 addOwnedComponents(std::move(argClone));
376 }
377 }
378
379 if (!_facList.empty()) {
380 oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables are " << _facList << std::endl ;
381 }
382
383
384 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
385 // * B) Check if list of dependents can be re-expressed in *
386 // * lvalues that are higher in the expression tree *
387 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
388
389
390 // Initial fill of list of LValue branches
391 RooArgSet exclLVBranches("exclLVBranches") ;
394 function.branchNodeServerList(&branchList) ;
395
396 for (auto branch: branchList) {
399 if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
400 exclLVBranches.add(*branch) ;
401 }
402 if (branch != &function && function.dependsOnValue(*branch)) {
403 branchListVD.add(*branch) ;
404 }
405 }
406 exclLVBranches.remove(depList,true,true) ;
407
408 // Initial fill of list of LValue leaf servers (put in intDepList, but the
409 // instances that are in the actual computation graph of the function)
410 RooArgSet exclLVServers("exclLVServers") ;
411 function.getObservables(&intDepList, exclLVServers);
412
413 // Obtain mutual exclusive dependence by iterative reduction
414 bool converged(false) ;
415 while(!converged) {
417
418 // Reduce exclLVServers to only those serving exclusively exclLVBranches
419 std::vector<RooAbsArg*> toBeRemoved;
420 for (auto server : exclLVServers) {
422 toBeRemoved.push_back(server);
424 }
425 }
427
428 // Reduce exclLVBranches to only those depending exclusively on exclLVservers
429 // Attention: counting loop, since erasing from container
430 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
431 const RooAbsArg* branch = exclLVBranches[i];
433 branch->getObservables(&intDepList, brDepList);
434 RooArgSet bsList(brDepList,"bsList") ;
435 bsList.remove(exclLVServers,true,true) ;
436 if (!bsList.empty()) {
437 exclLVBranches.remove(*branch,true,true) ;
438 --i;
440 }
441 }
442 }
443
444 // Eliminate exclLVBranches that do not depend on any LVServer
445 // Attention: Counting loop, since modifying container
446 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
447 const RooAbsArg* branch = exclLVBranches[i];
448 if (!branch->dependsOnValue(exclLVServers)) {
449 exclLVBranches.remove(*branch,true,true) ;
450 --i;
451 }
452 }
453
454 // Replace exclusive lvalue branch servers with lvalue branches
455 // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
456 if (!exclLVServers.empty() && !function.isBinnedDistribution(exclLVBranches)) {
457 intDepList.remove(exclLVServers) ;
459 }
460
461
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
463 // * C) Check for dependents that the PDF insists on integrating *
464 // analytically itself *
465 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
466
468 for (auto arg : intDepList) {
469 if (function.forceAnalyticalInt(*arg)) {
470 anIntOKDepList.add(*arg) ;
471 }
472 }
473
474 if (!anIntOKDepList.empty()) {
475 oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << std::endl ;
476 }
477
478
479 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
480 // * D) Make list of servers that can be integrated analytically *
481 // Add all parameters/dependents as value/shape servers *
482 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
483
486 // We will not add the servers just now, because it makes only sense to add
487 // them once we have made sure that this integral is not operating in
488 // pass-through mode. It will be done at the end of this constructor.
489
490 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
491 // * E) interact with function to make list of objects actually integrated analytically *
492 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
493
494 _mode = function.getAnalyticalIntegralWN(anIntOKDepList,_anaList,_funcNormSet.get(),RooNameReg::str(_rangeName)) ;
495
496 // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaList
497 if (_mode==0) {
499 }
500
501 if (_mode!=0) {
502 oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << std::endl ;
503 }
504
505 // WVE kludge: synchronize dset for use in analyticalIntegral
506 // LM : I think this is needed only if _funcNormSet is not an empty set
507 if (_funcNormSet && !_funcNormSet->empty()) {
508 function.getVal(_funcNormSet.get()) ;
509 }
510
511 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
512 // * F) Make list of numerical integration variables consisting of: *
513 // * - Category dependents of RealLValues in analytical integration *
514 // * - Expanded server lists of server that are not analytically integrated *
515 // * Make Jacobian list with analytically integrated RealLValues *
516 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
517
518 // Loop over actually analytically integrated dependents
519 for (const auto arg : _anaList) {
520
521 // Process only derived RealLValues
522 if (dynamic_cast<RooAbsRealLValue const *>(arg) && arg->isDerived() && !arg->isFundamental()) {
523
524 // Add to list of Jacobians to calculate
525 _jacList.add(*arg) ;
526
527 // Add category dependent of LValueReal used in integration
528 std::unique_ptr<RooArgSet> argDepList{arg->getObservables(&intDepList)};
529 for (const auto argDep : *argDepList) {
530 if (dynamic_cast<RooAbsCategoryLValue const *>(argDep) && intDepList.contains(*argDep)) {
532 }
533 }
534 }
535 }
536
537
538 // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
539 if (_anaList.empty()) {
540 if (!exclLVServers.empty()) {
541 //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << std::endl ;
542 intDepList.remove(exclLVBranches) ;
544 }
545 }
546 //cout << "NUMINT intDepList = " << intDepList << std::endl ;
547
548 // Loop again over function servers to add remaining numeric integrations
549 for (const auto arg : function.servers()) {
550
551 // Process only servers that are not treated analytically
552 if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
553
554 // Process only derived RealLValues
555 if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
556 addNumIntDep(*arg) ;
557 } else {
558
559 // WVE this will only get the observables, but not l-value transformations
560 // Expand server in final dependents
561 auto argDeps = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
562
563 // Add final dependents, that are not forcibly integrated analytically,
564 // to numerical integration list
565 for (const auto dep : *argDeps) {
566 if (!_anaList.find(dep->GetName())) {
568 }
569 }
570 }
571 }
572 }
573
574 if (!_anaList.empty()) {
575 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << std::endl ;
576 }
577 if (!_intList.empty()) {
578 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << std::endl ;
579 }
580 if (!_sumList.empty()) {
581 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << std::endl ;
582 }
583
584
585 // Determine operating mode
586 if (!_intList.empty() || !_sumList.empty()) {
587 // Numerical and optional Analytical integration
589 } else if (!_anaList.empty()) {
590 // Purely analytical integration
592 } else {
593 // No integration performed, where the function is a direct value server
595 _function._valueServer = true;
596 }
597 // We are only setting the function proxy now that it's clear if it's a value
598 // server or not.
599 _function.setArg(const_cast<RooAbsReal&>(function));
600
601 // Determine auto-dirty status
603
604 // Create value caches for _intList and _sumList
607
608
609 if (!_sumList.empty()) {
610 _sumCat.addOwned(std::make_unique<RooSuperCategory>(Form("%s_sumCat",GetName()),"sumCat",_sumList));
611 }
612
613 // Only if we are not in pass-through mode we need to add the shape and value
614 // servers separately.
616 for(auto const& toAdd : serversToAdd) {
617 addServer(*toAdd.arg, !toAdd.isShapeServer, toAdd.isShapeServer);
618 }
619 }
620
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// Set appropriate cache operation mode for integral depending on cache operation
626/// mode of server objects
627
629{
630 // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
631 for (const auto server : _serverList) {
632 if (server->isValueServer(*this)) {
634 server->leafNodeServerList(&leafSet) ;
635 for (const auto leaf : leafSet) {
636 if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
638 break ;
639 }
640 if (leaf->getAttribute("projectedDependent")) {
642 break ;
643 }
644 }
645 }
646 }
647}
648
649////////////////////////////////////////////////////////////////////////////////
650/// (Re)Initialize numerical integration engine if necessary. Return true if
651/// successful, or otherwise false.
652
654{
655 // if we already have an engine, check if it still works for the present limits.
656 if(_numIntEngine) {
657 if(_numIntEngine->isValid() && _numIntEngine->checkLimits() && !_restartNumIntEngine ) return true;
658 // otherwise, cleanup the old engine
659 _numIntEngine.reset();
660 _numIntegrand.reset();
661 }
662
663 // All done if there are no arguments to integrate numerically
664 if(_intList.empty()) return true;
665
666 // Bind the appropriate analytic integral of our RooRealVar object to
667 // those of its arguments that will be integrated out numerically.
668 if(_mode != 0) {
670 _numIntegrand = std::make_unique<RooRealBinding>(*analyticalPart,_intList,nullptr,false,_rangeName);
671 const_cast<RooRealIntegral*>(this)->addOwnedComponents(std::move(analyticalPart));
672 }
673 else {
674 _numIntegrand = std::make_unique<RooRealBinding>(*_function,_intList,actualFuncNormSet(),false,_rangeName);
675 }
676 if(nullptr == _numIntegrand || !_numIntegrand->isValid()) {
677 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << std::endl;
678 return false;
679 }
680
681 // Create appropriate numeric integrator using factory
683 std::string integratorName = RooNumIntFactory::instance().getIntegratorName(*_numIntegrand,*_iconfig,0,isBinned);
685
686 if(_numIntEngine == nullptr || !_numIntEngine->isValid()) {
687 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << std::endl;
688 return false;
689 }
690
691 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
692 << integratorName << " to calculate Int" << _intList << std::endl ;
693
694 if (_intList.size()>3) {
695 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 ;
696 }
697
699 return true;
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Copy constructor
704
707 _valid(other._valid),
708 _respectCompSelect(other._respectCompSelect),
709 _sumList("!sumList", this, other._sumList),
710 _intList("!intList", this, other._intList),
711 _anaList("!anaList", this, other._anaList),
712 _jacList("!jacList", this, other._jacList),
713 _facList("!facList", this, other._facList),
714 _function("!func", this, other._function),
715 _iconfig(other._iconfig),
716 _sumCat("!sumCat", this, other._sumCat),
717 _mode(other._mode),
718 _intOperMode(other._intOperMode),
719 _rangeName(other._rangeName)
720{
721 if(other._funcNormSet) {
722 _funcNormSet = std::make_unique<RooArgSet>();
723 other._funcNormSet->snapshot(*_funcNormSet, false);
724 }
725
726 other._intList.snapshot(_saveInt) ;
727 other._sumList.snapshot(_saveSum) ;
728
730}
731
732////////////////////////////////////////////////////////////////////////////////
733
738
739////////////////////////////////////////////////////////////////////////////////
740
742{
743 // Handle special case of no integration with default algorithm
744 if (iset.empty()) {
745 return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
746 }
747
748 // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
750 isetAll.add(_sumList) ;
751 isetAll.add(_intList) ;
752 isetAll.add(_anaList) ;
753 isetAll.add(_facList) ;
754
755 const RooArgSet* newNormSet(nullptr) ;
756 std::unique_ptr<RooArgSet> tmp;
757 if (nset && !_funcNormSet) {
758 newNormSet = nset ;
759 } else if (!nset && _funcNormSet) {
760 newNormSet = _funcNormSet.get();
761 } else if (nset && _funcNormSet) {
762 tmp = std::make_unique<RooArgSet>();
763 tmp->add(*nset) ;
764 tmp->add(*_funcNormSet,true) ;
765 newNormSet = tmp.get();
766 }
768}
769
770////////////////////////////////////////////////////////////////////////////////
771/// Return value of object. If the cache is clean, return the
772/// cached value, otherwise recalculate on the fly and refill
773/// the cache
774
775double RooRealIntegral::getValV(const RooArgSet* nset) const
776{
777// // fast-track clean-cache processing
778// if (_operMode==AClean) {
779// return _value ;
780// }
781
782 if (nset && nset->uniqueId().value() != _lastNormSetId) {
783 const_cast<RooRealIntegral*>(this)->setProxyNormSet(nset);
784 _lastNormSetId = nset->uniqueId().value();
785 }
786
788 _value = traceEval(nset) ;
789 }
790
791 return _value ;
792}
793
794////////////////////////////////////////////////////////////////////////////////
795/// Perform the integration and return the result
796
798{
800
801 double retVal(0) ;
802 switch (_intOperMode) {
803
804 case Hybrid:
805 {
806 // Cache numeric integrals in >1d expensive object cache
807 RooDouble const* cacheVal(nullptr) ;
808 if ((_cacheNum && !_intList.empty()) || int(_intList.size())>=_cacheAllNDim) {
810 }
811
812 if (cacheVal) {
813 retVal = *cacheVal ;
814 // std::cout << "using cached value of integral" << GetName() << std::endl ;
815 } else {
816
817
818 // Find any function dependents that are AClean
819 // and switch them temporarily to ADirty
820 bool origState = inhibitDirty() ;
821 setDirtyInhibit(true) ;
822
823 // try to initialize our numerical integration engine
824 if(!(_valid= initNumIntegrator())) {
825 coutE(Integration) << ClassName() << "::" << GetName()
826 << ":evaluate: cannot initialize numerical integrator" << std::endl;
827 return 0;
828 }
829
830 // Save current integral dependent values
833
834 // Evaluate sum/integral
835 retVal = sum() ;
836
837
838 // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
840
841 // Restore integral dependent values
844
845 // Cache numeric integrals in >1d expensive object cache
846 if ((_cacheNum && !_intList.empty()) || int(_intList.size())>=_cacheAllNDim) {
847 RooDouble* val = new RooDouble(retVal) ;
849 // std::cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << std::endl ;
850 }
851
852 }
853 break ;
854 }
855 case Analytic:
856 {
858 cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
859 << ")func = " << _function->ClassName() << "::" << _function->GetName()
860 << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << std::endl ;
861
862
863 break ;
864 }
865
866 case PassThrough:
867 {
868 // In pass through mode, the RooRealIntegral should have registered the
869 // function as a value server, because we directly depend on its value.
871 // There should be no other servers besides the actual function and the
872 // factorized observables that the function doesn't depend on but are
873 // integrated over later.
874 assert(servers().size() == _facList.size() + 1);
875
876 //setDirtyInhibit(true) ;
878 //setDirtyInhibit(false) ;
879 break ;
880 }
881 }
882
883
884 // Multiply answer with integration ranges of factorized variables
885 for (const auto arg : _facList) {
886 // Multiply by fit range for 'real' dependents
887 if (auto argLV = dynamic_cast<RooAbsRealLValue *>(arg)) {
888 retVal *= (argLV->getMax(intRange()) - argLV->getMin(intRange())) ;
889 }
890 // Multiply by number of states for category dependents
891 if (auto argLV = dynamic_cast<RooAbsCategoryLValue *>(arg)) {
892 retVal *= argLV->numTypes() ;
893 }
894 }
895
896
897 if (dologD(Tracing)) {
898 cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
899 switch(_intOperMode) {
900 case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
901 case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
902 case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
903 }
904
905 ccxcoutD(Tracing) << "raw*fact = " << retVal << std::endl ;
906 }
907
908 return retVal ;
909}
910
911////////////////////////////////////////////////////////////////////////////////
912/// Return product of jacobian terms originating from analytical integration
913
915{
916 if (_jacList.empty()) {
917 return 1 ;
918 }
919
920 double jacProd(1) ;
921 for (const auto elm : _jacList) {
922 auto arg = static_cast<const RooAbsRealLValue*>(elm);
923 jacProd *= arg->jacobian() ;
924 }
925
926 // Take std::abs() here: if jacobian is negative, min and max are swapped and analytical integral
927 // will be positive, so must multiply with positive jacobian.
928 return std::abs(jacProd) ;
929}
930
931////////////////////////////////////////////////////////////////////////////////
932/// Perform summation of list of category dependents to be integrated
933
935{
936 if (!_sumList.empty()) {
937 // Add integrals for all permutations of categories summed over
938 double total(0) ;
939
941 for (const auto& nameIdx : *sumCat) {
942 sumCat->setIndex(nameIdx);
943 if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
945 }
946 }
947
948 return total ;
949
950 } else {
951 // Simply return integral
952 double ret = integrate() / jacobianProduct() ;
953 return ret ;
954 }
955}
956
957////////////////////////////////////////////////////////////////////////////////
958/// Perform hybrid numerical/analytical integration over all real-valued dependents
959
961{
962 if (!_numIntEngine) {
963 // Trivial case, fully analytical integration
965 } else {
966 return _numIntEngine->calculate() ;
967 }
968}
969
970////////////////////////////////////////////////////////////////////////////////
971/// Intercept server redirects and reconfigure internal object accordingly
972
974 bool mustReplaceAll, bool nameChange, bool isRecursive)
975{
977
979
980 // Update contents value caches for _intList and _sumList
985
986 // Delete parameters cache if we have one
987 _params.reset();
988
990}
991
992////////////////////////////////////////////////////////////////////////////////
993
995{
996 if (!_params) {
997 _params = std::make_unique<RooArgSet>("params") ;
998
999 RooArgSet params ;
1000 for (const auto server : _serverList) {
1001 if (server->isValueServer(*this)) _params->add(*server) ;
1002 }
1003 }
1004
1005 return *_params ;
1006}
1007
1008////////////////////////////////////////////////////////////////////////////////
1009/// Check if current value is valid
1010
1011bool RooRealIntegral::isValidReal(double /*value*/, bool /*printError*/) const
1012{
1013 return true ;
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// Check if component selection is allowed
1018
1022
1023////////////////////////////////////////////////////////////////////////////////
1024/// Set component selection to be allowed/forbidden
1025
1029
1030////////////////////////////////////////////////////////////////////////////////
1031/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1032/// integration operation
1033
1034void RooRealIntegral::printMetaArgs(std::ostream& os) const
1035{
1036 if (!intVars().empty()) {
1037 os << "Int " ;
1038 }
1039 os << _function->GetName() ;
1040 if (_funcNormSet) {
1041 os << "_Norm" << *_funcNormSet << " " ;
1042 }
1043
1044 // List internally integrated observables and factorizing observables as analytically integrated
1046 tmp.add(_facList) ;
1047 if (!tmp.empty()) {
1048 os << "d[Ana]" << tmp << " ";
1049 }
1050
1051 // List numerically integrated and summed observables as numerically integrated
1053 tmp2.add(_sumList) ;
1054 if (!tmp2.empty()) {
1055 os << " d[Num]" << tmp2 << " ";
1056 }
1057}
1058
1059////////////////////////////////////////////////////////////////////////////////
1060/// Print the state of this object to the specified output stream.
1061
1062void RooRealIntegral::printMultiline(std::ostream& os, Int_t contents, bool verbose, TString indent) const
1063{
1064 RooAbsReal::printMultiline(os,contents,verbose,indent) ;
1065 os << indent << "--- RooRealIntegral ---" << std::endl;
1066 os << indent << " Integrates ";
1069 deeper.Append(" ");
1070 os << indent << " operating mode is "
1071 << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << std::endl ;
1072 os << indent << " Summed discrete args are " << _sumList << std::endl ;
1073 os << indent << " Numerically integrated args are " << _intList << std::endl;
1074 os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << std::endl ;
1075 os << indent << " Arguments included in Jacobian are " << _jacList << std::endl ;
1076 os << indent << " Factorized arguments are " << _facList << std::endl ;
1077 os << indent << " Function normalization set " ;
1078 if (_funcNormSet) {
1079 _funcNormSet->Print("1") ;
1080 } else {
1081 os << "<none>";
1082 }
1083
1084 os << std::endl ;
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1089
1091{
1092 _cacheAllNDim = ndim;
1093}
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Return minimum dimensions of numeric integration for which values are cached.
1097
1102
1103std::unique_ptr<RooAbsArg>
1108
1109/// Sort numeric integration variables in summation and integration lists.
1110/// To be used during construction.
1112{
1113 if (dynamic_cast<RooAbsRealLValue const *>(&arg)) {
1114 _intList.add(arg, true);
1115 } else if (dynamic_cast<RooAbsCategoryLValue const *>(&arg)) {
1116 _sumList.add(arg, true);
1117 }
1118}
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
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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:2489
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
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.
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition RooAbsArg.h:444
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:149
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.
virtual bool isDerived() const
Does value or shape of this arg depend on any other arg?
Definition RooAbsArg.h:97
bool isValueOrShapeDirtyAndClear() const
Definition RooAbsArg.h:396
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 _serverList
Definition RooAbsArg.h:573
Abstract base class for objects that represent a discrete value that can be set from the outside,...
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.
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
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for objects that are lvalues, i.e.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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:535
double traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
RooFit::UniqueId< RooArgSet >::Value_t _lastNormSetId
Component selection flag for RooAbsPdf::plotCompOn.
Definition RooAbsReal.h:542
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:348
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 ...
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:80
bool isValueServer() const
Returns true of contents is value server of owner.
Definition RooArgProxy.h:60
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
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...
Represents a constant real-valued object.
Definition RooConstVar.h:23
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,...
Registry for const char* names.
Definition RooNameReg.h:26
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:39
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
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,...
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 setAllowComponentSelection(bool allow)
Set component selection to be allowed/forbidden.
RooRealProxy _function
Function being integrated.
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
!
void addNumIntDep(RooAbsArg const &arg)
Sort numeric integration variables in summation and integration lists.
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.
double getValV(const RooArgSet *set=nullptr) const override
Return value of object.
RooSetProxy _anaList
Set of observables over which is integrated/summed analytically.
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
!
virtual double integrate() const
Perform hybrid numerical/analytical integration over all real-valued dependents.
RooListProxy _sumCat
!
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.
static Int_t _cacheAllNDim
! Cache all integrals with given numeric dimension
RooArgSet const * actualFuncNormSet() const
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...
const char * intRange() const
bool getAllowComponentSelection() const
Check if component selection is allowed.
Joins several RooAbsCategoryLValue objects into a single category.
bool setArg(T &newRef)
Change object held in proxy into newRef.
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:49
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
Basic string class.
Definition TString.h:139
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:35
constexpr Value_t value() const
Return numerical value of ID.
Definition UniqueId.h:59