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