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