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