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