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