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