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
51
52namespace {
53
54/// Utility function that returns true if 'object server' is a server
55/// to exactly one of the RooAbsArgs in 'exclLVBranches'
56bool servesExclusively(const RooAbsArg *server, const RooArgSet &exclLVBranches, const RooArgSet &allBranches)
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
72 Int_t numLVServ(0);
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) {
77 if (!servesExclusively(client, exclLVBranches, allBranches)) {
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
153 RooArgList allArgsList;
154 function.treeNodeServerList(&allArgsList, nullptr, true, true, /*valueOnly=*/false, false);
155 RooArgSet allArgs{allArgsList};
156 allArgs.sortTopologically();
157
158 // Figure out what are all the value servers only
159 RooArgList allValueArgsList;
160 function.treeNodeServerList(&allValueArgsList, nullptr, true, true, /*valueOnly=*/true, false);
161 RooArgSet allValueArgs{allValueArgsList};
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())) {
173 unmarkDepValueClients(*depInArgs, allArgs, marked);
174 addObservableToServers(function, *depInArgs, serversToAdd, rangeName);
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) {
185 addParameterToServers(function, *server, serversToAdd, !allValueArgs.find(*server));
186 marked[index] = MarkedState::AlreadyAdded;
187 }
188 }
189 }
190 }
191
192 return serversToAdd;
193}
194
195void fillAnIntOKDepList(RooAbsReal const &function, RooArgSet const &intDeps, RooArgSet &anIntOKDepList,
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.
204 RooArgSet filteredIntDeps;
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) {
215 filteredIntDeps.add(*intDep);
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
279{
281}
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// 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
360 RooArgSet intDepList(depList) ;
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())};
375 _facList.add(*argClone);
376 addOwnedComponents(std::move(argClone));
377 }
378 }
379
380 if (!_facList.empty()) {
381 oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables 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") ;
393 RooArgSet branchList;
394 RooArgSet branchListVD;
395 function.branchNodeServerList(&branchList) ;
396
397 for (auto branch: branchList) {
398 RooAbsRealLValue *realArgLV = dynamic_cast<RooAbsRealLValue*>(branch) ;
399 RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(branch) ;
400 if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
401 exclLVBranches.add(*branch) ;
402 }
403 if (branch != &function && function.dependsOnValue(*branch)) {
404 branchListVD.add(*branch) ;
405 }
406 }
407 exclLVBranches.remove(depList,true,true) ;
408
409 // Initial fill of list of LValue leaf servers (put in intDepList, but the
410 // instances that are in the actual computation graph of the function)
411 RooArgSet exclLVServers("exclLVServers") ;
412 function.getObservables(&intDepList, exclLVServers);
413
414 // Obtain mutual exclusive dependence by iterative reduction
415 bool converged(false) ;
416 while(!converged) {
417 converged=true ;
418
419 // Reduce exclLVServers to only those serving exclusively exclLVBranches
420 std::vector<RooAbsArg*> toBeRemoved;
421 for (auto server : exclLVServers) {
422 if (!servesExclusively(server,exclLVBranches,branchListVD)) {
423 toBeRemoved.push_back(server);
424 converged=false ;
425 }
426 }
427 exclLVServers.remove(toBeRemoved.begin(), toBeRemoved.end());
428
429 // Reduce exclLVBranches to only those depending exclusively on exclLVservers
430 // Attention: counting loop, since erasing from container
431 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
432 const RooAbsArg* branch = exclLVBranches[i];
433 RooArgSet brDepList;
434 branch->getObservables(&intDepList, brDepList);
435 RooArgSet bsList(brDepList,"bsList") ;
436 bsList.remove(exclLVServers,true,true) ;
437 if (!bsList.empty()) {
438 exclLVBranches.remove(*branch,true,true) ;
439 --i;
440 converged=false ;
441 }
442 }
443 }
444
445 // Eliminate exclLVBranches that do not depend on any LVServer
446 // Attention: Counting loop, since modifying container
447 for (std::size_t i=0; i < exclLVBranches.size(); ++i) {
448 const RooAbsArg* branch = exclLVBranches[i];
449 if (!branch->dependsOnValue(exclLVServers)) {
450 exclLVBranches.remove(*branch,true,true) ;
451 --i;
452 }
453 }
454
455 // Replace exclusive lvalue branch servers with lvalue branches
456 // WVE Don't do this for binned distributions - deal with this using numeric integration with transformed bin boundaroes
457 if (!exclLVServers.empty() && !function.isBinnedDistribution(exclLVBranches)) {
458 intDepList.remove(exclLVServers) ;
459 intDepList.add(exclLVBranches) ;
460 }
461
462
463 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
464 // * C) Check for dependents that the PDF insists on integrating *
465 // analytically itself *
466 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
467
468 RooArgSet anIntOKDepList ;
469 for (auto arg : intDepList) {
470 if (function.forceAnalyticalInt(*arg)) {
471 anIntOKDepList.add(*arg) ;
472 }
473 }
474
475 if (!anIntOKDepList.empty()) {
476 oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << std::endl ;
477 }
478
479
480 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
481 // * D) Make list of servers that can be integrated analytically *
482 // Add all parameters/dependents as value/shape servers *
483 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
484
485 auto serversToAdd = getValueAndShapeServers(function, depList, rangeName);
486 fillAnIntOKDepList(function, intDepList, anIntOKDepList, branchListVD);
487 // We will not add the servers just now, because it makes only sense to add
488 // them once we have made sure that this integral is not operating in
489 // pass-through mode. It will be done at the end of this constructor.
490
491 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
492 // * E) interact with function to make list of objects actually integrated analytically *
493 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
494
495 _mode = function.getAnalyticalIntegralWN(anIntOKDepList,_anaList,_funcNormSet.get(),RooNameReg::str(_rangeName)) ;
496
497 // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaList
498 if (_mode==0) {
500 }
501
502 if (_mode!=0) {
503 oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << std::endl ;
504 }
505
506 // WVE kludge: synchronize dset for use in analyticalIntegral
507 // LM : I think this is needed only if _funcNormSet is not an empty set
508 if (_funcNormSet && !_funcNormSet->empty()) {
509 function.getVal(_funcNormSet.get()) ;
510 }
511
512 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
513 // * F) Make list of numerical integration variables consisting of: *
514 // * - Category dependents of RealLValues in analytical integration *
515 // * - Expanded server lists of server that are not analytically integrated *
516 // * Make Jacobian list with analytically integrated RealLValues *
517 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
518
519 // Loop over actually analytically integrated dependents
520 for (const auto arg : _anaList) {
521
522 // Process only derived RealLValues
523 if (dynamic_cast<RooAbsRealLValue const *>(arg) && arg->isDerived() && !arg->isFundamental()) {
524
525 // Add to list of Jacobians to calculate
526 _jacList.add(*arg) ;
527
528 // Add category dependent of LValueReal used in integration
529 std::unique_ptr<RooArgSet> argDepList{arg->getObservables(&intDepList)};
530 for (const auto argDep : *argDepList) {
531 if (dynamic_cast<RooAbsCategoryLValue const *>(argDep) && intDepList.contains(*argDep)) {
532 addNumIntDep(*argDep);
533 }
534 }
535 }
536 }
537
538
539 // If nothing was integrated analytically, swap back LVbranches for LVservers for subsequent numeric integration
540 if (_anaList.empty()) {
541 if (!exclLVServers.empty()) {
542 //cout << "NUMINT phase analList is empty. exclLVServers = " << exclLVServers << std::endl ;
543 intDepList.remove(exclLVBranches) ;
544 intDepList.add(exclLVServers) ;
545 }
546 }
547 //cout << "NUMINT intDepList = " << intDepList << std::endl ;
548
549 // Loop again over function servers to add remaining numeric integrations
550 for (const auto arg : function.servers()) {
551
552 // Process only servers that are not treated analytically
553 if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
554
555 // Process only derived RealLValues
556 if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
557 addNumIntDep(*arg) ;
558 } else {
559
560 // WVE this will only get the observables, but not l-value transformations
561 // Expand server in final dependents
562 auto argDeps = std::unique_ptr<RooArgSet>(arg->getObservables(&intDepList));
563
564 // Add final dependents, that are not forcibly integrated analytically,
565 // to numerical integration list
566 for (const auto dep : *argDeps) {
567 if (!_anaList.find(dep->GetName())) {
568 addNumIntDep(*dep);
569 }
570 }
571 }
572 }
573 }
574
575 if (!_anaList.empty()) {
576 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << std::endl ;
577 }
578 if (!_intList.empty()) {
579 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << std::endl ;
580 }
581 if (!_sumList.empty()) {
582 oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << std::endl ;
583 }
584
585
586 // Determine operating mode
587 if (!_intList.empty() || !_sumList.empty()) {
588 // Numerical and optional Analytical integration
590 } else if (!_anaList.empty()) {
591 // Purely analytical integration
593 } else {
594 // No integration performed, where the function is a direct value server
596 _function._valueServer = true;
597 }
598 // We are only setting the function proxy now that it's clear if it's a value
599 // server or not.
600 _function.setArg(const_cast<RooAbsReal&>(function));
601
602 // Determine auto-dirty status
604
605 // Create value caches for _intList and _sumList
608
609
610 if (!_sumList.empty()) {
611 _sumCat.addOwned(std::make_unique<RooSuperCategory>(Form("%s_sumCat",GetName()),"sumCat",_sumList));
612 }
613
614 // Only if we are not in pass-through mode we need to add the shape and value
615 // servers separately.
617 for(auto const& toAdd : serversToAdd) {
618 addServer(*toAdd.arg, !toAdd.isShapeServer, toAdd.isShapeServer);
619 }
620 }
621
623}
624
625////////////////////////////////////////////////////////////////////////////////
626/// Set appropriate cache operation mode for integral depending on cache operation
627/// mode of server objects
628
630{
631 // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
632 for (const auto server : _serverList) {
633 if (server->isValueServer(*this)) {
634 RooArgSet leafSet ;
635 server->leafNodeServerList(&leafSet) ;
636 for (const auto leaf : leafSet) {
637 if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {
639 break ;
640 }
641 if (leaf->getAttribute("projectedDependent")) {
643 break ;
644 }
645 }
646 }
647 }
648}
649
650////////////////////////////////////////////////////////////////////////////////
651/// (Re)Initialize numerical integration engine if necessary. Return true if
652/// successful, or otherwise false.
653
655{
656 // if we already have an engine, check if it still works for the present limits.
657 if(_numIntEngine) {
658 if(_numIntEngine->isValid() && _numIntEngine->checkLimits() && !_restartNumIntEngine ) return true;
659 // otherwise, cleanup the old engine
660 _numIntEngine.reset();
661 _numIntegrand.reset();
662 }
663
664 // All done if there are no arguments to integrate numerically
665 if(_intList.empty()) return true;
666
667 // Bind the appropriate analytic integral of our RooRealVar object to
668 // those of its arguments that will be integrated out numerically.
669 if(_mode != 0) {
670 std::unique_ptr<RooAbsReal> analyticalPart{_function->createIntegral(_anaList, *actualFuncNormSet(), RooNameReg::str(_rangeName))};
671 _numIntegrand = std::make_unique<RooRealBinding>(*analyticalPart,_intList,nullptr,false,_rangeName);
672 const_cast<RooRealIntegral*>(this)->addOwnedComponents(std::move(analyticalPart));
673 }
674 else {
675 _numIntegrand = std::make_unique<RooRealBinding>(*_function,_intList,actualFuncNormSet(),false,_rangeName);
676 }
677 if(nullptr == _numIntegrand || !_numIntegrand->isValid()) {
678 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << std::endl;
679 return false;
680 }
681
682 // Create appropriate numeric integrator using factory
683 bool isBinned = _function->isBinnedDistribution(_intList) ;
684 std::string integratorName = RooNumIntFactory::instance().getIntegratorName(*_numIntegrand,*_iconfig,0,isBinned);
686
687 if(_numIntEngine == nullptr || !_numIntEngine->isValid()) {
688 coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << std::endl;
689 return false;
690 }
691
692 cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator "
693 << integratorName << " to calculate Int" << _intList << std::endl ;
694
695 if (_intList.size()>3) {
696 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 ;
697 }
698
699 _restartNumIntEngine = false ;
700 return true;
701}
702
703////////////////////////////////////////////////////////////////////////////////
704/// Copy constructor
705
707 : RooAbsReal(other, name),
708 _valid(other._valid),
709 _respectCompSelect(other._respectCompSelect),
710 _sumList("!sumList", this, other._sumList),
711 _intList("!intList", this, other._intList),
712 _anaList("!anaList", this, other._anaList),
713 _jacList("!jacList", this, other._jacList),
714 _facList("!facList", this, other._facList),
715 _function("!func", this, other._function),
716 _iconfig(other._iconfig),
717 _sumCat("!sumCat", this, other._sumCat),
718 _mode(other._mode),
719 _intOperMode(other._intOperMode),
720 _rangeName(other._rangeName)
721{
722 if(other._funcNormSet) {
723 _funcNormSet = std::make_unique<RooArgSet>();
724 other._funcNormSet->snapshot(*_funcNormSet, false);
725 }
726
727 other._intList.snapshot(_saveInt) ;
728 other._sumList.snapshot(_saveSum) ;
729
731}
732
733////////////////////////////////////////////////////////////////////////////////
734
736{
738}
739
740////////////////////////////////////////////////////////////////////////////////
741
742RooFit::OwningPtr<RooAbsReal> RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const
743{
744 // Handle special case of no integration with default algorithm
745 if (iset.empty()) {
746 return RooAbsReal::createIntegral(iset,nset,cfg,rangeName) ;
747 }
748
749 // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
750 RooArgSet isetAll(iset) ;
751 isetAll.add(_sumList) ;
752 isetAll.add(_intList) ;
753 isetAll.add(_anaList) ;
754 isetAll.add(_facList) ;
755
756 const RooArgSet* newNormSet(nullptr) ;
757 std::unique_ptr<RooArgSet> tmp;
758 if (nset && !_funcNormSet) {
759 newNormSet = nset ;
760 } else if (!nset && _funcNormSet) {
761 newNormSet = _funcNormSet.get();
762 } else if (nset && _funcNormSet) {
763 tmp = std::make_unique<RooArgSet>();
764 tmp->add(*nset) ;
765 tmp->add(*_funcNormSet,true) ;
766 newNormSet = tmp.get();
767 }
768 return _function->createIntegral(isetAll,newNormSet,cfg,rangeName);
769}
770
771////////////////////////////////////////////////////////////////////////////////
772/// Return value of object. If the cache is clean, return the
773/// cached value, otherwise recalculate on the fly and refill
774/// the cache
775
776double RooRealIntegral::getValV(const RooArgSet* nset) const
777{
778// // fast-track clean-cache processing
779// if (_operMode==AClean) {
780// return _value ;
781// }
782
783 if (nset && nset->uniqueId().value() != _lastNormSetId) {
784 const_cast<RooRealIntegral*>(this)->setProxyNormSet(nset);
785 _lastNormSetId = nset->uniqueId().value();
786 }
787
789 _value = traceEval(nset) ;
790 }
791
792 return _value ;
793}
794
795////////////////////////////////////////////////////////////////////////////////
796/// Perform the integration and return the result
797
799{
801
802 double retVal(0) ;
803 switch (_intOperMode) {
804
805 case Hybrid:
806 {
807 // Cache numeric integrals in >1d expensive object cache
808 RooDouble const* cacheVal(nullptr) ;
809 if ((_cacheNum && !_intList.empty()) || int(_intList.size())>=_cacheAllNDim) {
810 cacheVal = static_cast<RooDouble const*>(expensiveObjectCache().retrieveObject(GetName(),RooDouble::Class(),parameters())) ;
811 }
812
813 if (cacheVal) {
814 retVal = *cacheVal ;
815 // cout << "using cached value of integral" << GetName() << std::endl ;
816 } else {
817
818
819 // Find any function dependents that are AClean
820 // and switch them temporarily to ADirty
821 bool origState = inhibitDirty() ;
822 setDirtyInhibit(true) ;
823
824 // try to initialize our numerical integration engine
825 if(!(_valid= initNumIntegrator())) {
826 coutE(Integration) << ClassName() << "::" << GetName()
827 << ":evaluate: cannot initialize numerical integrator" << std::endl;
828 return 0;
829 }
830
831 // Save current integral dependent values
834
835 // Evaluate sum/integral
836 retVal = sum() ;
837
838
839 // This must happen BEFORE restoring dependents, otherwise no dirty state propagation in restore step
840 setDirtyInhibit(origState) ;
841
842 // Restore integral dependent values
845
846 // Cache numeric integrals in >1d expensive object cache
847 if ((_cacheNum && !_intList.empty()) || int(_intList.size())>=_cacheAllNDim) {
848 RooDouble* val = new RooDouble(retVal) ;
850 // cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << std::endl ;
851 }
852
853 }
854 break ;
855 }
856 case Analytic:
857 {
859 cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName()
860 << ")func = " << _function->ClassName() << "::" << _function->GetName()
861 << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << std::endl ;
862
863
864 break ;
865 }
866
867 case PassThrough:
868 {
869 // In pass through mode, the RooRealIntegral should have registered the
870 // function as a value server, because we directly depend on its value.
871 assert(_function.isValueServer());
872 // There should be no other servers besides the actual function and the
873 // factorized observables that the function doesn't depend on but are
874 // integrated over later.
875 assert(servers().size() == _facList.size() + 1);
876
877 //setDirtyInhibit(true) ;
879 //setDirtyInhibit(false) ;
880 break ;
881 }
882 }
883
884
885 // Multiply answer with integration ranges of factorized variables
886 for (const auto arg : _facList) {
887 // Multiply by fit range for 'real' dependents
888 if (auto argLV = dynamic_cast<RooAbsRealLValue *>(arg)) {
889 retVal *= (argLV->getMax(intRange()) - argLV->getMin(intRange())) ;
890 }
891 // Multiply by number of states for category dependents
892 if (auto argLV = dynamic_cast<RooAbsCategoryLValue *>(arg)) {
893 retVal *= argLV->numTypes() ;
894 }
895 }
896
897
898 if (dologD(Tracing)) {
899 cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
900 switch(_intOperMode) {
901 case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
902 case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
903 case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
904 }
905
906 ccxcoutD(Tracing) << "raw*fact = " << retVal << std::endl ;
907 }
908
909 return retVal ;
910}
911
912////////////////////////////////////////////////////////////////////////////////
913/// Return product of jacobian terms originating from analytical integration
914
916{
917 if (_jacList.empty()) {
918 return 1 ;
919 }
920
921 double jacProd(1) ;
922 for (const auto elm : _jacList) {
923 auto arg = static_cast<const RooAbsRealLValue*>(elm);
924 jacProd *= arg->jacobian() ;
925 }
926
927 // Take std::abs() here: if jacobian is negative, min and max are swapped and analytical integral
928 // will be positive, so must multiply with positive jacobian.
929 return std::abs(jacProd) ;
930}
931
932////////////////////////////////////////////////////////////////////////////////
933/// Perform summation of list of category dependents to be integrated
934
936{
937 if (!_sumList.empty()) {
938 // Add integrals for all permutations of categories summed over
939 double total(0) ;
940
941 RooSuperCategory* sumCat = static_cast<RooSuperCategory*>(_sumCat.first()) ;
942 for (const auto& nameIdx : *sumCat) {
943 sumCat->setIndex(nameIdx);
944 if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
946 }
947 }
948
949 return total ;
950
951 } else {
952 // Simply return integral
953 double ret = integrate() / jacobianProduct() ;
954 return ret ;
955 }
956}
957
958////////////////////////////////////////////////////////////////////////////////
959/// Perform hybrid numerical/analytical integration over all real-valued dependents
960
962{
963 if (!_numIntEngine) {
964 // Trivial case, fully analytical integration
966 } else {
967 return _numIntEngine->calculate() ;
968 }
969}
970
971////////////////////////////////////////////////////////////////////////////////
972/// Intercept server redirects and reconfigure internal object accordingly
973
975 bool mustReplaceAll, bool nameChange, bool isRecursive)
976{
977 _restartNumIntEngine = true ;
978
980
981 // Update contents value caches for _intList and _sumList
986
987 // Delete parameters cache if we have one
988 _params.reset();
989
990 return RooAbsReal::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
991}
992
993////////////////////////////////////////////////////////////////////////////////
994
996{
997 if (!_params) {
998 _params = std::make_unique<RooArgSet>("params") ;
999
1000 RooArgSet params ;
1001 for (const auto server : _serverList) {
1002 if (server->isValueServer(*this)) _params->add(*server) ;
1003 }
1004 }
1005
1006 return *_params ;
1007}
1008
1009////////////////////////////////////////////////////////////////////////////////
1010/// Check if current value is valid
1011
1012bool RooRealIntegral::isValidReal(double /*value*/, bool /*printError*/) const
1013{
1014 return true ;
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Check if component selection is allowed
1019
1021 return _respectCompSelect;
1022}
1023
1024////////////////////////////////////////////////////////////////////////////////
1025/// Set component selection to be allowed/forbidden
1026
1028 _respectCompSelect = allow;
1029}
1030
1031////////////////////////////////////////////////////////////////////////////////
1032/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
1033/// integration operation
1034
1035void RooRealIntegral::printMetaArgs(std::ostream& os) const
1036{
1037 if (!intVars().empty()) {
1038 os << "Int " ;
1039 }
1040 os << _function->GetName() ;
1041 if (_funcNormSet) {
1042 os << "_Norm" << *_funcNormSet << " " ;
1043 }
1044
1045 // List internally integrated observables and factorizing observables as analytically integrated
1046 RooArgSet tmp(_anaList) ;
1047 tmp.add(_facList) ;
1048 if (!tmp.empty()) {
1049 os << "d[Ana]" << tmp << " ";
1050 }
1051
1052 // List numerically integrated and summed observables as numerically integrated
1053 RooArgSet tmp2(_intList) ;
1054 tmp2.add(_sumList) ;
1055 if (!tmp2.empty()) {
1056 os << " d[Num]" << tmp2 << " ";
1057 }
1058}
1059
1060////////////////////////////////////////////////////////////////////////////////
1061/// Print the state of this object to the specified output stream.
1062
1063void RooRealIntegral::printMultiline(std::ostream& os, Int_t contents, bool verbose, TString indent) const
1064{
1065 RooAbsReal::printMultiline(os,contents,verbose,indent) ;
1066 os << indent << "--- RooRealIntegral ---" << std::endl;
1067 os << indent << " Integrates ";
1069 TString deeper(indent);
1070 deeper.Append(" ");
1071 os << indent << " operating mode is "
1072 << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << std::endl ;
1073 os << indent << " Summed discrete args are " << _sumList << std::endl ;
1074 os << indent << " Numerically integrated args are " << _intList << std::endl;
1075 os << indent << " Analytically integrated args using mode " << _mode << " are " << _anaList << std::endl ;
1076 os << indent << " Arguments included in Jacobian are " << _jacList << std::endl ;
1077 os << indent << " Factorized arguments are " << _facList << std::endl ;
1078 os << indent << " Function normalization set " ;
1079 if (_funcNormSet) {
1080 _funcNormSet->Print("1") ;
1081 } else {
1082 os << "<none>";
1083 }
1084
1085 os << std::endl ;
1086}
1087
1088////////////////////////////////////////////////////////////////////////////////
1089/// Global switch to cache all integral values that integrate at least ndim dimensions numerically
1090
1092{
1093 _cacheAllNDim = ndim;
1094}
1095
1096////////////////////////////////////////////////////////////////////////////////
1097/// Return minimum dimensions of numeric integration for which values are cached.
1098
1100{
1101 return _cacheAllNDim;
1102}
1103
1104std::unique_ptr<RooAbsArg>
1106{
1108}
1109
1110/// Sort numeric integration variables in summation and integration lists.
1111/// To be used during construction.
1113{
1114 if (dynamic_cast<RooAbsRealLValue const *>(&arg)) {
1115 _intList.add(arg, true);
1116 } else if (dynamic_cast<RooAbsCategoryLValue const *>(&arg)) {
1117 _sumList.add(arg, true);
1118 }
1119}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define cxcoutI(a)
#define cxcoutD(a)
#define oocxcoutD(o, a)
#define dologD(a)
#define coutE(a)
#define ccxcoutD(a)
#define ccoutD(a)
#define oocxcoutI(o, a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
int Int_t
Definition RtypesCore.h:45
#define ClassImp(name)
Definition Rtypes.h:382
static void indent(ostringstream &buf, int indent_level)
static unsigned int total
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
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.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
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
bool dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr) const
Check whether this object depends on values from an element in the serverList.
Definition RooAbsArg.h:106
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 getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
bool isValueOrShapeDirtyAndClear() const
Definition RooAbsArg.h:396
bool inhibitDirty() const
Delete watch flag.
bool hasClients() const
Definition RooAbsArg.h:126
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
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=nullptr, bool recurseNonDerived=false) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
bool isValueServer(const RooAbsArg &arg) const
Check if this is serving values to arg.
Definition RooAbsArg.h:166
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:425
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.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
Storage_t const & get() const
Const access to the underlying stl container.
void sortTopologically()
Sort collection topologically: the servers of any RooAbsArg will be before that RooAbsArg in the coll...
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
RooAbsArg * first() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract base class for objects that are lvalues, i.e.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
virtual bool isJacobianOK(const RooArgSet &depList) const
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:533
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:540
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:346
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:154
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.
std::unique_ptr< RooAbsIntegrator > createIntegrator(RooAbsFunc &func, const RooNumIntConfig &config, int ndim=0, bool isBinned=false) const
Construct a numeric integrator instance that operates on function 'func' and is configured with 'conf...
std::string getIntegratorName(RooAbsFunc &func, const RooNumIntConfig &config, int ndim=0, bool isBinned=false) const
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 setIndex(value_type index, bool printError=true) override
Set the value of the super category to the specified index.
bool inRange(const char *rangeName) const override
Check that all input category states are in the given range.
bool setArg(T &newRef)
Change object held in proxy into newRef.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:572
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