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