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