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