Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsArg.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/** \class RooAbsArg
19 \ingroup Roofitcore
20
21RooAbsArg is the common abstract base class for objects that
22represent a value and a "shape" in RooFit. Values or shapes usually depend on values
23or shapes of other RooAbsArg instances. Connecting several RooAbsArg in
24a computation graph models an expression tree that can be evaluated.
25
26### Building a computation graph of RooFit objects
27Therefore, RooAbsArg provides functionality to connect objects of type RooAbsArg into
28a computation graph to pass values between those objects.
29A value can e.g. be a real-valued number, (instances of RooAbsReal), or an integer, that is,
30catgory index (instances of RooAbsCategory). The third subclass of RooAbsArg is RooStringVar,
31but it is rarely used.
32
33The "shapes" that a RooAbsArg can possess can e.g. be the definition
34range of an observable, or how many states a category object has. In computations,
35values are expected to change often, while shapes remain mostly constant
36(unless e.g. a new range is set for an observable).
37
38Nodes of a computation graph are connected using instances of RooAbsProxy.
39If Node B declares a member `RooTemplateProxy<TypeOfNodeA>`, Node A will be
40registered as a server of values to Node B, and Node B will know that it is
41a client of node A. Using functions like dependsOn(), or getObservables()
42/ getParameters(), the relation of `A --> B` can be queried. Using graphVizTree(),
43one can create a visualisation of the expression tree.
44
45
46An instance of RooAbsArg can have named attributes. It also has flags
47to indicate that either its value or its shape were changed (= it is dirty).
48RooAbsArg provides functionality to manage client/server relations in
49a computation graph (\ref clientServerInterface), and helps propagating
50value/shape changes through the graph. RooAbsArg implements interfaces
51for inspecting client/server relationships (\ref clientServerInterface) and
52setting/clearing/querying named attributes.
53
54### Caching of values
55The values of nodes in the computation graph are cached in RooFit. If
56a value is used in two nodes of a graph, it doesn't need to be recomputed. If
57a node acquires a new value, it notifies its consumers ("clients") that
58their cached values are dirty. See the functions in \ref optimisationInterface
59for details.
60A node uses its isValueDirty() and isShapeDirty() functions to decide if a
61computation is necessary. Caching can be vetoed globally by setting a
62bit using setDirtyInhibit(). This will make computations slower, but all the
63nodes of the computation graph will be evaluated irrespective of whether their
64state is clean or dirty. Using setOperMode(), caching can also be enabled/disabled
65for single nodes.
66
67*/
68
69#include "RooFit.h"
70
71#include "TBuffer.h"
72#include "TClass.h"
74#include "strlcpy.h"
75
76#include "RooSecondMoment.h"
77#include "RooWorkspace.h"
78
79#include "RooMsgService.h"
80#include "RooAbsArg.h"
81#include "RooArgSet.h"
82#include "RooArgProxy.h"
83#include "RooSetProxy.h"
84#include "RooListProxy.h"
85#include "RooAbsData.h"
87#include "RooAbsRealLValue.h"
88#include "RooTrace.h"
89#include "RooRealIntegral.h"
90#include "RooConstVar.h"
92#include "RooAbsDataStore.h"
93#include "RooResolutionModel.h"
94#include "RooVectorDataStore.h"
95#include "RooTreeDataStore.h"
96#include "RooHelpers.h"
97
98#include <algorithm>
99#include <cstring>
100#include <fstream>
101#include <iostream>
102#include <memory>
103#include <sstream>
104
105using namespace std;
106
108;
109
113
114std::map<RooAbsArg*,std::unique_ptr<TRefArray>> RooAbsArg::_ioEvoList;
115std::stack<RooAbsArg*> RooAbsArg::_ioReadStack ;
116
117
118////////////////////////////////////////////////////////////////////////////////
119/// Default constructor
120
122 : TNamed(), _deleteWatch(kFALSE), _valueDirty(kTRUE), _shapeDirty(kTRUE), _operMode(Auto), _fast(kFALSE), _ownedComponents(nullptr),
123 _prohibitServerRedirect(kFALSE), _namePtr(0), _isConstant(kFALSE), _localNoInhibitDirty(kFALSE),
124 _myws(0)
125{
127
128}
129
130////////////////////////////////////////////////////////////////////////////////
131/// Create an object with the specified name and descriptive title.
132/// The newly created object has no clients or servers and has its
133/// dirty flags set.
134
135RooAbsArg::RooAbsArg(const char *name, const char *title)
136 : TNamed(name, title), _deleteWatch(kFALSE), _valueDirty(kTRUE), _shapeDirty(kTRUE), _operMode(Auto), _fast(kFALSE),
137 _ownedComponents(0), _prohibitServerRedirect(kFALSE), _namePtr(0), _isConstant(kFALSE),
138 _localNoInhibitDirty(kFALSE), _myws(0)
139{
140 if (name == nullptr || strlen(name) == 0) {
141 throw std::logic_error("Each RooFit object needs a name. "
142 "Objects representing the same entity (e.g. an observable 'x') are identified using their name.");
143 }
145}
146
147////////////////////////////////////////////////////////////////////////////////
148/// Copy constructor transfers all boolean and string properties of the original
149/// object. Transient properties and client-server links are not copied
150
151RooAbsArg::RooAbsArg(const RooAbsArg &other, const char *name)
152 : TNamed(name ? name : other.GetName(), other.GetTitle()), RooPrintable(other),
153 _boolAttrib(other._boolAttrib),
154 _stringAttrib(other._stringAttrib), _deleteWatch(other._deleteWatch), _operMode(Auto), _fast(kFALSE),
155 _ownedComponents(0), _prohibitServerRedirect(kFALSE),
156 _namePtr(name ? RooNameReg::instance().constPtr(name) : other._namePtr),
157 _isConstant(other._isConstant), _localNoInhibitDirty(other._localNoInhibitDirty), _myws(0)
158{
159
160 // Copy server list by hand
161 Bool_t valueProp, shapeProp ;
162 for (const auto server : other._serverList) {
163 valueProp = server->_clientListValue.containsByNamePtr(&other);
164 shapeProp = server->_clientListShape.containsByNamePtr(&other);
165 addServer(*server,valueProp,shapeProp) ;
166 }
167
168 setValueDirty() ;
169 setShapeDirty() ;
170
171 //setAttribute(Form("CloneOf(%08x)",&other)) ;
172 //cout << "RooAbsArg::cctor(" << this << ") #bools = " << _boolAttrib.size() << " #strings = " << _stringAttrib.size() << endl ;
173
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// Assign all boolean and string properties of the original
178/// object. Transient properties and client-server links are not assigned.
180 TNamed::operator=(other);
181 RooPrintable::operator=(other);
182 _boolAttrib = other._boolAttrib;
185 _operMode = other._operMode;
186 _fast = other._fast;
187 _ownedComponents = nullptr;
189 _namePtr = other._namePtr;
190 _isConstant = other._isConstant;
192 _myws = nullptr;
193
194 bool valueProp, shapeProp;
195 for (const auto server : other._serverList) {
196 valueProp = server->_clientListValue.containsByNamePtr(&other);
197 shapeProp = server->_clientListShape.containsByNamePtr(&other);
198 addServer(*server,valueProp,shapeProp) ;
199 }
200
203
204 return *this;
205}
206
207
208////////////////////////////////////////////////////////////////////////////////
209/// Destructor.
210
212{
213 // Notify all servers that they no longer need to serve us
214 while (!_serverList.empty()) {
216 }
217
218 // Notify all clients that they are in limbo
219 std::vector<RooAbsArg*> clientListTmp(_clientList.begin(), _clientList.end()); // have to copy, as we invalidate iterators
221 for (auto client : clientListTmp) {
222 client->setAttribute("ServerDied") ;
223 TString attr("ServerDied:");
224 attr.Append(GetName());
225 attr.Append(Form("(%zx)",(size_t)this)) ;
226 client->setAttribute(attr.Data());
227 client->removeServer(*this,kTRUE);
228
229 if (_verboseDirty) {
230
231 if (first) {
232 cxcoutD(Tracing) << "RooAbsArg::dtor(" << GetName() << "," << this << ") DeleteWatch: object is being destroyed" << endl ;
233 first = kFALSE ;
234 }
235
236 cxcoutD(Tracing) << fName << "::" << ClassName() << ":~RooAbsArg: dependent \""
237 << client->GetName() << "\" should have been deleted first" << endl ;
238 }
239 }
240
241 if (_ownedComponents) {
242 delete _ownedComponents ;
243 _ownedComponents = 0 ;
244 }
245
246}
247
248
249////////////////////////////////////////////////////////////////////////////////
250/// Control global dirty inhibit mode. When set to true no value or shape dirty
251/// flags are propagated and cache is always considered to be dirty.
252
254{
255 _inhibitDirty = flag ;
256}
257
258
259////////////////////////////////////////////////////////////////////////////////
260/// Activate verbose messaging related to dirty flag propagation
261
263{
264 _verboseDirty = flag ;
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// Check if this object was created as a clone of 'other'
269
271{
272 return (getAttribute(Form("CloneOf(%zx)",(size_t)&other)) ||
273 other.getAttribute(Form("CloneOf(%zx)",(size_t)this))) ;
274}
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// Set (default) or clear a named boolean attribute of this object.
279
281{
282 // Preserve backward compatibility - any strong
283 if(string("Constant")==name) {
284 _isConstant = value ;
285 }
286
287 if (value) {
288 _boolAttrib.insert(name) ;
289 } else {
290 set<string>::iterator iter = _boolAttrib.find(name) ;
291 if (iter != _boolAttrib.end()) {
292 _boolAttrib.erase(iter) ;
293 }
294
295 }
296
297}
298
299
300////////////////////////////////////////////////////////////////////////////////
301/// Check if a named attribute is set. By default, all attributes are unset.
302
304{
305 return (_boolAttrib.find(name) != _boolAttrib.end()) ;
306}
307
308
309////////////////////////////////////////////////////////////////////////////////
310/// Associate string 'value' to this object under key 'key'
311
312void RooAbsArg::setStringAttribute(const Text_t* key, const Text_t* value)
313{
314 if (value) {
315 _stringAttrib[key] = value ;
316 } else {
317 _stringAttrib.erase(key) ;
318 }
319}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Get string attribute mapped under key 'key'. Returns null pointer
323/// if no attribute exists under that key
324
326{
327 map<string,string>::const_iterator iter = _stringAttrib.find(key) ;
328 if (iter!=_stringAttrib.end()) {
329 return iter->second.c_str() ;
330 } else {
331 return 0 ;
332 }
333}
334
335
336////////////////////////////////////////////////////////////////////////////////
337/// Set (default) or clear a named boolean attribute of this object.
338
340{
341 if (value) {
342
343 _boolAttribTransient.insert(name) ;
344
345 } else {
346
347 set<string>::iterator iter = _boolAttribTransient.find(name) ;
348 if (iter != _boolAttribTransient.end()) {
349 _boolAttribTransient.erase(iter) ;
350 }
351
352 }
353
354}
355
356
357////////////////////////////////////////////////////////////////////////////////
358/// Check if a named attribute is set. By default, all attributes
359/// are unset.
360
362{
363 return (_boolAttribTransient.find(name) != _boolAttribTransient.end()) ;
364}
365
366
367
368
369////////////////////////////////////////////////////////////////////////////////
370/// Register another RooAbsArg as a server to us, ie, declare that
371/// we depend on it.
372/// \param server The server to be registered.
373/// \param valueProp In addition to the basic client-server relationship, declare dependence on the server's value.
374/// \param valueProp In addition to the basic client-server relationship, declare dependence on the server's shape.
375/// \param refCount Optionally add with higher reference count (if multiple components depend on it)
376
377void RooAbsArg::addServer(RooAbsArg& server, Bool_t valueProp, Bool_t shapeProp, std::size_t refCount)
378{
380 cxcoutF(LinkStateMgmt) << "RooAbsArg::addServer(" << this << "," << GetName()
381 << "): PROHIBITED SERVER ADDITION REQUESTED: adding server " << server.GetName()
382 << "(" << &server << ") for " << (valueProp?"value ":"") << (shapeProp?"shape":"") << endl ;
383 throw std::logic_error("PROHIBITED SERVER ADDITION REQUESTED in RooAbsArg::addServer");
384 }
385
386 cxcoutD(LinkStateMgmt) << "RooAbsArg::addServer(" << this << "," << GetName() << "): adding server " << server.GetName()
387 << "(" << &server << ") for " << (valueProp?"value ":"") << (shapeProp?"shape":"") << endl ;
388
389 if (server.operMode()==ADirty && operMode()!=ADirty && valueProp) {
391 }
392
393
394 // LM: use hash tables for larger lists
395// if (_serverList.GetSize() > 999 && _serverList.getHashTableSize() == 0) _serverList.setHashTableSize(1000);
396// if (server._clientList.GetSize() > 999 && server._clientList.getHashTableSize() == 0) server._clientList.setHashTableSize(1000);
397// if (server._clientListValue.GetSize() > 999 && server._clientListValue.getHashTableSize() == 0) server._clientListValue.setHashTableSize(1000);
398
399 // Add server link to given server
400 _serverList.Add(&server, refCount) ;
401
402 server._clientList.Add(this, refCount);
403 if (valueProp) server._clientListValue.Add(this, refCount);
404 if (shapeProp) server._clientListShape.Add(this, refCount);
405}
406
407
408
409////////////////////////////////////////////////////////////////////////////////
410/// Register a list of RooAbsArg as servers to us by calling
411/// addServer() for each arg in the list
412
413void RooAbsArg::addServerList(RooAbsCollection& serverList, Bool_t valueProp, Bool_t shapeProp)
414{
415 _serverList.reserve(_serverList.size() + serverList.size());
416
417 for (const auto arg : serverList) {
418 addServer(*arg,valueProp,shapeProp) ;
419 }
420}
421
422
423
424////////////////////////////////////////////////////////////////////////////////
425/// Unregister another RooAbsArg as a server to us, ie, declare that
426/// we no longer depend on its value and shape.
427
429{
431 cxcoutF(LinkStateMgmt) << "RooAbsArg::addServer(" << this << "," << GetName() << "): PROHIBITED SERVER REMOVAL REQUESTED: removing server "
432 << server.GetName() << "(" << &server << ")" << endl ;
433 assert(0) ;
434 }
435
436 if (_verboseDirty) {
437 cxcoutD(LinkStateMgmt) << "RooAbsArg::removeServer(" << GetName() << "): removing server "
438 << server.GetName() << "(" << &server << ")" << endl ;
439 }
440
441 // Remove server link to given server
442 _serverList.Remove(&server, force) ;
443
444 server._clientList.Remove(this, force) ;
445 server._clientListValue.Remove(this, force) ;
446 server._clientListShape.Remove(this, force) ;
447}
448
449
450////////////////////////////////////////////////////////////////////////////////
451/// Replace 'oldServer' with 'newServer'
452
453void RooAbsArg::replaceServer(RooAbsArg& oldServer, RooAbsArg& newServer, Bool_t propValue, Bool_t propShape)
454{
455 Int_t count = _serverList.refCount(&oldServer);
456 removeServer(oldServer, kTRUE);
457 addServer(newServer, propValue, propShape, count);
458}
459
460
461////////////////////////////////////////////////////////////////////////////////
462/// Change dirty flag propagation mask for specified server
463
464void RooAbsArg::changeServer(RooAbsArg& server, Bool_t valueProp, Bool_t shapeProp)
465{
466 if (!_serverList.containsByNamePtr(&server)) {
467 coutE(LinkStateMgmt) << "RooAbsArg::changeServer(" << GetName() << "): Server "
468 << server.GetName() << " not registered" << endl ;
469 return ;
470 }
471
472 // This condition should not happen, but check anyway
473 if (!server._clientList.containsByNamePtr(this)) {
474 coutE(LinkStateMgmt) << "RooAbsArg::changeServer(" << GetName() << "): Server "
475 << server.GetName() << " doesn't have us registered as client" << endl ;
476 return ;
477 }
478
479 // Remove all propagation links, then reinstall requested ones ;
480 Int_t vcount = server._clientListValue.refCount(this) ;
481 Int_t scount = server._clientListShape.refCount(this) ;
482 server._clientListValue.RemoveAll(this) ;
483 server._clientListShape.RemoveAll(this) ;
484 if (valueProp) {
485 server._clientListValue.Add(this, vcount) ;
486 }
487 if (shapeProp) {
488 server._clientListShape.Add(this, scount) ;
489 }
490}
491
492
493
494////////////////////////////////////////////////////////////////////////////////
495/// Fill supplied list with all leaf nodes of the arg tree, starting with
496/// ourself as top node. A leaf node is node that has no servers declared.
497
498void RooAbsArg::leafNodeServerList(RooAbsCollection* list, const RooAbsArg* arg, Bool_t recurseNonDerived) const
499{
500 treeNodeServerList(list,arg,kFALSE,kTRUE,kFALSE,recurseNonDerived) ;
501}
502
503
504
505////////////////////////////////////////////////////////////////////////////////
506/// Fill supplied list with all branch nodes of the arg tree starting with
507/// ourself as top node. A branch node is node that has one or more servers declared.
508
509void RooAbsArg::branchNodeServerList(RooAbsCollection* list, const RooAbsArg* arg, Bool_t recurseNonDerived) const
510{
511 treeNodeServerList(list,arg,kTRUE,kFALSE,kFALSE,recurseNonDerived) ;
512}
513
514
515////////////////////////////////////////////////////////////////////////////////
516/// Fill supplied list with nodes of the arg tree, following all server links,
517/// starting with ourself as top node.
518/// \param[in] list Output list
519/// \param[in] arg Start searching at this element of the tree.
520/// \param[in] doBranch Add branch nodes to the list.
521/// \param[in] doLeaf Add leaf nodes to the list.
522/// \param[in] valueOnly Only check if an element is a value server (no shape server).
523/// \param[in] recurseFundamental
524
525void RooAbsArg::treeNodeServerList(RooAbsCollection* list, const RooAbsArg* arg, Bool_t doBranch, Bool_t doLeaf, Bool_t valueOnly, Bool_t recurseFundamental) const
526{
527// if (arg==0) {
528// cout << "treeNodeServerList(" << GetName() << ") doBranch=" << (doBranch?"T":"F") << " doLeaf = " << (doLeaf?"T":"F") << " valueOnly=" << (valueOnly?"T":"F") << endl ;
529// }
530
531 if (!arg) {
532 list->reserve(10);
533 arg=this ;
534 }
535
536 // Decide if to add current node
537 if ((doBranch&&doLeaf) ||
538 (doBranch&&arg->isDerived()) ||
539 (doLeaf&&arg->isFundamental()&&(!(recurseFundamental&&arg->isDerived()))) ||
540 (doLeaf && !arg->isFundamental() && !arg->isDerived())) {
541
542 list->add(*arg,kTRUE) ;
543 }
544
545 // Recurse if current node is derived
546 if (arg->isDerived() && (!arg->isFundamental() || recurseFundamental)) {
547 for (const auto server : arg->_serverList) {
548
549 // Skip non-value server nodes if requested
550 Bool_t isValueSrv = server->_clientListValue.containsByNamePtr(arg);
551 if (valueOnly && !isValueSrv) {
552 continue ;
553 }
554 treeNodeServerList(list,server,doBranch,doLeaf,valueOnly,recurseFundamental) ;
555 }
556 }
557}
558
559
560////////////////////////////////////////////////////////////////////////////////
561/// Create a list of leaf nodes in the arg tree starting with
562/// ourself as top node that don't match any of the names of the variable list
563/// of the supplied data set (the dependents). The caller of this
564/// function is responsible for deleting the returned argset.
565/// The complement of this function is getObservables()
566
567RooArgSet* RooAbsArg::getParameters(const RooAbsData* set, Bool_t stripDisconnected) const
568{
569 return getParameters(set?set->get():0,stripDisconnected) ;
570}
571
572
573////////////////////////////////////////////////////////////////////////////////
574/// Add all parameters of the function and its daughters to `params`.
575/// \param[in] params Collection that stores all parameters. Add all new parameters to this.
576/// \param[in] nset Normalisation set (optional). If a value depends on this set, it's not a parameter.
577/// \param[in] stripDisconnected Passed on to getParametersHook().
578
579void RooAbsArg::addParameters(RooAbsCollection& params, const RooArgSet* nset, bool stripDisconnected) const
580{
581
582 RooArgSet nodeParamServers;
583 std::vector<RooAbsArg*> branchList;
584 for (const auto server : _serverList) {
585 if (server->isValueServer(*this)) {
586 if (server->isFundamental()) {
587 if (!nset || !server->dependsOn(*nset)) {
588 nodeParamServers.add(*server);
589 }
590 } else {
591 branchList.push_back(server);
592 }
593 }
594 }
595
596 // Allow pdf to strip parameters from list before adding it
597 getParametersHook(nset,&nodeParamServers,stripDisconnected) ;
598
599 // Add parameters of this node to the combined list
600 params.add(nodeParamServers,kTRUE) ;
601
602 // Now recurse into branch servers
603 std::sort(branchList.begin(), branchList.end());
604 const auto last = std::unique(branchList.begin(), branchList.end());
605 for (auto serverIt = branchList.begin(); serverIt < last; ++serverIt) {
606 (*serverIt)->addParameters(params, nset);
607 }
608}
609
610////////////////////////////////////////////////////////////////////////////////
611/// Obtain an estimate of the number of parameters of the function and its daughters.
612/// Calling `addParamters` for large functions (NLL) can cause many reallocations of
613/// `params` due to the recursive behaviour. This utility function aims to pre-compute
614/// the total number of parameters, so that enough memory is reserved.
615/// The estimate is not fully accurate (overestimate) as there is no equivalent to `getParametersHook`.
616/// \param[in] nset Normalisation set (optional). If a value depends on this set, it's not a parameter.
617
619{
620
621 std::size_t res = 0;
622 std::vector<RooAbsArg*> branchList;
623 for (const auto server : _serverList) {
624 if (server->isValueServer(*this)) {
625 if (server->isFundamental()) {
626 if (!nset || !server->dependsOn(*nset)) {
627 res++;
628 }
629 } else {
630 branchList.push_back(server);
631 }
632 }
633 }
634
635 // Now recurse into branch servers
636 std::sort(branchList.begin(), branchList.end());
637 const auto last = std::unique(branchList.begin(), branchList.end());
638 for (auto serverIt = branchList.begin(); serverIt < last; ++serverIt) {
639 res += (*serverIt)->getParametersSizeEstimate(nset);
640 }
641
642 return res;
643}
644
645////////////////////////////////////////////////////////////////////////////////
646/// Create a list of leaf nodes in the arg tree starting with
647/// ourself as top node that don't match any of the names the args in the
648/// supplied argset. The caller of this function is responsible
649/// for deleting the returned argset. The complement of this function
650/// is getObservables().
651
652RooArgSet* RooAbsArg::getParameters(const RooArgSet* observables, bool stripDisconnected) const {
653 auto * outputSet = new RooArgSet;
654 getParameters(observables, *outputSet, stripDisconnected);
655 return outputSet;
656}
657
658
659////////////////////////////////////////////////////////////////////////////////
660/// Fills a list with leaf nodes in the arg tree starting with
661/// ourself as top node that don't match any of the names the args in the
662/// supplied argset. Returns `true` only if something went wrong.
663/// The complement of this function is getObservables().
664/// \param[in] observables Set of leafs to ignore because they are observables and not parameters.
665/// \param[out] outputSet Output set.
666/// \param[in] stripDisconnected Allow pdf to strip parameters from list before adding it.
667
668bool RooAbsArg::getParameters(const RooArgSet* observables, RooArgSet& outputSet, bool stripDisconnected) const
669{
671
672 // Check for cached parameter set
673 if (_myws) {
674 auto nsetObs = getColonSeparatedNameString(observables ? *observables : RooArgSet());
675 const RooArgSet *paramSet = _myws->set(Form("CACHE_PARAMS_OF_PDF_%s_FOR_OBS_%s", GetName(), nsetObs.c_str()));
676 if (paramSet) {
677 outputSet.add(*paramSet);
678 return false;
679 }
680 }
681
682 outputSet.clear();
683 outputSet.setName("parameters");
684
685 RooArgList tempList;
686 // reserve all memory needed in one go
687 tempList.reserve(getParametersSizeEstimate(observables));
688
689 addParameters(tempList, observables, stripDisconnected);
690
691 // The adding from the list to the set has to be silent to not complain
692 // about duplicate parameters. After all, it's normal that parameters can
693 // appear in sifferent components of the model.
694 outputSet.add(tempList, /*silent=*/true);
695 outputSet.sort();
696
697 // Cache parameter set
698 if (_myws && outputSet.size() > 10) {
699 auto nsetObs = getColonSeparatedNameString(observables ? *observables : RooArgSet());
700 _myws->defineSetInternal(Form("CACHE_PARAMS_OF_PDF_%s_FOR_OBS_%s", GetName(), nsetObs.c_str()), outputSet);
701 }
702
703 return false;
704}
705
706
707////////////////////////////////////////////////////////////////////////////////
708/// Create a list of leaf nodes in the arg tree starting with
709/// ourself as top node that match any of the names of the variable list
710/// of the supplied data set (the dependents). The caller of this
711/// function is responsible for deleting the returned argset.
712/// The complement of this function is getParameters().
713
715{
716 if (!set) return new RooArgSet ;
717
718 return getObservables(set->get()) ;
719}
720
721
722////////////////////////////////////////////////////////////////////////////////
723/// Create a list of leaf nodes in the arg tree starting with
724/// ourself as top node that match any of the names the args in the
725/// supplied argset. The caller of this function is responsible
726/// for deleting the returned argset. The complement of this function
727/// is getParameters().
728
729RooArgSet* RooAbsArg::getObservables(const RooArgSet* dataList, bool valueOnly) const
730{
731 auto depList = new RooArgSet;
732 getObservables(dataList, *depList, valueOnly);
733 return depList;
734}
735
736
737////////////////////////////////////////////////////////////////////////////////
738/// Create a list of leaf nodes in the arg tree starting with
739/// ourself as top node that match any of the names the args in the
740/// supplied argset.
741/// Returns `true` only if something went wrong.
742/// The complement of this function is getParameters().
743/// \param[in] dataList Set of leaf nodes to match.
744/// \param[out] outputSet Output set.
745/// \param[in] valueOnly If this parameter is true, we only match leafs that
746/// depend on the value of any arg in `dataList`.
747
748bool RooAbsArg::getObservables(const RooAbsCollection* dataList, RooArgSet& outputSet, bool valueOnly) const
749{
750 outputSet.clear();
751 outputSet.setName("dependents");
752
753 if (!dataList) return false;
754
755 // Make iterator over tree leaf node list
756 RooArgSet leafList("leafNodeServerList") ;
757 treeNodeServerList(&leafList,0,kFALSE,kTRUE,valueOnly) ;
758
759 if (valueOnly) {
760 for (const auto arg : leafList) {
761 if (arg->dependsOnValue(*dataList) && arg->isLValue()) {
762 outputSet.add(*arg) ;
763 }
764 }
765 } else {
766 for (const auto arg : leafList) {
767 if (arg->dependsOn(*dataList) && arg->isLValue()) {
768 outputSet.add(*arg) ;
769 }
770 }
771 }
772
773 return false;
774}
775
776
777////////////////////////////////////////////////////////////////////////////////
778/// Create a RooArgSet with all components (branch nodes) of the
779/// expression tree headed by this object.
781{
782 TString name(GetName()) ;
783 name.Append("_components") ;
784
785 RooArgSet* set = new RooArgSet(name) ;
787
788 return set ;
789}
790
791
792
793////////////////////////////////////////////////////////////////////////////////
794/// Overloadable function in which derived classes can implement
795/// consistency checks of the variables. If this function returns
796/// true, indicating an error, the fitter or generator will abort.
797
799{
800 return kFALSE ;
801}
802
803
804////////////////////////////////////////////////////////////////////////////////
805/// Recursively call checkObservables on all nodes in the expression tree
806
808{
809 RooArgSet nodeList ;
810 treeNodeServerList(&nodeList) ;
811 RooFIter iter = nodeList.fwdIterator() ;
812
813 RooAbsArg* arg ;
814 Bool_t ret(kFALSE) ;
815 while((arg=iter.next())) {
816 if (arg->getAttribute("ServerDied")) {
817 coutE(LinkStateMgmt) << "RooAbsArg::recursiveCheckObservables(" << GetName() << "): ERROR: one or more servers of node "
818 << arg->GetName() << " no longer exists!" << endl ;
819 arg->Print("v") ;
820 ret = kTRUE ;
821 }
822 ret |= arg->checkObservables(nset) ;
823 }
824
825 return ret ;
826}
827
828
829////////////////////////////////////////////////////////////////////////////////
830/// Test whether we depend on (ie, are served by) any object in the
831/// specified collection. Uses the dependsOn(RooAbsArg&) member function.
832
833Bool_t RooAbsArg::dependsOn(const RooAbsCollection& serverList, const RooAbsArg* ignoreArg, Bool_t valueOnly) const
834{
835 // Test whether we depend on (ie, are served by) any object in the
836 // specified collection. Uses the dependsOn(RooAbsArg&) member function.
837
838 for (auto server : serverList) {
839 if (dependsOn(*server,ignoreArg,valueOnly)) {
840 return kTRUE;
841 }
842 }
843 return kFALSE;
844}
845
846
847////////////////////////////////////////////////////////////////////////////////
848/// Test whether we depend on (ie, are served by) the specified object.
849/// Note that RooAbsArg objects are considered equivalent if they have
850/// the same name.
851
852Bool_t RooAbsArg::dependsOn(const RooAbsArg& testArg, const RooAbsArg* ignoreArg, Bool_t valueOnly) const
853{
854 if (this==ignoreArg) return kFALSE ;
855
856 // First check if testArg is self
857 //if (!TString(testArg.GetName()).CompareTo(GetName())) return kTRUE ;
858 if (testArg.namePtr()==namePtr()) return kTRUE ;
859
860
861 // Next test direct dependence
862 RooAbsArg* foundServer = findServer(testArg) ;
863 if (foundServer) {
864
865 // Return true if valueOnly is FALSE or if server is value server, otherwise keep looking
866 if ( !valueOnly || foundServer->isValueServer(*this)) {
867 return kTRUE ;
868 }
869 }
870
871 // If not, recurse
872 for (const auto server : _serverList) {
873 if ( !valueOnly || server->isValueServer(*this)) {
874 if (server->dependsOn(testArg,ignoreArg,valueOnly)) {
875 return kTRUE ;
876 }
877 }
878 }
879
880 return kFALSE ;
881}
882
883
884
885////////////////////////////////////////////////////////////////////////////////
886/// Test if any of the nodes of tree are shared with that of the given tree
887
888Bool_t RooAbsArg::overlaps(const RooAbsArg& testArg, Bool_t valueOnly) const
889{
890 RooArgSet list("treeNodeList") ;
891 treeNodeServerList(&list) ;
892
893 return valueOnly ? testArg.dependsOnValue(list) : testArg.dependsOn(list) ;
894}
895
896
897
898////////////////////////////////////////////////////////////////////////////////
899/// Test if any of the dependents of the arg tree (as determined by getObservables)
900/// overlaps with those of the testArg.
901
903{
904 return observableOverlaps(dset->get(),testArg) ;
905}
906
907
908////////////////////////////////////////////////////////////////////////////////
909/// Test if any of the dependents of the arg tree (as determined by getObservables)
910/// overlaps with those of the testArg.
911
912Bool_t RooAbsArg::observableOverlaps(const RooArgSet* nset, const RooAbsArg& testArg) const
913{
914 RooArgSet* depList = getObservables(nset) ;
915 Bool_t ret = testArg.dependsOn(*depList) ;
916 delete depList ;
917 return ret ;
918}
919
920
921
922////////////////////////////////////////////////////////////////////////////////
923/// Mark this object as having changed its value, and propagate this status
924/// change to all of our clients. If the object is not in automatic dirty
925/// state propagation mode, this call has no effect.
926
928{
929 if (_operMode!=Auto || _inhibitDirty) return ;
930
931 // Handle no-propagation scenarios first
932 if (_clientListValue.size() == 0) {
934 return ;
935 }
936
937 // Cyclical dependency interception
938 if (source==0) {
939 source=this ;
940 } else if (source==this) {
941 // Cyclical dependency, abort
942 coutE(LinkStateMgmt) << "RooAbsArg::setValueDirty(" << GetName()
943 << "): cyclical dependency detected, source = " << source->GetName() << endl ;
944 //assert(0) ;
945 return ;
946 }
947
948 // Propagate dirty flag to all clients if this is a down->up transition
949 if (_verboseDirty) {
950 cxcoutD(LinkStateMgmt) << "RooAbsArg::setValueDirty(" << (source?source->GetName():"self") << "->" << GetName() << "," << this
951 << "): dirty flag " << (_valueDirty?"already ":"") << "raised" << endl ;
952 }
953
955
956
957 for (auto client : _clientListValue) {
958 client->setValueDirty(source) ;
959 }
960
961
962}
963
964
965////////////////////////////////////////////////////////////////////////////////
966/// Mark this object as having changed its shape, and propagate this status
967/// change to all of our clients.
968
970{
971 if (_verboseDirty) {
972 cxcoutD(LinkStateMgmt) << "RooAbsArg::setShapeDirty(" << GetName()
973 << "): dirty flag " << (_shapeDirty?"already ":"") << "raised" << endl ;
974 }
975
976 if (_clientListShape.empty()) {
978 return ;
979 }
980
981 // Set 'dirty' shape state for this object and propagate flag to all its clients
982 if (source==0) {
983 source=this ;
984 } else if (source==this) {
985 // Cyclical dependency, abort
986 coutE(LinkStateMgmt) << "RooAbsArg::setShapeDirty(" << GetName()
987 << "): cyclical dependency detected" << endl ;
988 return ;
989 }
990
991 // Propagate dirty flag to all clients if this is a down->up transition
993
994 for (auto client : _clientListShape) {
995 client->setShapeDirty(source) ;
996 client->setValueDirty(source) ;
997 }
998
999}
1000
1001
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// Replace all direct servers of this object with the new servers in `newServerList`.
1005/// This substitutes objects that we receive values from with new objects that have the same name.
1006/// \see recursiveRedirectServers() Use recursive version if servers that are only indirectly serving this object should be replaced as well.
1007/// \see redirectServers() If only the direct servers of an object need to be replaced.
1008///
1009/// Note that changing the types of objects is generally allowed, but can be wrong if the interface of an object changes.
1010/// For example, one can reparametrise a model by substituting a variable with a function:
1011/// \f[
1012/// f(x\, |\, a) = a \cdot x \rightarrow f(x\, |\, b) = (2.1 \cdot b) \cdot x
1013/// \f]
1014/// If an object, however, expects a PDF, and this is substituted with a function that isn't normalised, wrong results might be obtained
1015/// or it might even crash the program. The types of the objects being substituted are not checked.
1016///
1017/// \param[in] newSetOrig Set of new servers that should be used instead of the current servers.
1018/// \param[in] mustReplaceAll A warning is printed and error status is returned if not all servers could be
1019/// substituted successfully.
1020/// \param[in] nameChange If false, an object named "x" is only replaced with an object also named "x" in `newSetOrig`.
1021/// If the object in `newSet` is called differently, set `nameChange` to true and use setAttribute() on the x object:
1022/// ```
1023/// objectToReplaceX.setAttribute("ORIGNAME:x")
1024/// ```
1025/// Now, the renamed object will be selected based on the attribute "ORIGNAME:<name>".
1026/// \param[in] isRecursionStep Internal switch used when called from recursiveRedirectServers().
1027Bool_t RooAbsArg::redirectServers(const RooAbsCollection& newSetOrig, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursionStep)
1028{
1029 // Trivial case, no servers
1030 if (_serverList.empty()) return kFALSE ;
1031 if (newSetOrig.getSize()==0) return kFALSE ;
1032
1033 // Strip any non-matching removal nodes from newSetOrig
1034 RooAbsCollection* newSet ;
1035
1036 if (nameChange) {
1037 newSet = new RooArgSet ;
1038 for (auto arg : newSetOrig) {
1039
1040 if (string("REMOVAL_DUMMY")==arg->GetName()) {
1041
1042 if (arg->getAttribute("REMOVE_ALL")) {
1043 newSet->add(*arg) ;
1044 } else if (arg->getAttribute(Form("REMOVE_FROM_%s",getStringAttribute("ORIGNAME")))) {
1045 newSet->add(*arg) ;
1046 }
1047 } else {
1048 newSet->add(*arg) ;
1049 }
1050 }
1051 } else {
1052 newSet = (RooAbsCollection*) &newSetOrig ;
1053 }
1054
1055 // Replace current servers with new servers with the same name from the given list
1056 Bool_t ret(kFALSE) ;
1057
1058 //Copy original server list to not confuse the iterator while deleting
1059 std::vector<RooAbsArg*> origServerList, origServerValue, origServerShape;
1060 auto origSize = _serverList.size();
1061 origServerList.reserve(origSize);
1062 origServerValue.reserve(origSize);
1063
1064 for (const auto oldServer : _serverList) {
1065 origServerList.push_back(oldServer) ;
1066
1067 // Retrieve server side link state information
1068 if (oldServer->_clientListValue.containsByNamePtr(this)) {
1069 origServerValue.push_back(oldServer) ;
1070 }
1071 if (oldServer->_clientListShape.containsByNamePtr(this)) {
1072 origServerShape.push_back(oldServer) ;
1073 }
1074 }
1075
1076 // Delete all previously registered servers
1077 for (auto oldServer : origServerList) {
1078
1079 RooAbsArg * newServer= oldServer->findNewServer(*newSet, nameChange);
1080
1081 if (newServer && _verboseDirty) {
1082 cxcoutD(LinkStateMgmt) << "RooAbsArg::redirectServers(" << (void*)this << "," << GetName() << "): server " << oldServer->GetName()
1083 << " redirected from " << oldServer << " to " << newServer << endl ;
1084 }
1085
1086 if (!newServer) {
1087 if (mustReplaceAll) {
1088 coutE(LinkStateMgmt) << "RooAbsArg::redirectServers(" << (void*)this << "," << GetName() << "): server " << oldServer->GetName()
1089 << " (" << (void*)oldServer << ") not redirected" << (nameChange?"[nameChange]":"") << endl ;
1090 ret = kTRUE ;
1091 }
1092 continue ;
1093 }
1094
1095 auto findByNamePtr = [&oldServer](const RooAbsArg * item) {
1096 return oldServer->namePtr() == item->namePtr();
1097 };
1098 bool propValue = std::any_of(origServerValue.begin(), origServerValue.end(), findByNamePtr);
1099 bool propShape = std::any_of(origServerShape.begin(), origServerShape.end(), findByNamePtr);
1100
1101 if (newServer != this) {
1102 replaceServer(*oldServer,*newServer,propValue,propShape) ;
1103 }
1104 }
1105
1106
1107 setValueDirty() ;
1108 setShapeDirty() ;
1109
1110 // Process the proxies
1111 for (int i=0 ; i<numProxies() ; i++) {
1112 RooAbsProxy* p = getProxy(i) ;
1113 if (!p) continue ;
1114 Bool_t ret2 = p->changePointer(*newSet,nameChange,kFALSE) ;
1115
1116 if (mustReplaceAll && !ret2) {
1117 auto ap = dynamic_cast<const RooArgProxy*>(p);
1118 coutE(LinkStateMgmt) << "RooAbsArg::redirectServers(" << GetName()
1119 << "): ERROR, proxy '" << p->name()
1120 << "' with arg '" << (ap ? ap->absArg()->GetName() : "<could not cast>") << "' could not be adjusted" << endl;
1121 ret = kTRUE ;
1122 }
1123 }
1124
1125
1126 // Optional subclass post-processing
1127 for (Int_t i=0 ;i<numCaches() ; i++) {
1128 ret |= getCache(i)->redirectServersHook(*newSet,mustReplaceAll,nameChange,isRecursionStep) ;
1129 }
1130 ret |= redirectServersHook(*newSet,mustReplaceAll,nameChange,isRecursionStep) ;
1131
1132 if (nameChange) {
1133 delete newSet ;
1134 }
1135
1136 return ret ;
1137}
1138
1139////////////////////////////////////////////////////////////////////////////////
1140/// Find the new server in the specified set that matches the old server.
1141///
1142/// \param[in] newSet Search this set by name for a new server.
1143/// \param[in] If true, search for an item with the bool attribute "ORIGNAME:<oldName>" set.
1144/// Use `<object>.setAttribute("ORIGNAME:<oldName>")` to set this attribute.
1145/// \return Pointer to the new server or `nullptr` if there's no unique match.
1147{
1148 RooAbsArg *newServer = 0;
1149 if (!nameChange) {
1150 newServer = newSet.find(*this) ;
1151 }
1152 else {
1153 // Name changing server redirect:
1154 // use 'ORIGNAME:<oldName>' attribute instead of name of new server
1155 TString nameAttrib("ORIGNAME:") ;
1156 nameAttrib.Append(GetName()) ;
1157
1158 RooArgSet* tmp = (RooArgSet*) newSet.selectByAttrib(nameAttrib,kTRUE) ;
1159 if(0 != tmp) {
1160
1161 // Check if any match was found
1162 if (tmp->getSize()==0) {
1163 delete tmp ;
1164 return 0 ;
1165 }
1166
1167 // Check if match is unique
1168 if(tmp->getSize()>1) {
1169 coutF(LinkStateMgmt) << "RooAbsArg::redirectServers(" << GetName() << "): FATAL Error, " << tmp->getSize() << " servers with "
1170 << nameAttrib << " attribute" << endl ;
1171 tmp->Print("v") ;
1172 assert(0) ;
1173 }
1174
1175 // use the unique element in the set
1176 newServer= tmp->first();
1177 delete tmp ;
1178 }
1179 }
1180 return newServer;
1181}
1182
1183
1184////////////////////////////////////////////////////////////////////////////////
1185/// Recursively replace all servers with the new servers in `newSet`.
1186/// This substitutes objects that we receive values from (also indirectly through other objects) with new objects that have the same name.
1187///
1188/// *Copied from redirectServers:*
1189///
1190/// \copydetails RooAbsArg::redirectServers
1191Bool_t RooAbsArg::recursiveRedirectServers(const RooAbsCollection& newSet, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t recurseInNewSet)
1192{
1193 // Cyclic recursion protection
1194 static std::set<const RooAbsArg*> callStack;
1195 {
1196 std::set<const RooAbsArg*>::iterator it = callStack.lower_bound(this);
1197 if (it != callStack.end() && this == *it) {
1198 return kFALSE;
1199 } else {
1200 callStack.insert(it, this);
1201 }
1202 }
1203
1204 // Do not recurse into newset if not so specified
1205// if (!recurseInNewSet && newSet.contains(*this)) {
1206// return kFALSE ;
1207// }
1208
1209
1210 // Apply the redirectServers function recursively on all branch nodes in this argument tree.
1211 Bool_t ret(kFALSE) ;
1212
1213 cxcoutD(LinkStateMgmt) << "RooAbsArg::recursiveRedirectServers(" << this << "," << GetName() << ") newSet = " << newSet << " mustReplaceAll = "
1214 << (mustReplaceAll?"T":"F") << " nameChange = " << (nameChange?"T":"F") << " recurseInNewSet = " << (recurseInNewSet?"T":"F") << endl ;
1215
1216 // Do redirect on self (identify operation as recursion step)
1217 ret |= redirectServers(newSet,mustReplaceAll,nameChange,kTRUE) ;
1218
1219 // Do redirect on servers
1220 for (const auto server : _serverList){
1221 ret |= server->recursiveRedirectServers(newSet,mustReplaceAll,nameChange,recurseInNewSet) ;
1222 }
1223
1224 callStack.erase(this);
1225 return ret ;
1226}
1227
1228
1229
1230////////////////////////////////////////////////////////////////////////////////
1231/// Register an RooArgProxy in the proxy list. This function is called by owned
1232/// proxies upon creation. After registration, this arg wil forward pointer
1233/// changes from serverRedirects and updates in cached normalization sets
1234/// to the proxies immediately after they occur. The proxied argument is
1235/// also added as value and/or shape server
1236
1238{
1239 // Every proxy can be registered only once
1240 if (_proxyList.FindObject(&proxy)) {
1241 coutE(LinkStateMgmt) << "RooAbsArg::registerProxy(" << GetName() << "): proxy named "
1242 << proxy.GetName() << " for arg " << proxy.absArg()->GetName()
1243 << " already registered" << endl ;
1244 return ;
1245 }
1246
1247// cout << (void*)this << " " << GetName() << ": registering proxy "
1248// << (void*)&proxy << " with name " << proxy.name() << " in mode "
1249// << (proxy.isValueServer()?"V":"-") << (proxy.isShapeServer()?"S":"-") << endl ;
1250
1251 // Register proxied object as server
1252 if (proxy.absArg()) {
1253 addServer(*proxy.absArg(),proxy.isValueServer(),proxy.isShapeServer()) ;
1254 }
1255
1256 // Register proxy itself
1257 _proxyList.Add(&proxy) ;
1258 _proxyListCache.isDirty = true;
1259}
1260
1261
1262////////////////////////////////////////////////////////////////////////////////
1263/// Remove proxy from proxy list. This functions is called by owned proxies
1264/// upon their destruction.
1265
1267{
1268 _proxyList.Remove(&proxy) ;
1270 _proxyListCache.isDirty = true;
1271}
1272
1273
1274
1275////////////////////////////////////////////////////////////////////////////////
1276/// Register an RooSetProxy in the proxy list. This function is called by owned
1277/// proxies upon creation. After registration, this arg wil forward pointer
1278/// changes from serverRedirects and updates in cached normalization sets
1279/// to the proxies immediately after they occur.
1280
1282{
1283 // Every proxy can be registered only once
1284 if (_proxyList.FindObject(&proxy)) {
1285 coutE(LinkStateMgmt) << "RooAbsArg::registerProxy(" << GetName() << "): proxy named "
1286 << proxy.GetName() << " already registered" << endl ;
1287 return ;
1288 }
1289
1290 // Register proxy itself
1291 _proxyList.Add(&proxy) ;
1292 _proxyListCache.isDirty = true;
1293}
1294
1295
1296
1297////////////////////////////////////////////////////////////////////////////////
1298/// Remove proxy from proxy list. This functions is called by owned proxies
1299/// upon their destruction.
1300
1302{
1303 _proxyList.Remove(&proxy) ;
1305 _proxyListCache.isDirty = true;
1306}
1307
1308
1309
1310////////////////////////////////////////////////////////////////////////////////
1311/// Register an RooListProxy in the proxy list. This function is called by owned
1312/// proxies upon creation. After registration, this arg wil forward pointer
1313/// changes from serverRedirects and updates in cached normalization sets
1314/// to the proxies immediately after they occur.
1315
1317{
1318 // Every proxy can be registered only once
1319 if (_proxyList.FindObject(&proxy)) {
1320 coutE(LinkStateMgmt) << "RooAbsArg::registerProxy(" << GetName() << "): proxy named "
1321 << proxy.GetName() << " already registered" << endl ;
1322 return ;
1323 }
1324
1325 // Register proxy itself
1326 Int_t nProxyOld = _proxyList.GetEntries() ;
1327 _proxyList.Add(&proxy) ;
1328 _proxyListCache.isDirty = true;
1329 if (_proxyList.GetEntries()!=nProxyOld+1) {
1330 cout << "RooAbsArg::registerProxy(" << GetName() << ") proxy registration failure! nold=" << nProxyOld << " nnew=" << _proxyList.GetEntries() << endl ;
1331 }
1332}
1333
1334
1335
1336////////////////////////////////////////////////////////////////////////////////
1337/// Remove proxy from proxy list. This functions is called by owned proxies
1338/// upon their destruction.
1339
1341{
1342 _proxyList.Remove(&proxy) ;
1344 _proxyListCache.isDirty = true;
1345}
1346
1347
1348
1349////////////////////////////////////////////////////////////////////////////////
1350/// Return the nth proxy from the proxy list.
1351
1353{
1354 // Cross cast: proxy list returns TObject base pointer, we need
1355 // a RooAbsProxy base pointer. C++ standard requires
1356 // a dynamic_cast for this.
1357 return dynamic_cast<RooAbsProxy*> (_proxyList.At(index)) ;
1358}
1359
1360
1361
1362////////////////////////////////////////////////////////////////////////////////
1363/// Return the number of registered proxies.
1364
1366{
1367 return _proxyList.GetEntriesFast();
1368}
1369
1370
1371
1372////////////////////////////////////////////////////////////////////////////////
1373/// Forward a change in the cached normalization argset
1374/// to all the registered proxies.
1375
1377{
1379 // First time we loop over proxies: cache the results to avoid future
1380 // costly dynamic_casts
1381 _proxyListCache.cache.clear();
1382 for (int i=0 ; i<numProxies() ; i++) {
1383 RooAbsProxy* p = getProxy(i) ;
1384 if (!p) continue ;
1385 _proxyListCache.cache.push_back(p);
1386 }
1387 _proxyListCache.isDirty = false;
1388 }
1389
1390 for ( auto& p : _proxyListCache.cache ) {
1391 p->changeNormSet(nset);
1392 }
1393}
1394
1395
1396
1397////////////////////////////////////////////////////////////////////////////////
1398/// Overloadable function for derived classes to implement
1399/// attachment as branch to a TTree
1400
1402{
1403 coutE(Contents) << "RooAbsArg::attachToTree(" << GetName()
1404 << "): Cannot be attached to a TTree" << endl ;
1405}
1406
1407
1408
1409////////////////////////////////////////////////////////////////////////////////
1410/// WVE (08/21/01) Probably obsolete now
1411
1413{
1414 return kTRUE ;
1415}
1416
1417
1418
1419
1420////////////////////////////////////////////////////////////////////////////////
1421/// Print object name
1422
1423void RooAbsArg::printName(ostream& os) const
1424{
1425 os << GetName() ;
1426}
1427
1428
1429
1430////////////////////////////////////////////////////////////////////////////////
1431/// Print object title
1432
1433void RooAbsArg::printTitle(ostream& os) const
1434{
1435 os << GetTitle() ;
1436}
1437
1438
1439
1440////////////////////////////////////////////////////////////////////////////////
1441/// Print object class name
1442
1443void RooAbsArg::printClassName(ostream& os) const
1444{
1445 os << IsA()->GetName() ;
1446}
1447
1448
1449void RooAbsArg::printAddress(ostream& os) const
1450{
1451 // Print addrss of this RooAbsArg
1452 os << this ;
1453}
1454
1455
1456
1457////////////////////////////////////////////////////////////////////////////////
1458/// Print object arguments, ie its proxies
1459
1460void RooAbsArg::printArgs(ostream& os) const
1461{
1462 // Print nothing if there are no dependencies
1463 if (numProxies()==0) return ;
1464
1465 os << "[ " ;
1466 for (Int_t i=0 ; i<numProxies() ; i++) {
1467 RooAbsProxy* p = getProxy(i) ;
1468 if (p==0) continue ;
1469 if (!TString(p->name()).BeginsWith("!")) {
1470 p->print(os) ;
1471 os << " " ;
1472 }
1473 }
1474 printMetaArgs(os) ;
1475 os << "]" ;
1476}
1477
1478
1479
1480////////////////////////////////////////////////////////////////////////////////
1481/// Define default contents to print
1482
1484{
1485 return kName|kClassName|kValue|kArgs ;
1486}
1487
1488
1489
1490////////////////////////////////////////////////////////////////////////////////
1491/// Implement multi-line detailed printing
1492
1493void RooAbsArg::printMultiline(ostream& os, Int_t /*contents*/, Bool_t /*verbose*/, TString indent) const
1494{
1495 os << indent << "--- RooAbsArg ---" << endl;
1496 // dirty state flags
1497 os << indent << " Value State: " ;
1498 switch(_operMode) {
1499 case ADirty: os << "FORCED DIRTY" ; break ;
1500 case AClean: os << "FORCED clean" ; break ;
1501 case Auto: os << (isValueDirty() ? "DIRTY":"clean") ; break ;
1502 }
1503 os << endl
1504 << indent << " Shape State: " << (isShapeDirty() ? "DIRTY":"clean") << endl;
1505 // attribute list
1506 os << indent << " Attributes: " ;
1507 printAttribList(os) ;
1508 os << endl ;
1509 // our memory address (for x-referencing with client addresses of other args)
1510 os << indent << " Address: " << (void*)this << endl;
1511 // client list
1512 os << indent << " Clients: " << endl;
1513 for (const auto client : _clientList) {
1514 os << indent << " (" << (void*)client << ","
1515 << (_clientListValue.containsByNamePtr(client)?"V":"-")
1516 << (_clientListShape.containsByNamePtr(client)?"S":"-")
1517 << ") " ;
1518 client->printStream(os,kClassName|kTitle|kName,kSingleLine);
1519 }
1520
1521 // server list
1522 os << indent << " Servers: " << endl;
1523 for (const auto server : _serverList) {
1524 os << indent << " (" << (void*)server << ","
1525 << (server->_clientListValue.containsByNamePtr(this)?"V":"-")
1526 << (server->_clientListShape.containsByNamePtr(this)?"S":"-")
1527 << ") " ;
1528 server->printStream(os,kClassName|kName|kTitle,kSingleLine);
1529 }
1530
1531 // proxy list
1532 os << indent << " Proxies: " << endl ;
1533 for (int i=0 ; i<numProxies() ; i++) {
1534 RooAbsProxy* proxy=getProxy(i) ;
1535 if (!proxy) continue ;
1536 if (proxy->IsA()->InheritsFrom(RooArgProxy::Class())) {
1537 os << indent << " " << proxy->name() << " -> " ;
1538 RooAbsArg* parg = ((RooArgProxy*)proxy)->absArg() ;
1539 if (parg) {
1540 parg->printStream(os,kName,kSingleLine) ;
1541 } else {
1542 os << " (empty)" << endl ; ;
1543 }
1544 } else if(dynamic_cast<RooAbsCollection*>(proxy)) {
1545 os << indent << " " << proxy->name() << " -> " ;
1546 os << endl ;
1547 TString moreIndent(indent) ;
1548 moreIndent.Append(" ") ;
1549 dynamic_cast<RooAbsCollection*>(proxy)->printStream(os,kName,kStandard,moreIndent.Data()) ;
1550 }
1551 }
1552}
1553
1554
1555////////////////////////////////////////////////////////////////////////////////
1556/// Print object tree structure
1557
1558void RooAbsArg::printTree(ostream& os, TString /*indent*/) const
1559{
1560 const_cast<RooAbsArg*>(this)->printCompactTree(os) ;
1561}
1562
1563
1564////////////////////////////////////////////////////////////////////////////////
1565/// Ostream operator
1566
1567ostream& operator<<(ostream& os, RooAbsArg const& arg)
1568{
1569 arg.writeToStream(os,kTRUE) ;
1570 return os ;
1571}
1572
1573////////////////////////////////////////////////////////////////////////////////
1574/// Istream operator
1575
1576istream& operator>>(istream& is, RooAbsArg &arg)
1577{
1578 arg.readFromStream(is,kTRUE,kFALSE) ;
1579 return is ;
1580}
1581
1582////////////////////////////////////////////////////////////////////////////////
1583/// Print the attribute list
1584
1585void RooAbsArg::printAttribList(ostream& os) const
1586{
1587 set<string>::const_iterator iter = _boolAttrib.begin() ;
1588 Bool_t first(kTRUE) ;
1589 while (iter != _boolAttrib.end()) {
1590 os << (first?" [":",") << *iter ;
1591 first=kFALSE ;
1592 ++iter ;
1593 }
1594 if (!first) os << "] " ;
1595}
1596
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Bind this node to objects in `set`.
1600/// Search the set for objects that have the same name as our servers, and
1601/// attach ourselves to those. After this operation, this node is computing its
1602/// values based on the new servers. This can be used to e.g. read values from
1603// a dataset.
1604
1605
1607{
1608 RooArgSet branches;
1609 branchNodeServerList(&branches,0,true);
1610
1611 for(auto const& branch : branches) {
1612 branch->redirectServers(set,false,false);
1613 }
1614}
1615
1616
1617////////////////////////////////////////////////////////////////////////////////
1618/// Replace server nodes with names matching the dataset variable names
1619/// with those data set variables, making this PDF directly dependent on the dataset.
1620
1622{
1623 attachArgs(*data.get());
1624}
1625
1626
1627////////////////////////////////////////////////////////////////////////////////
1628/// Replace server nodes with names matching the dataset variable names
1629/// with those data set variables, making this PDF directly dependent on the dataset
1630
1632{
1633 attachArgs(*dstore.get());
1634}
1635
1636
1637////////////////////////////////////////////////////////////////////////////////
1638/// Utility function used by TCollection::Sort to compare contained TObjects
1639/// We implement comparison by name, resulting in alphabetical sorting by object name.
1640
1642{
1643 return strcmp(GetName(),other->GetName()) ;
1644}
1645
1646
1647
1648////////////////////////////////////////////////////////////////////////////////
1649/// Print information about current value dirty state information.
1650/// If depth flag is true, information is recursively printed for
1651/// all nodes in this arg tree.
1652
1654{
1655 if (depth) {
1656
1657 RooArgSet branchList ;
1658 branchNodeServerList(&branchList) ;
1659 RooFIter bIter = branchList.fwdIterator() ;
1660 RooAbsArg* branch ;
1661 while((branch=bIter.next())) {
1662 branch->printDirty(kFALSE) ;
1663 }
1664
1665 } else {
1666 cout << GetName() << " : " ;
1667 switch (_operMode) {
1668 case AClean: cout << "FORCED clean" ; break ;
1669 case ADirty: cout << "FORCED DIRTY" ; break ;
1670 case Auto: cout << "Auto " << (isValueDirty()?"DIRTY":"clean") ;
1671 }
1672 cout << endl ;
1673 }
1674}
1675
1676
1677////////////////////////////////////////////////////////////////////////////////
1678/// Activate cache mode optimization with given definition of observables.
1679/// The cache operation mode of all objects in the expression tree will
1680/// modified such that all nodes that depend directly or indirectly on
1681/// any of the listed observables will be set to ADirty, as they are
1682/// expected to change every time. This save change tracking overhead for
1683/// nodes that are a priori known to change every time
1684
1686{
1687 RooLinkedList proc;
1688 RooArgSet opt ;
1689 optimizeCacheMode(observables,opt,proc) ;
1690
1691 coutI(Optimization) << "RooAbsArg::optimizeCacheMode(" << GetName() << ") nodes " << opt << " depend on observables, "
1692 << "changing cache operation mode from change tracking to unconditional evaluation" << endl ;
1693}
1694
1695
1696////////////////////////////////////////////////////////////////////////////////
1697/// Activate cache mode optimization with given definition of observables.
1698/// The cache operation mode of all objects in the expression tree will
1699/// modified such that all nodes that depend directly or indirectly on
1700/// any of the listed observables will be set to ADirty, as they are
1701/// expected to change every time. This save change tracking overhead for
1702/// nodes that are a priori known to change every time
1703
1704void RooAbsArg::optimizeCacheMode(const RooArgSet& observables, RooArgSet& optimizedNodes, RooLinkedList& processedNodes)
1705{
1706 // Optimization applies only to branch nodes, not to leaf nodes
1707 if (!isDerived()) {
1708 return ;
1709 }
1710
1711
1712 // Terminate call if this node was already processed (tree structure may be cyclical)
1713 // LM : RooLinkedList::findArg looks by name and not but by object pointer,
1714 // should one use RooLinkedList::FindObject (look byt pointer) instead of findArg when
1715 // tree contains nodes with the same name ?
1716 // Add an info message if the require node does not exist but a different node already exists with same name
1717
1718 if (processedNodes.FindObject(this))
1719 return;
1720
1721 // check if findArgs returns something different (i.e. a different node with same name) when
1722 // this node has not been processed (FindObject returns a null pointer)
1723 auto obj = processedNodes.findArg(this);
1724 assert(obj != this); // obj == this cannot happen
1725 if (obj)
1726 // here for nodes with duplicate names
1727 cxcoutI(Optimization) << "RooAbsArg::optimizeCacheMode(" << GetName()
1728 << " node " << this << " exists already as " << obj << " but with the SAME name !" << endl;
1729
1730 processedNodes.Add(this);
1731
1732 // Set cache mode operator to 'AlwaysDirty' if we depend on any of the given observables
1733 if (dependsOnValue(observables)) {
1734
1735 if (dynamic_cast<RooRealIntegral*>(this)) {
1736 cxcoutI(Integration) << "RooAbsArg::optimizeCacheMode(" << GetName() << ") integral depends on value of one or more observables and will be evaluated for every event" << endl ;
1737 }
1738 optimizedNodes.add(*this,kTRUE) ;
1739 if (operMode()==AClean) {
1740 } else {
1741 setOperMode(ADirty,kTRUE) ; // WVE propagate flag recursively to top of tree
1742 }
1743 } else {
1744 }
1745 // Process any RooAbsArgs contained in any of the caches of this object
1746 for (Int_t i=0 ;i<numCaches() ; i++) {
1747 getCache(i)->optimizeCacheMode(observables,optimizedNodes,processedNodes) ;
1748 }
1749
1750 // Forward calls to all servers
1751 for (const auto server : _serverList) {
1752 server->optimizeCacheMode(observables,optimizedNodes,processedNodes) ;
1753 }
1754
1755}
1756
1757////////////////////////////////////////////////////////////////////////////////
1758/// Find branch nodes with all-constant parameters, and add them to the list of
1759/// nodes that can be cached with a dataset in a test statistic calculation
1760
1762{
1763 RooLinkedList proc ;
1764 Bool_t ret = findConstantNodes(observables,cacheList,proc) ;
1765
1766 // If node can be optimized and hasn't been identified yet, add it to the list
1767 coutI(Optimization) << "RooAbsArg::findConstantNodes(" << GetName() << "): components "
1768 << cacheList << " depend exclusively on constant parameters and will be precalculated and cached" << endl ;
1769
1770 return ret ;
1771}
1772
1773
1774
1775////////////////////////////////////////////////////////////////////////////////
1776/// Find branch nodes with all-constant parameters, and add them to the list of
1777/// nodes that can be cached with a dataset in a test statistic calculation
1778
1779Bool_t RooAbsArg::findConstantNodes(const RooArgSet& observables, RooArgSet& cacheList, RooLinkedList& processedNodes)
1780{
1781 // Caching only applies to branch nodes
1782 if (!isDerived()) {
1783 return kFALSE;
1784 }
1785
1786 // Terminate call if this node was already processed (tree structure may be cyclical)
1787 if (processedNodes.findArg(this)) {
1788 return kFALSE ;
1789 } else {
1790 processedNodes.Add(this) ;
1791 }
1792
1793 // Check if node depends on any non-constant parameter
1794 Bool_t canOpt(kTRUE) ;
1795 RooArgSet* paramSet = getParameters(observables) ;
1796 RooFIter iter = paramSet->fwdIterator() ;
1797 RooAbsArg* param ;
1798 while((param = iter.next())) {
1799 if (!param->isConstant()) {
1800 canOpt=kFALSE ;
1801 break ;
1802 }
1803 }
1804 delete paramSet ;
1805
1806
1807 if (getAttribute("NeverConstant")) {
1808 canOpt = kFALSE ;
1809 }
1810
1811 if (canOpt) {
1812 setAttribute("ConstantExpression") ;
1813 }
1814
1815 // If yes, list node eligible for caching, if not test nodes one level down
1816 if (canOpt||getAttribute("CacheAndTrack")) {
1817
1818 if (!cacheList.find(*this) && dependsOnValue(observables) && !observables.find(*this) ) {
1819
1820 // Add to cache list
1821 cxcoutD(Optimization) << "RooAbsArg::findConstantNodes(" << GetName() << ") adding self to list of constant nodes" << endl ;
1822
1823 if (canOpt) setAttribute("ConstantExpressionCached") ;
1824 cacheList.add(*this,kFALSE) ;
1825 }
1826 }
1827
1828 if (!canOpt) {
1829
1830 // If not, see if next level down can be cached
1831 for (const auto server : _serverList) {
1832 if (server->isDerived()) {
1833 server->findConstantNodes(observables,cacheList,processedNodes) ;
1834 }
1835 }
1836 }
1837
1838 // Forward call to all cached contained in current object
1839 for (Int_t i=0 ;i<numCaches() ; i++) {
1840 getCache(i)->findConstantNodes(observables,cacheList,processedNodes) ;
1841 }
1842
1843 return kFALSE ;
1844}
1845
1846
1847
1848
1849////////////////////////////////////////////////////////////////////////////////
1850/// Interface function signaling a request to perform constant term
1851/// optimization. This default implementation takes no action other than to
1852/// forward the calls to all servers
1853
1855{
1856 for (const auto server : _serverList) {
1857 server->constOptimizeTestStatistic(opcode,doAlsoTrackingOpt) ;
1858 }
1859}
1860
1861
1862////////////////////////////////////////////////////////////////////////////////
1863/// Change cache operation mode to given mode. If recurseAdirty
1864/// is true, then a mode change to AlwaysDirty will automatically
1865/// be propagated recursively to all client nodes
1866
1867void RooAbsArg::setOperMode(OperMode mode, Bool_t recurseADirty)
1868{
1869 // Prevent recursion loops
1870 if (mode==_operMode) return ;
1871
1872 _operMode = mode ;
1873 _fast = ((mode==AClean) || dynamic_cast<RooRealVar*>(this)!=0 || dynamic_cast<RooConstVar*>(this)!=0 ) ;
1874 for (Int_t i=0 ;i<numCaches() ; i++) {
1875 getCache(i)->operModeHook() ;
1876 }
1877 operModeHook() ;
1878
1879 // Propagate to all clients
1880 if (mode==ADirty && recurseADirty) {
1881 for (auto clientV : _clientListValue) {
1882 clientV->setOperMode(mode) ;
1883 }
1884 }
1885}
1886
1887
1888////////////////////////////////////////////////////////////////////////////////
1889/// Print tree structure of expression tree on stdout, or to file if filename is specified.
1890/// If namePat is not "*", only nodes with names matching the pattern will be printed.
1891/// The client argument is used in recursive calls to properly display the value or shape nature
1892/// of the client-server links. It should be zero in calls initiated by users.
1893
1894void RooAbsArg::printCompactTree(const char* indent, const char* filename, const char* namePat, RooAbsArg* client)
1895{
1896 if (filename) {
1897 ofstream ofs(filename) ;
1898 printCompactTree(ofs,indent,namePat,client) ;
1899 } else {
1900 printCompactTree(cout,indent,namePat,client) ;
1901 }
1902}
1903
1904
1905////////////////////////////////////////////////////////////////////////////////
1906/// Print tree structure of expression tree on given ostream.
1907/// If namePat is not "*", only nodes with names matching the pattern will be printed.
1908/// The client argument is used in recursive calls to properly display the value or shape nature
1909/// of the client-server links. It should be zero in calls initiated by users.
1910
1911void RooAbsArg::printCompactTree(ostream& os, const char* indent, const char* namePat, RooAbsArg* client)
1912{
1913 if ( !namePat || TString(GetName()).Contains(namePat)) {
1914 os << indent << this ;
1915 if (client) {
1916 os << "/" ;
1917 if (isValueServer(*client)) os << "V" ; else os << "-" ;
1918 if (isShapeServer(*client)) os << "S" ; else os << "-" ;
1919 }
1920 os << " " ;
1921
1922 os << IsA()->GetName() << "::" << GetName() << " = " ;
1923 printValue(os) ;
1924
1925 if (!_serverList.empty()) {
1926 switch(operMode()) {
1927 case Auto: os << " [Auto," << (isValueDirty()?"Dirty":"Clean") << "] " ; break ;
1928 case AClean: os << " [ACLEAN] " ; break ;
1929 case ADirty: os << " [ADIRTY] " ; break ;
1930 }
1931 }
1932 os << endl ;
1933
1934 for (Int_t i=0 ;i<numCaches() ; i++) {
1936 }
1938 }
1939
1940 TString indent2(indent) ;
1941 indent2 += " " ;
1942 for (const auto arg : _serverList) {
1943 arg->printCompactTree(os,indent2,namePat,this) ;
1944 }
1945}
1946
1947
1948////////////////////////////////////////////////////////////////////////////////
1949/// Print tree structure of expression tree on given ostream, only branch nodes are printed.
1950/// Lead nodes (variables) will not be shown
1951///
1952/// If namePat is not "*", only nodes with names matching the pattern will be printed.
1953
1954void RooAbsArg::printComponentTree(const char* indent, const char* namePat, Int_t nLevel)
1955{
1956 if (nLevel==0) return ;
1957 if (isFundamental()) return ;
1958 RooResolutionModel* rmodel = dynamic_cast<RooResolutionModel*>(this) ;
1959 if (rmodel && rmodel->isConvolved()) return ;
1960 if (InheritsFrom("RooConstVar")) return ;
1961
1962 if ( !namePat || TString(GetName()).Contains(namePat)) {
1963 cout << indent ;
1964 Print() ;
1965 }
1966
1967 TString indent2(indent) ;
1968 indent2 += " " ;
1969 for (const auto arg : _serverList) {
1970 arg->printComponentTree(indent2.Data(),namePat,nLevel-1) ;
1971 }
1972}
1973
1974
1975////////////////////////////////////////////////////////////////////////////////
1976/// Construct a mangled name from the actual name that
1977/// is free of any math symbols that might be interpreted by TTree
1978
1980{
1981 // Check for optional alternate name of branch for this argument
1982 TString rawBranchName = GetName() ;
1983 if (getStringAttribute("BranchName")) {
1984 rawBranchName = getStringAttribute("BranchName") ;
1985 }
1986
1987 TString cleanName(rawBranchName) ;
1988 cleanName.ReplaceAll("/","D") ;
1989 cleanName.ReplaceAll("-","M") ;
1990 cleanName.ReplaceAll("+","P") ;
1991 cleanName.ReplaceAll("*","X") ;
1992 cleanName.ReplaceAll("[","L") ;
1993 cleanName.ReplaceAll("]","R") ;
1994 cleanName.ReplaceAll("(","L") ;
1995 cleanName.ReplaceAll(")","R") ;
1996 cleanName.ReplaceAll("{","L") ;
1997 cleanName.ReplaceAll("}","R") ;
1998
1999 return cleanName;
2000}
2001
2002
2003////////////////////////////////////////////////////////////////////////////////
2004/// Hook function interface for object to insert additional information
2005/// when printed in the context of a tree structure. This default
2006/// implementation prints nothing
2007
2008void RooAbsArg::printCompactTreeHook(ostream&, const char *)
2009{
2010}
2011
2012
2013////////////////////////////////////////////////////////////////////////////////
2014/// Register RooAbsCache with this object. This function is called
2015/// by RooAbsCache constructors for objects that are a datamember
2016/// of this RooAbsArg. By registering itself the RooAbsArg is aware
2017/// of all its cache data members and will forward server change
2018/// and cache mode change calls to the cache objects, which in turn
2019/// can forward them their contents
2020
2022{
2023 _cacheList.push_back(&cache) ;
2024}
2025
2026
2027////////////////////////////////////////////////////////////////////////////////
2028/// Unregister a RooAbsCache. Called from the RooAbsCache destructor
2029
2031{
2032 _cacheList.erase(std::remove(_cacheList.begin(), _cacheList.end(), &cache),
2033 _cacheList.end());
2034}
2035
2036
2037////////////////////////////////////////////////////////////////////////////////
2038/// Return number of registered caches
2039
2041{
2042 return _cacheList.size() ;
2043}
2044
2045
2046////////////////////////////////////////////////////////////////////////////////
2047/// Return registered cache object by index
2048
2050{
2051 return _cacheList[index] ;
2052}
2053
2054
2055////////////////////////////////////////////////////////////////////////////////
2056/// Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
2057
2059{
2060 return getParameters(RooArgSet(),stripDisconnected) ;
2061}
2062
2063
2064////////////////////////////////////////////////////////////////////////////////
2065/// Return ancestors in cloning chain of this RooAbsArg. NOTE: Returned pointers
2066/// are not guaranteed to be 'live', so do not dereference without proper caution
2067
2069{
2070 RooLinkedList retVal ;
2071
2072 set<string>::const_iterator iter= _boolAttrib.begin() ;
2073 while(iter != _boolAttrib.end()) {
2074 if (TString(*iter).BeginsWith("CloneOf(")) {
2075 char buf[128] ;
2076 strlcpy(buf,iter->c_str(),128) ;
2077 strtok(buf,"(") ;
2078 char* ptrToken = strtok(0,")") ;
2079 RooAbsArg* ptr = (RooAbsArg*) strtoll(ptrToken,0,16) ;
2080 retVal.Add(ptr) ;
2081 }
2082 }
2083
2084 return retVal ;
2085}
2086
2087
2088////////////////////////////////////////////////////////////////////////////////
2089/// Create a GraphViz .dot file visualizing the expression tree headed by
2090/// this RooAbsArg object. Use the GraphViz tool suite to make e.g. a gif
2091/// or ps file from the .dot file.
2092/// If a node derives from RooAbsReal, its current (unnormalised) value is
2093/// printed as well.
2094///
2095/// Based on concept developed by Kyle Cranmer.
2096
2097void RooAbsArg::graphVizTree(const char* fileName, const char* delimiter, bool useTitle, bool useLatex)
2098{
2099 ofstream ofs(fileName) ;
2100 if (!ofs) {
2101 coutE(InputArguments) << "RooAbsArg::graphVizTree() ERROR: Cannot open graphViz output file with name " << fileName << endl ;
2102 return ;
2103 }
2104 graphVizTree(ofs, delimiter, useTitle, useLatex) ;
2105}
2106
2107////////////////////////////////////////////////////////////////////////////////
2108/// Write the GraphViz representation of the expression tree headed by
2109/// this RooAbsArg object to the given ostream.
2110/// If a node derives from RooAbsReal, its current (unnormalised) value is
2111/// printed as well.
2112///
2113/// Based on concept developed by Kyle Cranmer.
2114
2115void RooAbsArg::graphVizTree(ostream& os, const char* delimiter, bool useTitle, bool useLatex)
2116{
2117 if (!os) {
2118 coutE(InputArguments) << "RooAbsArg::graphVizTree() ERROR: output stream provided as input argument is in invalid state" << endl ;
2119 }
2120
2121 // silent warning messages coming when evaluating a RooAddPdf without a normalization set
2123
2124 // Write header
2125 os << "digraph \"" << GetName() << "\"{" << endl ;
2126
2127 // First list all the tree nodes
2128 RooArgSet nodeSet ;
2129 treeNodeServerList(&nodeSet) ;
2130 RooFIter iter = nodeSet.fwdIterator() ;
2131 RooAbsArg* node ;
2132
2133 // iterate over nodes
2134 while((node=iter.next())) {
2135 string nodeName = node->GetName();
2136 string nodeTitle = node->GetTitle();
2137 string nodeLabel = (useTitle && !nodeTitle.empty()) ? nodeTitle : nodeName;
2138
2139 // if using latex, replace ROOT's # with normal latex backslash
2140 string::size_type position = nodeLabel.find("#") ;
2141 while(useLatex && position!=nodeLabel.npos){
2142 nodeLabel.replace(position, 1, "\\");
2143 }
2144
2145 string typeFormat = "\\texttt{";
2146 string nodeType = (useLatex) ? typeFormat+node->IsA()->GetName()+"}" : node->IsA()->GetName();
2147
2148 if (auto realNode = dynamic_cast<RooAbsReal*>(node)) {
2149 nodeLabel += delimiter + std::to_string(realNode->getVal());
2150 }
2151
2152 os << "\"" << nodeName << "\" [ color=" << (node->isFundamental()?"blue":"red")
2153 << ", label=\"" << nodeType << delimiter << nodeLabel << "\"];" << endl ;
2154
2155 }
2156
2157 // Get set of all server links
2158 set<pair<RooAbsArg*,RooAbsArg*> > links ;
2159 graphVizAddConnections(links) ;
2160
2161 // And write them out
2162 for(auto const& link : links) {
2163 os << "\"" << link.first->GetName() << "\" -> \"" << link.second->GetName() << "\";" << endl ;
2164 }
2165
2166 // Write trailer
2167 os << "}" << endl ;
2168
2169}
2170
2171////////////////////////////////////////////////////////////////////////////////
2172/// Utility function that inserts all point-to-point client-server connections
2173/// between any two RooAbsArgs in the expression tree headed by this object
2174/// in the linkSet argument.
2175
2176void RooAbsArg::graphVizAddConnections(set<pair<RooAbsArg*,RooAbsArg*> >& linkSet)
2177{
2178 for (const auto server : _serverList) {
2179 linkSet.insert(make_pair(this,server)) ;
2180 server->graphVizAddConnections(linkSet) ;
2181 }
2182}
2183
2184
2185////////////////////////////////////////////////////////////////////////////////
2186/// Take ownership of the contents of 'comps'.
2187
2189{
2190 if (!_ownedComponents) {
2191 _ownedComponents = new RooArgSet("owned components") ;
2192 }
2193 return _ownedComponents->addOwned(comps) ;
2194}
2195
2196
2197////////////////////////////////////////////////////////////////////////////////
2198/// Take ownership of the contents of 'comps'. Different from the overload that
2199/// taked the RooArgSet by `const&`, this version can also take an owning
2200/// RooArgSet without error, because the ownership will not be ambiguous afterwards.
2201
2203{
2204 if (!_ownedComponents) {
2205 _ownedComponents = new RooArgSet("owned components") ;
2206 }
2207 return _ownedComponents->addOwned(std::move(comps)) ;
2208}
2209
2210
2211////////////////////////////////////////////////////////////////////////////////
2212/// \copydoc RooAbsArg::addOwnedComponents(RooAbsCollection&& comps)
2213
2215 return addOwnedComponents(static_cast<RooAbsCollection&&>(std::move(comps)));
2216}
2217
2218
2219////////////////////////////////////////////////////////////////////////////////
2220/// Clone tree expression of objects. All tree nodes will be owned by
2221/// the head node return by cloneTree()
2222
2223RooAbsArg* RooAbsArg::cloneTree(const char* newname) const
2224{
2225 // Clone tree using snapshot
2226 RooArgSet clonedNodes;
2227 RooArgSet(*this).snapshot(clonedNodes, true);
2228
2229 // Find the head node in the cloneSet
2230 RooAbsArg* head = clonedNodes.find(*this) ;
2231 assert(head);
2232
2233 // Remove the head node from the cloneSet
2234 // To release it from the set ownership
2235 clonedNodes.remove(*head) ;
2236
2237 // Add the set as owned component of the head
2238 head->addOwnedComponents(std::move(clonedNodes)) ;
2239
2240 // Adjust name of head node if requested
2241 if (newname) {
2242 head->TNamed::SetName(newname) ;
2243 head->_namePtr = RooNameReg::instance().constPtr(newname) ;
2244 }
2245
2246 // Return the head
2247 return head ;
2248}
2249
2250
2251
2252////////////////////////////////////////////////////////////////////////////////
2253
2255{
2256 if (dynamic_cast<RooTreeDataStore*>(&store)) {
2257 attachToTree(((RooTreeDataStore&)store).tree()) ;
2258 } else if (dynamic_cast<RooVectorDataStore*>(&store)) {
2260 }
2261}
2262
2263
2264
2265////////////////////////////////////////////////////////////////////////////////
2266
2268{
2269 if (_eocache) {
2270 return *_eocache ;
2271 } else {
2273 }
2274}
2275
2276
2277////////////////////////////////////////////////////////////////////////////////
2278
2280{
2281 string suffix ;
2282
2283 RooArgSet branches ;
2284 branchNodeServerList(&branches) ;
2285 RooFIter iter = branches.fwdIterator();
2286 RooAbsArg* arg ;
2287 while((arg=iter.next())) {
2288 const char* tmp = arg->cacheUniqueSuffix() ;
2289 if (tmp) suffix += tmp ;
2290 }
2291 return Form("%s",suffix.c_str()) ;
2292}
2293
2294
2295////////////////////////////////////////////////////////////////////////////////
2296
2298{
2299 RooArgSet branches ;
2300 branchNodeServerList(&branches) ;
2301 for(auto const& arg : branches) {
2302 for (auto const& arg2 : arg->_cacheList) {
2303 arg2->wireCache() ;
2304 }
2305 }
2306}
2307
2308
2309
2310////////////////////////////////////////////////////////////////////////////////
2311
2312void RooAbsArg::SetName(const char* name)
2313{
2315 auto newPtr = RooNameReg::instance().constPtr(GetName()) ;
2316 if (newPtr != _namePtr) {
2317 //cout << "Rename '" << _namePtr->GetName() << "' to '" << name << "' (set flag in new name)" << endl;
2318 _namePtr = newPtr;
2321 }
2322}
2323
2324
2325
2326
2327////////////////////////////////////////////////////////////////////////////////
2328
2329void RooAbsArg::SetNameTitle(const char *name, const char *title)
2330{
2331 TNamed::SetTitle(title) ;
2332 SetName(name);
2333}
2334
2335
2336////////////////////////////////////////////////////////////////////////////////
2337/// Stream an object of class RooAbsArg.
2338
2339void RooAbsArg::Streamer(TBuffer &R__b)
2340{
2341 if (R__b.IsReading()) {
2342 _ioReadStack.push(this) ;
2343 R__b.ReadClassBuffer(RooAbsArg::Class(),this);
2344 _ioReadStack.pop() ;
2346 _isConstant = getAttribute("Constant") ;
2347 } else {
2348 R__b.WriteClassBuffer(RooAbsArg::Class(),this);
2349 }
2350}
2351
2352////////////////////////////////////////////////////////////////////////////////
2353/// Method called by workspace container to finalize schema evolution issues
2354/// that cannot be handled in a single ioStreamer pass.
2355///
2356/// A second pass is typically needed when evolving data member of RooAbsArg-derived
2357/// classes that are container classes with references to other members, which may
2358/// not yet be 'live' in the first ioStreamer() evolution pass.
2359///
2360/// Classes may overload this function, but must call the base method in the
2361/// overloaded call to ensure base evolution is handled properly
2362
2364{
2365 // Handling of v5-v6 migration (TRefArray _proxyList --> RooRefArray _proxyList)
2366 auto iter = _ioEvoList.find(this);
2367 if (iter != _ioEvoList.end()) {
2368
2369 // Transfer contents of saved TRefArray to RooRefArray now
2371 _proxyList.Expand(iter->second->GetEntriesFast());
2372 for (int i = 0; i < iter->second->GetEntriesFast(); i++) {
2373 _proxyList.Add(iter->second->At(i));
2374 }
2375 // Delete TRefArray and remove from list
2376 _ioEvoList.erase(iter);
2377 }
2378}
2379
2380
2381
2382
2383////////////////////////////////////////////////////////////////////////////////
2384/// Method called by workspace container to finalize schema evolution issues
2385/// that cannot be handled in a single ioStreamer pass. This static finalize method
2386/// is called after ioStreamerPass2() is called on each directly listed object
2387/// in the workspace. It's purpose is to complete schema evolution of any
2388/// objects in the workspace that are not directly listed as content elements
2389/// (e.g. analytical convolution tokens )
2390
2392{
2393 // Handling of v5-v6 migration (TRefArray _proxyList --> RooRefArray _proxyList)
2394 for (const auto& iter : _ioEvoList) {
2395
2396 // Transfer contents of saved TRefArray to RooRefArray now
2397 if (!iter.first->_proxyList.GetEntriesFast())
2398 iter.first->_proxyList.Expand(iter.second->GetEntriesFast());
2399 for (int i = 0; i < iter.second->GetEntriesFast(); i++) {
2400 iter.first->_proxyList.Add(iter.second->At(i));
2401 }
2402 }
2403
2404 _ioEvoList.clear();
2405}
2406
2407
2411}
2412
2413////////////////////////////////////////////////////////////////////////////////
2414/// Stream an object of class RooRefArray.
2415
2416void RooRefArray::Streamer(TBuffer &R__b)
2417{
2418 UInt_t R__s, R__c;
2419 if (R__b.IsReading()) {
2420
2421 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
2422
2423 // Make temporary refArray and read that from the streamer
2424 auto refArray = std::make_unique<TRefArray>();
2425 refArray->Streamer(R__b) ;
2426 R__b.CheckByteCount(R__s, R__c, refArray->IsA());
2427
2428 // Schedule deferred processing of TRefArray into proxy list
2429 RooAbsArg::_ioEvoList[RooAbsArg::_ioReadStack.top()] = std::move(refArray);
2430
2431 } else {
2432
2433 R__c = R__b.WriteVersion(RooRefArray::IsA(), kTRUE);
2434
2435 // Make a temporary refArray and write that to the streamer
2436 TRefArray refArray(GetEntriesFast());
2437 TIterator* iter = MakeIterator() ;
2438 TObject* tmpObj ; while ((tmpObj = iter->Next())) {
2439 refArray.Add(tmpObj) ;
2440 }
2441 delete iter ;
2442
2443 refArray.Streamer(R__b) ;
2444 R__b.SetByteCount(R__c, kTRUE) ;
2445
2446 }
2447}
2448
2449/// Print at the prompt
2450namespace cling {
2451std::string printValue(RooAbsArg *raa)
2452{
2453 std::stringstream s;
2454 if (0 == *raa->GetName() && 0 == *raa->GetTitle()) {
2455 s << "An instance of " << raa->ClassName() << ".";
2456 return s.str();
2457 }
2458 raa->printStream(s, raa->defaultPrintContents(""), raa->defaultPrintStyle(""));
2459 return s.str();
2460}
2461} // namespace cling
2462
2463
2464/// Disables or enables the usage of squared weights. Needs to be overloaded in
2465/// the likelihood classes for which this is relevant.
2467 for(auto * server : servers()) {
2468 server->applyWeightSquared(flag);
2469 }
2470}
2471
2472
2473/// Fills a RooArgSet to be used as the normalization set for a server, given a
2474/// normalization set for this RooAbsArg. If the output is a `nullptr`, it
2475/// means that the normalization set doesn't change.
2476///
2477/// \param[in] normSet The normalization set for this RooAbsArg.
2478/// \param[in] server A server of this RooAbsArg that we determine the
2479/// normalization set for.
2480/// \param[out] serverNormSet Output parameter. Normalization set for the
2481/// server.
2482std::unique_ptr<RooArgSet> RooAbsArg::fillNormSetForServer(RooArgSet const& /*normSet*/, RooAbsArg const& /*server*/) const {
2483 return nullptr;
2484}
istream & operator>>(istream &is, RooAbsArg &arg)
Istream operator.
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define coutF(a)
#define coutE(a)
#define cxcoutF(a)
char Text_t
Definition RtypesCore.h:62
short Version_t
Definition RtypesCore.h:65
const Bool_t kFALSE
Definition RtypesCore.h:101
bool Bool_t
Definition RtypesCore.h:63
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
TBuffer & operator<<(TBuffer &buf, const Tmpl *obj)
Definition TBuffer.h:399
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooExpensiveObjectCache & expensiveObjectCache() const
virtual Bool_t checkObservables(const RooArgSet *nset) const
Overloadable function in which derived classes can implement consistency checks of the variables.
virtual RooAbsArg * cloneTree(const char *newname=0) const
Clone tree expression of objects.
RooRefArray _proxyList
Definition RooAbsArg.h:655
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
Bool_t _fast
Definition RooAbsArg.h:735
virtual Bool_t readFromStream(std::istream &is, Bool_t compact, Bool_t verbose=kFALSE)=0
virtual Bool_t isValid() const
WVE (08/21/01) Probably obsolete now.
static void verboseDirty(Bool_t flag)
Activate verbose messaging related to dirty flag propagation.
void attachToStore(RooAbsDataStore &store)
Attach this argument to the data store such that it reads data from there.
virtual Int_t defaultPrintContents(Option_t *opt) const
Define default contents to print.
const TNamed * namePtr() const
De-duplicated pointer to this object's name.
Definition RooAbsArg.h:577
const char * aggregateCacheUniqueSuffix() const
Bool_t redirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t isRecursionStep=kFALSE)
Replace all direct servers of this object with the new servers in newServerList.
void changeServer(RooAbsArg &server, Bool_t valueProp, Bool_t shapeProp)
Change dirty flag propagation mask for specified server.
Bool_t isValueServer(const RooAbsArg &arg) const
Check if this is serving values to arg.
Definition RooAbsArg.h:213
ProxyListCache _proxyListCache
Definition RooAbsArg.h:710
RooWorkspace * _myws
Prevent 'AlwaysDirty' mode for this node.
Definition RooAbsArg.h:752
void printCompactTree(const char *indent="", const char *fileName=0, const char *namePat=0, RooAbsArg *client=0)
Print tree structure of expression tree on stdout, or to file if filename is specified.
void attachDataStore(const RooAbsDataStore &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
virtual void printClassName(std::ostream &os) const
Print object class name.
RooArgSet * _ownedComponents
Definition RooAbsArg.h:738
void setShapeDirty()
Notify that a shape-like property (e.g. binning) has changed.
Definition RooAbsArg.h:510
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
void registerProxy(RooArgProxy &proxy)
Register an RooArgProxy in the proxy list.
void attachArgs(const RooAbsCollection &set)
Bind this node to objects in set.
Bool_t isShapeServer(const RooAbsArg &arg) const
Check if this is serving shape to arg.
Definition RooAbsArg.h:221
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
static void ioStreamerPass2Finalize()
Method called by workspace container to finalize schema evolution issues that cannot be handled in a ...
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Implement multi-line detailed printing.
void addParameters(RooAbsCollection &params, const RooArgSet *nset=nullptr, bool stripDisconnected=true) const
Add all parameters of the function and its daughters to params.
Bool_t isCloneOf(const RooAbsArg &other) const
Check if this object was created as a clone of 'other'.
Bool_t isShapeDirty() const
Definition RooAbsArg.h:431
static void setDirtyInhibit(Bool_t flag)
Control global dirty inhibit mode.
Bool_t findConstantNodes(const RooArgSet &observables, RooArgSet &cacheList)
Find branch nodes with all-constant parameters, and add them to the list of nodes that can be cached ...
void graphVizAddConnections(std::set< std::pair< RooAbsArg *, RooAbsArg * > > &)
Utility function that inserts all point-to-point client-server connections between any two RooAbsArgs...
void unRegisterProxy(RooArgProxy &proxy)
Remove proxy from proxy list.
void printComponentTree(const char *indent="", const char *namePat=0, Int_t nLevel=999)
Print tree structure of expression tree on given ostream, only branch nodes are printed.
friend class RooArgSet
Definition RooAbsArg.h:642
Bool_t inhibitDirty() const
Delete watch flag.
std::set< std::string > _boolAttrib
Definition RooAbsArg.h:677
void addServer(RooAbsArg &server, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
virtual void printTree(std::ostream &os, TString indent="") const
Print object tree structure.
void unRegisterCache(RooAbsCache &cache)
Unregister a RooAbsCache. Called from the RooAbsCache destructor.
virtual Bool_t isFundamental() const
Is this object a fundamental type that can be added to a dataset? Fundamental-type subclasses overrid...
Definition RooAbsArg.h:239
virtual void printTitle(std::ostream &os) const
Print object title.
RefCountList_t _clientListValue
Definition RooAbsArg.h:653
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
void printAttribList(std::ostream &os) const
Transient boolean attributes (not copied in ctor)
virtual void Print(Option_t *options=0) const
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:337
friend void RooRefArray::Streamer(TBuffer &)
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
void treeNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t doBranch=kTRUE, Bool_t doLeaf=kTRUE, Bool_t valueOnly=kFALSE, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with nodes of the arg tree, following all server links, starting with ourself as t...
void setAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
virtual void writeToStream(std::ostream &os, Bool_t compact) const =0
virtual void printAddress(std::ostream &os) const
Print class name of object.
void setTransientAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
virtual void printCompactTreeHook(std::ostream &os, const char *ind="")
Hook function interface for object to insert additional information when printed in the context of a ...
const TNamed * _namePtr
Definition RooAbsArg.h:744
Bool_t _isConstant
De-duplicated name pointer. This will be equal for all objects with the same name.
Definition RooAbsArg.h:745
Bool_t isValueDirty() const
Definition RooAbsArg.h:436
Bool_t overlaps(const RooAbsArg &testArg, Bool_t valueOnly=kFALSE) const
Test if any of the nodes of tree are shared with that of the given tree.
Int_t Compare(const TObject *other) const
Utility function used by TCollection::Sort to compare contained TObjects We implement comparison by n...
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
virtual void ioStreamerPass2()
Method called by workspace container to finalize schema evolution issues that cannot be handled in a ...
void wireAllCaches()
Bool_t _valueDirty
Definition RooAbsArg.h:731
virtual ~RooAbsArg()
Destructor.
virtual Bool_t isDerived() const
Does value or shape of this arg depend on any other arg?
Definition RooAbsArg.h:89
static Bool_t _inhibitDirty
Definition RooAbsArg.h:714
virtual const char * cacheUniqueSuffix() const
Definition RooAbsArg.h:513
RefCountListLegacyIterator_t * makeLegacyIterator(const RefCountList_t &list) const
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:196
Bool_t _localNoInhibitDirty
Cached isConstant status.
Definition RooAbsArg.h:747
virtual void constOptimizeTestStatistic(ConstOpCode opcode, Bool_t doAlsoTrackingOpt=kTRUE)
Interface function signaling a request to perform constant term optimization.
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
virtual void attachToTree(TTree &t, Int_t bufSize=32000)=0
Overloadable function for derived classes to implement attachment as branch to a TTree.
OperMode _operMode
Definition RooAbsArg.h:734
RooLinkedList getCloningAncestors() const
Return ancestors in cloning chain of this RooAbsArg.
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:505
Bool_t recursiveCheckObservables(const RooArgSet *nset) const
Recursively call checkObservables on all nodes in the expression tree.
RooAbsCache * getCache(Int_t index) const
Return registered cache object by index.
void printDirty(Bool_t depth=kTRUE) const
Print information about current value dirty state information.
static Bool_t _verboseDirty
cache of the list of proxies. Avoids type casting.
Definition RooAbsArg.h:713
void registerCache(RooAbsCache &cache)
Register RooAbsCache with this object.
virtual void optimizeCacheMode(const RooArgSet &observables)
Activate cache mode optimization with given definition of observables.
RefCountList_t _clientListShape
Definition RooAbsArg.h:652
virtual void attachToVStore(RooVectorDataStore &vstore)=0
TString cleanBranchName() const
Construct a mangled name from the actual name that is free of any math symbols that might be interpre...
void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Bool_t _prohibitServerRedirect
Set of owned component.
Definition RooAbsArg.h:740
void removeServer(RooAbsArg &server, Bool_t force=kFALSE)
Unregister another RooAbsArg as a server to us, ie, declare that we no longer depend on its value and...
Int_t numProxies() const
Return the number of registered proxies.
Bool_t getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
virtual void printName(std::ostream &os) const
Print object name.
void replaceServer(RooAbsArg &oldServer, RooAbsArg &newServer, Bool_t valueProp, Bool_t shapeProp)
Replace 'oldServer' with 'newServer'.
virtual void getParametersHook(const RooArgSet *, RooArgSet *, Bool_t) const
Definition RooAbsArg.h:608
virtual void printMetaArgs(std::ostream &) const
Definition RooAbsArg.h:347
void addServerList(RooAbsCollection &serverList, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE)
Register a list of RooAbsArg as servers to us by calling addServer() for each arg in the list.
virtual void applyWeightSquared(bool flag)
Disables or enables the usage of squared weights.
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Check whether this object depends on values from an element in the serverList.
Definition RooAbsArg.h:99
void setProxyNormSet(const RooArgSet *nset)
Forward a change in the cached normalization argset to all the registered proxies.
RefCountList_t _clientList
Definition RooAbsArg.h:651
RooAbsProxy * getProxy(Int_t index) const
Return the nth proxy from the proxy list.
RooAbsArg & operator=(const RooAbsArg &other)
Assign all boolean and string properties of the original object.
virtual bool redirectServersHook(const RooAbsCollection &, bool, bool, bool)
Function that is called at the end of redirectServers().
Definition RooAbsArg.h:271
RefCountList_t _serverList
Definition RooAbsArg.h:650
RooExpensiveObjectCache * _eocache
Prohibit server redirects – Debugging tool.
Definition RooAbsArg.h:742
RooArgSet * getComponents() const
Create a RooArgSet with all components (branch nodes) of the expression tree headed by this object.
std::set< std::string > _boolAttribTransient
Definition RooAbsArg.h:679
Bool_t _shapeDirty
Definition RooAbsArg.h:732
RooArgSet * getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooAbsArg * findNewServer(const RooAbsCollection &newSet, Bool_t nameChange) const
Find the new server in the specified set that matches the old server.
Bool_t isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:377
void branchNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all branch nodes of the arg tree starting with ourself as top node.
std::size_t getParametersSizeEstimate(const RooArgSet *nset=nullptr) const
Obtain an estimate of the number of parameters of the function and its daughters.
void graphVizTree(const char *fileName, const char *delimiter="\n", bool useTitle=false, bool useLatex=false)
Create a GraphViz .dot file visualizing the expression tree headed by this RooAbsArg object.
virtual void operModeHook()
Definition RooAbsArg.h:602
Bool_t getTransientAttribute(const Text_t *name) const
Check if a named attribute is set.
std::map< std::string, std::string > _stringAttrib
Definition RooAbsArg.h:678
Int_t numCaches() const
Return number of registered caches.
RooAbsArg()
Default constructor.
virtual std::unique_ptr< RooArgSet > fillNormSetForServer(RooArgSet const &normSet, RooAbsArg const &server) const
Fills a RooArgSet to be used as the normalization set for a server, given a normalization set for thi...
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Recursively replace all servers with the new servers in newSet.
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
virtual void printArgs(std::ostream &os) const
Print object arguments, ie its proxies.
TIteratorToSTLInterface< RefCountList_t::Container_t > RefCountListLegacyIterator_t
Definition RooAbsArg.h:72
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition RooAbsArg.h:200
void SetName(const char *name)
Set the name of the TNamed.
Bool_t _deleteWatch
Definition RooAbsArg.h:715
std::vector< RooAbsCache * > _cacheList
Definition RooAbsArg.h:657
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:499
Bool_t observableOverlaps(const RooAbsData *dset, const RooAbsArg &testArg) const
Test if any of the dependents of the arg tree (as determined by getObservables) overlaps with those o...
RooAbsCache is the abstract base class for data members of RooAbsArgs that cache other (composite) Ro...
Definition RooAbsCache.h:27
virtual void operModeHook()
Interface for operation mode changes.
Definition RooAbsCache.h:46
virtual void findConstantNodes(const RooArgSet &, RooArgSet &, RooLinkedList &)
Interface for constant term node finding calls.
Definition RooAbsCache.h:52
virtual Bool_t redirectServersHook(const RooAbsCollection &, bool, bool, bool)
Interface for server redirect calls.
Definition RooAbsCache.h:40
virtual void printCompactTreeHook(std::ostream &, const char *)
Interface for printing of cache guts in tree mode printing.
Definition RooAbsCache.h:55
virtual void optimizeCacheMode(const RooArgSet &, RooArgSet &, RooLinkedList &)
Interface for processing of cache mode optimization calls.
Definition RooAbsCache.h:49
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
void sort(Bool_t reverse=false)
Sort collection using std::sort and name comparison.
RooFIter fwdIterator() const
One-time forward iterator.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add an argument and transfer the ownership to the collection.
Storage_t::size_type size() const
RooAbsArg * first() const
void reserve(Storage_t::size_type count)
void clear()
Clear contents. If the collection is owning, it will also delete the contents.
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
RooAbsCollection * selectByAttrib(const char *name, Bool_t value) const
Create a subset of the current collection, consisting only of those elements with the specified attri...
const char * GetName() const
Returns name of object.
void setName(const char *name)
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsDataStore is the abstract base class for data collection that use a TTree as internal storage m...
virtual const RooArgSet * get(Int_t index) const =0
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:82
virtual const RooArgSet * get() const
Definition RooAbsData.h:128
RooAbsProxy is the abstact interface for proxy classes.
Definition RooAbsProxy.h:30
virtual void print(std::ostream &os, Bool_t addContents=kFALSE) const
Print proxy name.
virtual const char * name() const
Definition RooAbsProxy.h:40
virtual Bool_t changePointer(const RooAbsCollection &newServerSet, Bool_t nameChange=kFALSE, Bool_t factoryInitMode=kFALSE)=0
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgProxy is the abstract interface for RooAbsArg proxy classes.
Definition RooArgProxy.h:24
Bool_t isShapeServer() const
Definition RooArgProxy.h:55
Bool_t isValueServer() const
Definition RooArgProxy.h:51
RooAbsArg * absArg() const
Definition RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:158
RooConstVar represent a constant real-valued object.
Definition RooConstVar.h:26
RooExpensiveObjectCache is a singleton class that serves as repository for objects that are expensive...
static RooExpensiveObjectCache & instance()
Return reference to singleton instance.
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
Switches the message service to a different level while the instance is alive.
Definition RooHelpers.h:42
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
TObject * FindObject(const char *name) const
Return pointer to obejct with given name.
RooAbsArg * findArg(const RooAbsArg *) const
Return pointer to object with given name in collection.
virtual void Add(TObject *arg)
RooLinkedListIterImpl begin() const
RooListProxy is the concrete proxy for RooArgList objects.
RooNameReg is a registry for const char* names.
Definition RooNameReg.h:25
const TNamed * constPtr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
static RooNameReg & instance()
Return reference to singleton instance.
static void incrementRenameCounter()
The renaming counter has to be incremented every time a RooAbsArg is renamed.
RooPlotable is a 'mix-in' base class that define the standard RooFit plotting and printing methods.
virtual StyleOption defaultPrintStyle(Option_t *opt) const
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,...
virtual void printValue(std::ostream &os) const
Interface to print value of object.
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
Bool_t isConvolved() const
std::size_t refCount(typename Container_t::const_iterator item) const
Return ref count of item that iterator points to.
const Container_t & containedObjects() const
Direct reference to container of objects held by this list.
void Remove(const T *obj, bool force=false)
Decrease ref count of given object.
Container_t::const_iterator begin() const
Iterator over contained objects.
void Add(T *obj, std::size_t initialCount=1)
Add an object or increase refCount if it is already present.
Container_t::const_iterator end() const
End of contained objects.
bool empty() const
Check if empty.
std::size_t size() const
Number of contained objects (neglecting the ref count).
void RemoveAll(const T *obj)
Remove from list irrespective of ref count.
void reserve(std::size_t amount)
bool containsByNamePtr(const T *obj) const
Check if list contains an item using findByNamePointer().
RooSetProxy is the concrete proxy for RooArgSet objects.
Definition RooSetProxy.h:23
RooTreeDataStore is a TTree-backed data storage.
RooVectorDataStore uses std::vectors to store data columns.
Bool_t defineSetInternal(const char *name, const RooArgSet &aset)
const RooArgSet * set(const char *name)
Return pointer to previously defined named set with given nmame If no such set is found a null pointe...
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
TIterator and GenericRooFIter front end with STL back end.
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
TString fName
Definition TNamed.h:32
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
virtual const char * GetTitle() const
Returns title of object.
Definition TNamed.h:48
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Int_t GetEntriesFast() const
Definition TObjArray.h:58
virtual void Expand(Int_t newSize)
Expand or shrink the array to newSize elements.
void Add(TObject *obj)
Definition TObjArray.h:68
virtual void Compress()
Remove empty slots from array.
Int_t GetEntries() const
Return the number of objects in array (i.e.
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
TIterator * MakeIterator(Bool_t dir=kIterForward) const
Returns an array iterator.
virtual TObject * Remove(TObject *obj)
Remove object from array.
TObject * At(Int_t idx) const
Definition TObjArray.h:164
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:429
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:515
An array of references to TObjects.
Definition TRefArray.h:33
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:692
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Definition TString.h:615
TString & Append(const char *cs)
Definition TString.h:564
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
A TTree represents a columnar dataset.
Definition TTree.h:79
std::string getColonSeparatedNameString(RooArgSet const &argSet)
Create a string with all sorted names of RooArgSet elements separated by colons.
Definition first.py:1
Definition tree.py:1
std::vector< RooAbsProxy * > cache
Definition RooAbsArg.h:707