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