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