Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooWorkspace.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  * *
8  * Copyright (c) 2000-2005, Regents of the University of California *
9  * and Stanford University. All rights reserved. *
10  * *
11  * Redistribution and use in source and binary forms, *
12  * with or without modification, are permitted according to the terms *
13  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14  *****************************************************************************/
15 
16 /**
17 \file RooWorkspace.cxx
18 \class RooWorkspace
19 \ingroup Roofitcore
20 
21 The RooWorkspace is a persistable container for RooFit projects. A workspace
22 can contain and own variables, p.d.f.s, functions and datasets. All objects
23 that live in the workspace are owned by the workspace. The `import()` method
24 enforces consistency of objects upon insertion into the workspace (e.g. no
25 duplicate object with the same name are allowed) and makes sure all objects
26 in the workspace are connected to each other. Easy accessor methods like
27 `pdf()`, `var()` and `data()` allow to refer to the contents of the workspace by
28 object name. The entire RooWorkspace can be saved into a ROOT TFile and organises
29 the consistent streaming of its contents without duplication.
30 If a RooWorkspace contains custom classes, i.e. classes not in the
31 ROOT distribution, portability of workspaces can be enhanced by
32 storing the source code of those classes in the workspace as well.
33 This process is also organized by the workspace through the
34 `importClassCode()` method.
35 **/
36 
37 #include "RooWorkspace.h"
38 #include "RooHelpers.h"
39 #include "RooFit.h"
40 #include "RooWorkspace.h"
41 #include "RooWorkspaceHandle.h"
42 #include "RooFit.h"
43 #include "RooAbsPdf.h"
44 #include "RooRealVar.h"
45 #include "RooCategory.h"
46 #include "RooAbsData.h"
47 #include "RooCmdConfig.h"
48 #include "RooMsgService.h"
49 #include "RooConstVar.h"
50 #include "RooResolutionModel.h"
51 #include "RooPlot.h"
52 #include "RooRandom.h"
53 #include "TInterpreter.h"
54 #include "TClassTable.h"
55 #include "TBaseClass.h"
56 #include "TSystem.h"
57 #include "TRegexp.h"
58 #include "RooFactoryWSTool.h"
59 #include "RooAbsStudy.h"
60 #include "RooTObjWrap.h"
61 #include "RooAbsOptTestStatistic.h"
62 #include "TROOT.h"
63 #include "TFile.h"
64 #include "TH1.h"
65 #include <map>
66 #include <sstream>
67 #include <string>
68 
69 using namespace std ;
70 
71 #include "TClass.h"
72 #include "Riostream.h"
73 #include <string.h>
74 
75 ClassImp(RooWorkspace);
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 
79 ClassImp(RooWorkspace::CodeRepo);
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 
83 ClassImp(RooWorkspace::WSDir);
84 
85 list<string> RooWorkspace::_classDeclDirList ;
86 list<string> RooWorkspace::_classImplDirList ;
87 string RooWorkspace::_classFileExportDir = ".wscode.%s.%s" ;
88 Bool_t RooWorkspace::_autoClass = kFALSE ;
89 
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Add `dir` to search path for class declaration (header) files. This is needed
93 /// to find class headers custom classes are imported into the workspace.
94 void RooWorkspace::addClassDeclImportDir(const char* dir)
95 {
96  _classDeclDirList.push_back(dir) ;
97 }
98 
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Add `dir` to search path for class implementation (.cxx) files. This is needed
102 /// to find class headers custom classes are imported into the workspace.
103 void RooWorkspace::addClassImplImportDir(const char* dir)
104 {
105  _classImplDirList.push_back(dir) ;
106 }
107 
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Specify the name of the directory in which embedded source
111 /// code is unpacked and compiled. The specified string may contain
112 /// one '%s' token which will be substituted by the workspace name
113 
114 void RooWorkspace::setClassFileExportDir(const char* dir)
115 {
116  if (dir) {
117  _classFileExportDir = dir ;
118  } else {
119  _classFileExportDir = ".wscode.%s.%s" ;
120  }
121 }
122 
123 
124 ////////////////////////////////////////////////////////////////////////////////
125 /// If flag is true, source code of classes not the the ROOT distribution
126 /// is automatically imported if on object of such a class is imported
127 /// in the workspace
128 
129 void RooWorkspace::autoImportClassCode(Bool_t flag)
130 {
131  _autoClass = flag ;
132 }
133 
134 
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Default constructor
138 
139 RooWorkspace::RooWorkspace() : _classes(this), _dir(nullptr), _factory(nullptr), _doExport(kFALSE), _openTrans(kFALSE)
140 {
141 }
142 
143 
144 
145 ////////////////////////////////////////////////////////////////////////////////
146 /// Construct empty workspace with given name and title
147 
148 RooWorkspace::RooWorkspace(const char* name, const char* title) :
149  TNamed(name,title?title:name), _classes(this), _dir(nullptr), _factory(nullptr), _doExport(kFALSE), _openTrans(kFALSE)
150 {
151 }
152 
153 
154 RooWorkspace::RooWorkspace(const char* name, Bool_t doCINTExport) :
155  TNamed(name,name), _classes(this), _dir(nullptr), _factory(nullptr), _doExport(kFALSE), _openTrans(kFALSE)
156 {
157  // Construct empty workspace with given name and option to export reference to all workspace contents to a CINT namespace with the same name
158  if (doCINTExport) {
159  exportToCint(name) ;
160  }
161 }
162 
163 
164 ////////////////////////////////////////////////////////////////////////////////
165 /// Workspace copy constructor
166 
167 RooWorkspace::RooWorkspace(const RooWorkspace& other) :
168  TNamed(other), _uuid(other._uuid), _classes(other._classes,this), _dir(nullptr), _factory(nullptr), _doExport(kFALSE), _openTrans(kFALSE)
169 {
170  // Copy owned nodes
171  other._allOwnedNodes.snapshot(_allOwnedNodes,kTRUE) ;
172 
173  // Copy datasets
174  TIterator* iter = other._dataList.MakeIterator() ;
175  TObject* data2 ;
176  while((data2=iter->Next())) {
177  _dataList.Add(data2->Clone()) ;
178  }
179  delete iter ;
180 
181  // Copy snapshots
182  TIterator* iter2 = other._snapshots.MakeIterator() ;
183  RooArgSet* snap ;
184  while((snap=(RooArgSet*)iter2->Next())) {
185  RooArgSet* snapClone = (RooArgSet*) snap->snapshot() ;
186  snapClone->setName(snap->GetName()) ;
187  _snapshots.Add(snapClone) ;
188  }
189  delete iter2 ;
190 
191  // Copy named sets
192  for (map<string,RooArgSet>::const_iterator iter3 = other._namedSets.begin() ; iter3 != other._namedSets.end() ; ++iter3) {
193  // Make RooArgSet with equivalent content of this workspace
194  RooArgSet* tmp = (RooArgSet*) _allOwnedNodes.selectCommon(iter3->second) ;
195  _namedSets[iter3->first].add(*tmp) ;
196  delete tmp ;
197  }
198 
199  // Copy generic objects
200  TIterator* iter4 = other._genObjects.MakeIterator() ;
201  TObject* gobj ;
202  while((gobj=iter4->Next())) {
203  TObject *theClone = gobj->Clone();
204 
205  auto handle = dynamic_cast<RooWorkspaceHandle*>(theClone);
206  if (handle) {
207  handle->ReplaceWS(this);
208  }
209 
210  _genObjects.Add(theClone);
211  }
212  delete iter4 ;
213 
214 }
215 
216 
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 /// Workspace destructor
220 
221 RooWorkspace::~RooWorkspace()
222 {
223  // Delete references to variables that were declared in CINT
224  if (_doExport) {
225  unExport() ;
226  }
227 
228  // Delete contents
229  _dataList.Delete() ;
230  if (_dir) {
231  delete _dir ;
232  }
233  _snapshots.Delete() ;
234 
235  // WVE named sets too?
236 
237  _genObjects.Delete() ;
238 }
239 
240 
241 ////////////////////////////////////////////////////////////////////////////////
242 /// Import a RooAbsArg or RooAbsData set from a workspace in a file. Filespec should be constructed as "filename:wspacename:objectname"
243 /// The arguments will be passed to the relevant import() or import(RooAbsData&, ...) import calls
244 
245 Bool_t RooWorkspace::import(const char* fileSpec,
246  const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
247  const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
248  const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
249 {
250  // Parse file/workspace/objectname specification
251  std::vector<std::string> tokens = RooHelpers::tokenise(fileSpec, ":");
252 
253  // Check that parsing was successful
254  if (tokens.size() != 3) {
255  std::ostringstream stream;
256  for (const auto& token : tokens) {
257  stream << "\n\t" << token;
258  }
259  coutE(InputArguments) << "RooWorkspace(" << GetName() << ") ERROR in file specification, expecting 'filename:wsname:objname', but '" << fileSpec << "' given."
260  << "\nTokens read are:" << stream.str() << endl;
261  return kTRUE ;
262  }
263 
264  const std::string& filename = tokens[0];
265  const std::string& wsname = tokens[1];
266  const std::string& objname = tokens[2];
267 
268  // Check that file can be opened
269  TFile* f = TFile::Open(filename.c_str()) ;
270  if (f==0) {
271  coutE(InputArguments) << "RooWorkspace(" << GetName() << ") ERROR opening file " << filename << endl ;
272  return 0 ;
273  }
274 
275  // That that file contains workspace
276  RooWorkspace* w = dynamic_cast<RooWorkspace*>(f->Get(wsname.c_str())) ;
277  if (w==0) {
278  coutE(InputArguments) << "RooWorkspace(" << GetName() << ") ERROR: No object named " << wsname << " in file " << filename
279  << " or object is not a RooWorkspace" << endl ;
280  return 0 ;
281  }
282 
283  // Check that workspace contains object and forward to appropriate import method
284  RooAbsArg* warg = w->arg(objname.c_str()) ;
285  if (warg) {
286  Bool_t ret = import(*warg,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9) ;
287  delete f ;
288  return ret ;
289  }
290  RooAbsData* wdata = w->data(objname.c_str()) ;
291  if (wdata) {
292  Bool_t ret = import(*wdata,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9) ;
293  delete f ;
294  return ret ;
295  }
296 
297  coutE(InputArguments) << "RooWorkspace(" << GetName() << ") ERROR: No RooAbsArg or RooAbsData object named " << objname
298  << " in workspace " << wsname << " in file " << filename << endl ;
299  return kTRUE ;
300 }
301 
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// Import multiple RooAbsArg objects into workspace. For details on arguments see documentation
305 /// of import() method for single RooAbsArg
306 
307 Bool_t RooWorkspace::import(const RooArgSet& args,
308  const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
309  const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
310  const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
311 {
312  unique_ptr<TIterator> iter(args.createIterator()) ;
313  RooAbsArg* oneArg ;
314  Bool_t ret(kFALSE) ;
315  while((oneArg=(RooAbsArg*)iter->Next())) {
316  ret |= import(*oneArg,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9) ;
317  }
318  return ret ;
319 }
320 
321 
322 
323 ////////////////////////////////////////////////////////////////////////////////
324 /// Import a RooAbsArg object, e.g. function, p.d.f or variable into the workspace. This import function clones the input argument and will
325 /// own the clone. If a composite object is offered for import, e.g. a p.d.f with parameters and observables, the
326 /// complete tree of objects is imported. If any of the _variables_ of a composite object (parameters/observables) are already
327 /// in the workspace the imported p.d.f. is connected to the already existing variables. If any of the _function_ objects (p.d.f, formulas)
328 /// to be imported already exists in the workspace an error message is printed and the import of the entire tree of objects is cancelled.
329 /// Several optional arguments can be provided to modify the import procedure.
330 ///
331 /// <table>
332 /// <tr><th> Accepted arguments
333 /// <tr><td> `RenameConflictNodes(const char* suffix)` <td> Add suffix to branch node name if name conflicts with existing node in workspace
334 /// <tr><td> `RenameAllNodes(const char* suffix)` <td> Add suffix to all branch node names including top level node
335 /// <tr><td> `RenameAllVariables(const char* suffix)` <td> Add suffix to all variables names
336 /// <tr><td> `RenameAllVariablesExcept(const char* suffix, const char* exceptionList)` <td> Add suffix to all variables names, except ones listed
337 /// <tr><td> `RenameVariable(const char* inputName, const char* outputName)` <td> Rename variable as specified upon import.
338 /// <tr><td> `RecycleConflictNodes()` <td> If any of the function objects to be imported already exist in the name space, connect the
339 /// imported expression to the already existing nodes.
340 /// \attention Use with care! If function definitions do not match, this alters the definition of your function upon import
341 ///
342 /// <tr><td> `Silence()` <td> Do not issue any info message
343 /// </table>
344 ///
345 /// The RenameConflictNodes, RenameNodes and RecycleConflictNodes arguments are mutually exclusive. The RenameVariable argument can be repeated
346 /// as often as necessary to rename multiple variables. Alternatively, a single RenameVariable argument can be given with
347 /// two comma separated lists.
348 
349 Bool_t RooWorkspace::import(const RooAbsArg& inArg,
350  const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
351  const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
352  const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
353 {
354  RooLinkedList args ;
355  args.Add((TObject*)&arg1) ;
356  args.Add((TObject*)&arg2) ;
357  args.Add((TObject*)&arg3) ;
358  args.Add((TObject*)&arg4) ;
359  args.Add((TObject*)&arg5) ;
360  args.Add((TObject*)&arg6) ;
361  args.Add((TObject*)&arg7) ;
362  args.Add((TObject*)&arg8) ;
363  args.Add((TObject*)&arg9) ;
364 
365  // Select the pdf-specific commands
366  RooCmdConfig pc(Form("RooWorkspace::import(%s)",GetName())) ;
367 
368  pc.defineString("conflictSuffix","RenameConflictNodes",0) ;
369  pc.defineInt("renameConflictOrig","RenameConflictNodes",0,0) ;
370  pc.defineString("allSuffix","RenameAllNodes",0) ;
371  pc.defineString("allVarsSuffix","RenameAllVariables",0) ;
372  pc.defineString("allVarsExcept","RenameAllVariables",1) ;
373  pc.defineString("varChangeIn","RenameVar",0,"",kTRUE) ;
374  pc.defineString("varChangeOut","RenameVar",1,"",kTRUE) ;
375  pc.defineString("factoryTag","FactoryTag",0) ;
376  pc.defineInt("useExistingNodes","RecycleConflictNodes",0,0) ;
377  pc.defineInt("silence","Silence",0,0) ;
378  pc.defineInt("noRecursion","NoRecursion",0,0) ;
379  pc.defineMutex("RenameConflictNodes","RenameAllNodes") ;
380  pc.defineMutex("RenameConflictNodes","RecycleConflictNodes") ;
381  pc.defineMutex("RenameAllNodes","RecycleConflictNodes") ;
382  pc.defineMutex("RenameVariable","RenameAllVariables") ;
383 
384  // Process and check varargs
385  pc.process(args) ;
386  if (!pc.ok(kTRUE)) {
387  return kTRUE ;
388  }
389 
390  // Decode renaming logic into suffix string and boolean for conflictOnly mode
391  const char* suffixC = pc.getString("conflictSuffix") ;
392  const char* suffixA = pc.getString("allSuffix") ;
393  const char* suffixV = pc.getString("allVarsSuffix") ;
394  const char* exceptVars = pc.getString("allVarsExcept") ;
395  const char* varChangeIn = pc.getString("varChangeIn") ;
396  const char* varChangeOut = pc.getString("varChangeOut") ;
397  Bool_t renameConflictOrig = pc.getInt("renameConflictOrig") ;
398  Int_t useExistingNodes = pc.getInt("useExistingNodes") ;
399  Int_t silence = pc.getInt("silence") ;
400  Int_t noRecursion = pc.getInt("noRecursion") ;
401 
402 
403  // Turn zero length strings into null pointers
404  if (suffixC && strlen(suffixC)==0) suffixC = 0 ;
405  if (suffixA && strlen(suffixA)==0) suffixA = 0 ;
406 
407  Bool_t conflictOnly = suffixA ? kFALSE : kTRUE ;
408  const char* suffix = suffixA ? suffixA : suffixC ;
409 
410  // Process any change in variable names
411  map<string,string> varMap ;
412  if (strlen(varChangeIn)>0) {
413 
414  // Parse comma separated lists into map<string,string>
415  const std::vector<std::string> tokIn = RooHelpers::tokenise(varChangeIn, ", ");
416  const std::vector<std::string> tokOut = RooHelpers::tokenise(varChangeOut, ", ");
417  for (unsigned int i=0; i < tokIn.size(); ++i) {
418  varMap.insert(std::make_pair(tokIn[i], tokOut[i]));
419  }
420 
421  assert(tokIn.size() == tokOut.size());
422  }
423 
424  // Process RenameAllVariables argument if specified
425  // First convert exception list if provided
426  std::set<string> exceptVarNames ;
427  if (exceptVars && strlen(exceptVars)) {
428  const std::vector<std::string> toks = RooHelpers::tokenise(exceptVars, ", ");
429  exceptVarNames.insert(toks.begin(), toks.end());
430  }
431 
432  if (suffixV != 0 && strlen(suffixV)>0) {
433  RooArgSet* vars = inArg.getVariables() ;
434  TIterator* iter = vars->createIterator() ;
435  RooAbsArg* v ;
436  while((v=(RooAbsArg*)iter->Next())) {
437  if (exceptVarNames.find(v->GetName())==exceptVarNames.end()) {
438  varMap[v->GetName()] = Form("%s_%s",v->GetName(),suffixV) ;
439  }
440  }
441  delete iter ;
442  delete vars ;
443  }
444 
445  // Scan for overlaps with current contents
446  RooAbsArg* wsarg = _allOwnedNodes.find(inArg.GetName()) ;
447 
448  // Check for factory specification match
449  const char* tagIn = inArg.getStringAttribute("factory_tag") ;
450  const char* tagWs = wsarg ? wsarg->getStringAttribute("factory_tag") : 0 ;
451  Bool_t factoryMatch = (tagIn && tagWs && !strcmp(tagIn,tagWs)) ;
452  if (factoryMatch) {
453  ((RooAbsArg&)inArg).setAttribute("RooWorkspace::Recycle") ;
454  }
455 
456  if (!suffix && wsarg && !useExistingNodes && !(inArg.isFundamental() && varMap[inArg.GetName()]!="")) {
457  if (!factoryMatch) {
458  if (wsarg!=&inArg) {
459  coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR importing object named " << inArg.GetName()
460  << ": another instance with same name already in the workspace and no conflict resolution protocol specified" << endl ;
461  return kTRUE ;
462  } else {
463  if (!silence) {
464  coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") Object " << inArg.GetName() << " is already in workspace!" << endl ;
465  }
466  return kTRUE ;
467  }
468  } else {
469  coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") Recycling existing object " << inArg.GetName() << " created with identical factory specification" << endl ;
470  }
471  }
472 
473  // Make list of conflicting nodes
474  RooArgSet conflictNodes ;
475  RooArgSet branchSet ;
476  if (noRecursion) {
477  branchSet.add(inArg) ;
478  } else {
479  inArg.branchNodeServerList(&branchSet) ;
480  }
481  TIterator* iter = branchSet.createIterator() ;
482  RooAbsArg* branch ;
483  while ((branch=(RooAbsArg*)iter->Next())) {
484  RooAbsArg* wsbranch = _allOwnedNodes.find(branch->GetName()) ;
485  if (wsbranch && wsbranch!=branch && !branch->getAttribute("RooWorkspace::Recycle") && !useExistingNodes) {
486  conflictNodes.add(*branch) ;
487  }
488  }
489  delete iter ;
490 
491  // Terminate here if there are conflicts and no resolution protocol
492  if (conflictNodes.getSize()>0 && !suffix && !useExistingNodes) {
493  coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR object named " << inArg.GetName() << ": component(s) "
494  << conflictNodes << " already in the workspace and no conflict resolution protocol specified" << endl ;
495  return kTRUE ;
496  }
497 
498  // Now create a working copy of the incoming object tree
499  RooArgSet* cloneSet = (RooArgSet*) RooArgSet(inArg).snapshot(noRecursion==kFALSE) ;
500  RooAbsArg* cloneTop = cloneSet->find(inArg.GetName()) ;
501 
502  // Mark all nodes for renaming if we are not in conflictOnly mode
503  if (!conflictOnly) {
504  conflictNodes.removeAll() ;
505  conflictNodes.add(branchSet) ;
506  }
507 
508  // Mark nodes that are to be renamed with special attribute
509  string topName2 = cloneTop->GetName() ;
510  if (!renameConflictOrig) {
511  // Mark all nodes to be imported for renaming following conflict resolution protocol
512  TIterator* citer = conflictNodes.createIterator() ;
513  RooAbsArg* cnode ;
514  while ((cnode=(RooAbsArg*)citer->Next())) {
515  RooAbsArg* cnode2 = cloneSet->find(cnode->GetName()) ;
516  string origName = cnode2->GetName() ;
517  cnode2->SetName(Form("%s_%s",cnode2->GetName(),suffix)) ;
518  cnode2->SetTitle(Form("%s (%s)",cnode2->GetTitle(),suffix)) ;
519  string tag = Form("ORIGNAME:%s",origName.c_str()) ;
520  cnode2->setAttribute(tag.c_str()) ;
521  if (!cnode2->getStringAttribute("origName")) {
522  string tag2 = Form("%s",origName.c_str()) ;
523  cnode2->setStringAttribute("origName",tag2.c_str()) ;
524  }
525 
526  // Save name of new top level node for later use
527  if (cnode2==cloneTop) {
528  topName2 = cnode2->GetName() ;
529  }
530 
531  if (!silence) {
532  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName()
533  << ") Resolving name conflict in workspace by changing name of imported node "
534  << origName << " to " << cnode2->GetName() << endl ;
535  }
536  }
537  delete citer ;
538  } else {
539 
540  // Rename all nodes already in the workspace to 'clear the way' for the imported nodes
541  TIterator* citer = conflictNodes.createIterator() ;
542  RooAbsArg* cnode ;
543  while ((cnode=(RooAbsArg*)citer->Next())) {
544 
545  string origName = cnode->GetName() ;
546  RooAbsArg* wsnode = _allOwnedNodes.find(origName.c_str()) ;
547  if (wsnode) {
548 
549  if (!wsnode->getStringAttribute("origName")) {
550  wsnode->setStringAttribute("origName",wsnode->GetName()) ;
551  }
552 
553  if (!_allOwnedNodes.find(Form("%s_%s",cnode->GetName(),suffix))) {
554  wsnode->SetName(Form("%s_%s",cnode->GetName(),suffix)) ;
555  wsnode->SetTitle(Form("%s (%s)",cnode->GetTitle(),suffix)) ;
556  } else {
557  // Name with suffix already taken, add additional suffix
558  Int_t n=1 ;
559  while (true) {
560  string newname = Form("%s_%s_%d",cnode->GetName(),suffix,n) ;
561  if (!_allOwnedNodes.find(newname.c_str())) {
562  wsnode->SetName(newname.c_str()) ;
563  wsnode->SetTitle(Form("%s (%s %d)",cnode->GetTitle(),suffix,n)) ;
564  break ;
565  }
566  n++ ;
567  }
568  }
569  if (!silence) {
570  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName()
571  << ") Resolving name conflict in workspace by changing name of original node "
572  << origName << " to " << wsnode->GetName() << endl ;
573  }
574  } else {
575  coutW(ObjectHandling) << "RooWorkspce::import(" << GetName() << ") Internal error: expected to find existing node "
576  << origName << " to be renamed, but didn't find it..." << endl ;
577  }
578 
579  }
580  delete citer ;
581 
582  }
583 
584  // Process any change in variable names
585  if (strlen(varChangeIn)>0 || (suffixV && strlen(suffixV)>0)) {
586 
587  // Process all changes in variable names
588  TIterator* cliter = cloneSet->createIterator() ;
589  RooAbsArg* cnode ;
590  while ((cnode=(RooAbsArg*)cliter->Next())) {
591 
592  if (varMap.find(cnode->GetName())!=varMap.end()) {
593  string origName = cnode->GetName() ;
594  cnode->SetName(varMap[cnode->GetName()].c_str()) ;
595  string tag = Form("ORIGNAME:%s",origName.c_str()) ;
596  cnode->setAttribute(tag.c_str()) ;
597  if (!cnode->getStringAttribute("origName")) {
598  string tag2 = Form("%s",origName.c_str()) ;
599  cnode->setStringAttribute("origName",tag2.c_str()) ;
600  }
601 
602  if (!silence) {
603  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") Changing name of variable "
604  << origName << " to " << cnode->GetName() << " on request" << endl ;
605  }
606 
607  if (cnode==cloneTop) {
608  topName2 = cnode->GetName() ;
609  }
610 
611  }
612  }
613  delete cliter ;
614  }
615 
616  // Now clone again with renaming effective
617  RooArgSet* cloneSet2 = (RooArgSet*) RooArgSet(*cloneTop).snapshot(noRecursion==kFALSE) ;
618  RooAbsArg* cloneTop2 = cloneSet2->find(topName2.c_str()) ;
619 
620  // Make final check list of conflicting nodes
621  RooArgSet conflictNodes2 ;
622  RooArgSet branchSet2 ;
623  //inArg.branchNodeServerList(&branchSet) ; // WVE not needed
624  TIterator* iter2 = branchSet2.createIterator() ;
625  RooAbsArg* branch2 ;
626  while ((branch2=(RooAbsArg*)iter2->Next())) {
627  if (_allOwnedNodes.find(branch2->GetName())) {
628  conflictNodes2.add(*branch2) ;
629  }
630  }
631  delete iter2 ;
632 
633  // Terminate here if there are conflicts and no resolution protocol
634  if (conflictNodes2.getSize()) {
635  coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR object named " << inArg.GetName() << ": component(s) "
636  << conflictNodes2 << " cause naming conflict after conflict resolution protocol was executed" << endl ;
637  return kTRUE ;
638  }
639 
640  // Print a message for each imported node
641  iter = cloneSet2->createIterator() ;
642 
643  // Perform any auxiliary imports at this point
644  RooAbsArg* node ;
645  while((node=(RooAbsArg*)iter->Next())) {
646  if (node->importWorkspaceHook(*this)) {
647  coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR object named " << node->GetName()
648  << " has an error in importing in one or more of its auxiliary objects, aborting" << endl ;
649  return kTRUE ;
650  }
651  }
652  iter->Reset() ;
653 
654  if (cloneSet2->getSize()+_allOwnedNodes.getSize() > 999) _allOwnedNodes.setHashTableSize(1000);
655 
656  RooArgSet recycledNodes ;
657  RooArgSet nodesToBeDeleted ;
658  while((node=(RooAbsArg*)iter->Next())) {
659 
660  if (_autoClass) {
661  if (!_classes.autoImportClass(node->IsA())) {
662  coutW(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") WARNING: problems import class code of object "
663  << node->IsA()->GetName() << "::" << node->GetName() << ", reading of workspace will require external definition of class" << endl ;
664  }
665  }
666 
667  // Point expensiveObjectCache to copy in this workspace
668  RooExpensiveObjectCache& oldCache = node->expensiveObjectCache() ;
669  node->setExpensiveObjectCache(_eocache) ;
670  _eocache.importCacheObjects(oldCache,node->GetName(),kTRUE) ;
671 
672  // Check if node is already in workspace (can only happen for variables or identical instances, unless RecycleConflictNodes is specified)
673  RooAbsArg* wsnode = _allOwnedNodes.find(node->GetName()) ;
674 
675  if (wsnode) {
676  // Do not import node, add not to list of nodes that require reconnection
677  if (!silence && useExistingNodes) {
678  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") using existing copy of " << node->IsA()->GetName()
679  << "::" << node->GetName() << " for import of " << cloneTop2->IsA()->GetName() << "::"
680  << cloneTop2->GetName() << endl ;
681  }
682  recycledNodes.add(*_allOwnedNodes.find(node->GetName())) ;
683 
684  // Delete clone of incoming node
685  nodesToBeDeleted.addOwned(*node) ;
686 
687  //cout << "WV: recycling existing node " << existingNode << " = " << existingNode->GetName() << " for imported node " << node << endl ;
688 
689  } else {
690  // Import node
691  if (!silence) {
692  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") importing " << node->IsA()->GetName() << "::"
693  << node->GetName() << endl ;
694  }
695  _allOwnedNodes.addOwned(*node) ;
696  if (_openTrans) {
697  _sandboxNodes.add(*node) ;
698  } else {
699  if (_dir && node->IsA() != RooConstVar::Class()) {
700  _dir->InternalAppend(node) ;
701  }
702  if (_doExport && node->IsA() != RooConstVar::Class()) {
703  exportObj(node) ;
704  }
705  }
706  }
707  }
708 
709  // Release working copy
710  // no need to do a safe list since it was generated from a snapshot
711  // just take ownership and delte elements by hand
712  cloneSet->releaseOwnership() ;
713  RooFIter cloneSet_iter = cloneSet->fwdIterator() ;
714  RooAbsArg* cloneNode ;
715  while ((cloneNode=(RooAbsArg*)cloneSet_iter.next())) {
716  delete cloneNode;
717  }
718  delete cloneSet ;
719 
720  // Reconnect any nodes that need to be
721  if (recycledNodes.getSize()>0) {
722  iter->Reset() ;
723  while((node=(RooAbsArg*)iter->Next())) {
724  node->redirectServers(recycledNodes) ;
725  }
726  }
727  delete iter ;
728 
729  cloneSet2->releaseOwnership() ;
730  delete cloneSet2 ;
731 
732  return kFALSE ;
733 }
734 
735 
736 
737 ////////////////////////////////////////////////////////////////////////////////
738 /// Import a dataset (RooDataSet or RooDataHist) into the work space. The workspace will contain a copy of the data.
739 /// The dataset and its variables can be renamed upon insertion with the options below
740 ///
741 /// <table>
742 /// <tr><th> Accepted arguments
743 /// <tr><td> `Rename(const char* suffix)` <td> Rename dataset upon insertion
744 /// <tr><td> `RenameVariable(const char* inputName, const char* outputName)` <td> Change names of observables in dataset upon insertion
745 /// <tr><td> `Silence` <td> Be quiet, except in case of errors
746 
747 Bool_t RooWorkspace::import(RooAbsData& inData,
748  const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
749  const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
750  const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
751 
752 {
753 
754  RooLinkedList args ;
755  args.Add((TObject*)&arg1) ;
756  args.Add((TObject*)&arg2) ;
757  args.Add((TObject*)&arg3) ;
758  args.Add((TObject*)&arg4) ;
759  args.Add((TObject*)&arg5) ;
760  args.Add((TObject*)&arg6) ;
761  args.Add((TObject*)&arg7) ;
762  args.Add((TObject*)&arg8) ;
763  args.Add((TObject*)&arg9) ;
764 
765  // Select the pdf-specific commands
766  RooCmdConfig pc(Form("RooWorkspace::import(%s)",GetName())) ;
767 
768  pc.defineString("dsetName","Rename",0,"") ;
769  pc.defineString("varChangeIn","RenameVar",0,"",kTRUE) ;
770  pc.defineString("varChangeOut","RenameVar",1,"",kTRUE) ;
771  pc.defineInt("embedded","Embedded",0,0) ;
772  pc.defineInt("silence","Silence",0,0) ;
773 
774  // Process and check varargs
775  pc.process(args) ;
776  if (!pc.ok(kTRUE)) {
777  return kTRUE ;
778  }
779 
780  // Decode renaming logic into suffix string and boolean for conflictOnly mode
781  const char* dsetName = pc.getString("dsetName") ;
782  const char* varChangeIn = pc.getString("varChangeIn") ;
783  const char* varChangeOut = pc.getString("varChangeOut") ;
784  Bool_t embedded = pc.getInt("embedded") ;
785  Int_t silence = pc.getInt("silence") ;
786 
787  if (!silence)
788  coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") importing dataset " << inData.GetName() << endl ;
789 
790  // Transform emtpy string into null pointer
791  if (dsetName && strlen(dsetName)==0) {
792  dsetName=0 ;
793  }
794 
795  RooLinkedList& dataList = embedded ? _embeddedDataList : _dataList ;
796 
797  // Check that no dataset with target name already exists
798  if (dsetName && dataList.FindObject(dsetName)) {
799  coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << dsetName << " already exists in workspace, import aborted" << endl ;
800  return kTRUE ;
801  }
802  if (!dsetName && dataList.FindObject(inData.GetName())) {
803  coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << inData.GetName() << " already exists in workspace, import aborted" << endl ;
804  return kTRUE ;
805  }
806 
807  // Rename dataset if required
808  RooAbsData* clone ;
809  if (dsetName) {
810  if (!silence)
811  coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset from " << inData.GetName() << " to " << dsetName << endl ;
812  clone = (RooAbsData*) inData.Clone(dsetName) ;
813  } else {
814  clone = (RooAbsData*) inData.Clone(inData.GetName()) ;
815  }
816 
817 
818  // Process any change in variable names
819  if (strlen(varChangeIn)>0) {
820  // Parse comma separated lists of variable name changes
821  const std::vector<std::string> tokIn = RooHelpers::tokenise(varChangeIn, ",");
822  const std::vector<std::string> tokOut = RooHelpers::tokenise(varChangeOut, ",");
823  for (unsigned int i=0; i < tokIn.size(); ++i) {
824  if (!silence)
825  coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset observable " << tokIn[i] << " to " << tokOut[i] << endl ;
826  clone->changeObservableName(tokIn[i].c_str(), tokOut[i].c_str());
827  }
828  }
829 
830  // Now import the dataset observables, unless dataset is embedded
831  RooAbsArg* carg ;
832  if (!embedded) {
833  TIterator* iter = clone->get()->createIterator() ;
834  while((carg=(RooAbsArg*)iter->Next())) {
835  if (!arg(carg->GetName())) {
836  import(*carg) ;
837  }
838  }
839  delete iter ;
840  }
841 
842  dataList.Add(clone) ;
843  if (_dir) {
844  _dir->InternalAppend(clone) ;
845  }
846  if (_doExport) {
847  exportObj(clone) ;
848  }
849 
850  // Set expensive object cache of dataset internal buffers to that of workspace
851  RooFIter iter2 = clone->get()->fwdIterator() ;
852  while ((carg=iter2.next())) {
853  carg->setExpensiveObjectCache(expensiveObjectCache()) ;
854  }
855 
856 
857  return kFALSE ;
858 }
859 
860 
861 
862 
863 ////////////////////////////////////////////////////////////////////////////////
864 /// Define a named RooArgSet with given constituents. If importMissing is true, any constituents
865 /// of aset that are not in the workspace will be imported, otherwise an error is returned
866 /// for missing components
867 
868 Bool_t RooWorkspace::defineSet(const char* name, const RooArgSet& aset, Bool_t importMissing)
869 {
870  // Check if set was previously defined, if so print warning
871  map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
872  if (i!=_namedSets.end()) {
873  coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << endl ;
874  }
875 
876  RooArgSet wsargs ;
877 
878  // Check all constituents of provided set
879  TIterator* iter = aset.createIterator() ;
880  RooAbsArg* sarg ;
881  while((sarg=(RooAbsArg*)iter->Next())) {
882  // If missing, either import or report error
883  if (!arg(sarg->GetName())) {
884  if (importMissing) {
885  import(*sarg) ;
886  } else {
887  coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR set constituent \"" << sarg->GetName()
888  << "\" is not in workspace and importMissing option is disabled" << endl ;
889  return kTRUE ;
890  }
891  }
892  wsargs.add(*arg(sarg->GetName())) ;
893  }
894  delete iter ;
895 
896  // Install named set
897  _namedSets[name].removeAll() ;
898  _namedSets[name].add(wsargs) ;
899 
900  return kFALSE ;
901 }
902 
903 //_____________________________________________________________________________
904 Bool_t RooWorkspace::defineSetInternal(const char *name, const RooArgSet &aset)
905 {
906  // Define a named RooArgSet with given constituents. If importMissing is true, any constituents
907  // of aset that are not in the workspace will be imported, otherwise an error is returned
908  // for missing components
909 
910  // Check if set was previously defined, if so print warning
911  map<string, RooArgSet>::iterator i = _namedSets.find(name);
912  if (i != _namedSets.end()) {
913  coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName()
914  << ") WARNING redefining previously defined named set " << name << endl;
915  }
916 
917  // Install named set
918  _namedSets[name].removeAll();
919  _namedSets[name].add(aset);
920 
921  return kFALSE;
922 }
923 
924 ////////////////////////////////////////////////////////////////////////////////
925 /// Define a named set in the work space through a comma separated list of
926 /// names of objects already in the workspace
927 
928 Bool_t RooWorkspace::defineSet(const char* name, const char* contentList)
929 {
930  // Check if set was previously defined, if so print warning
931  map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
932  if (i!=_namedSets.end()) {
933  coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << endl ;
934  }
935 
936  RooArgSet wsargs ;
937 
938  // Check all constituents of provided set
939  for (const std::string& token : RooHelpers::tokenise(contentList, ",")) {
940  // If missing, either import or report error
941  if (!arg(token.c_str())) {
942  coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
943  << "\" is not in workspace" << endl ;
944  return kTRUE ;
945  }
946  wsargs.add(*arg(token.c_str())) ;
947  }
948 
949  // Install named set
950  _namedSets[name].removeAll() ;
951  _namedSets[name].add(wsargs) ;
952 
953  return kFALSE ;
954 }
955 
956 
957 
958 
959 ////////////////////////////////////////////////////////////////////////////////
960 /// Define a named set in the work space through a comma separated list of
961 /// names of objects already in the workspace
962 
963 Bool_t RooWorkspace::extendSet(const char* name, const char* newContents)
964 {
965  RooArgSet wsargs ;
966 
967  // Check all constituents of provided set
968  for (const std::string& token : RooHelpers::tokenise(newContents, ",")) {
969  // If missing, either import or report error
970  if (!arg(token.c_str())) {
971  coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
972  << "\" is not in workspace" << endl ;
973  return kTRUE ;
974  }
975  wsargs.add(*arg(token.c_str())) ;
976  }
977 
978  // Extend named set
979  _namedSets[name].add(wsargs,kTRUE) ;
980 
981  return kFALSE ;
982 }
983 
984 
985 
986 ////////////////////////////////////////////////////////////////////////////////
987 /// Return pointer to previously defined named set with given nmame
988 /// If no such set is found a null pointer is returned
989 
990 const RooArgSet* RooWorkspace::set(const char* name)
991 {
992  map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
993  return (i!=_namedSets.end()) ? &(i->second) : 0 ;
994 }
995 
996 
997 
998 
999 ////////////////////////////////////////////////////////////////////////////////
1000 /// Rename set to a new name
1001 
1002 Bool_t RooWorkspace::renameSet(const char* name, const char* newName)
1003 {
1004  // First check if set exists
1005  if (!set(name)) {
1006  coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << name
1007  << " does not exist" << endl ;
1008  return kTRUE ;
1009  }
1010 
1011  // Check if no set exists with new name
1012  if (set(newName)) {
1013  coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << newName
1014  << " already exists" << endl ;
1015  return kTRUE ;
1016  }
1017 
1018  // Copy entry under 'name' to 'newName'
1019  _namedSets[newName].add(_namedSets[name]) ;
1020 
1021  // Remove entry under old name
1022  _namedSets.erase(name) ;
1023 
1024  return kFALSE ;
1025 }
1026 
1027 
1028 
1029 
1030 ////////////////////////////////////////////////////////////////////////////////
1031 /// Remove a named set from the workspace
1032 
1033 Bool_t RooWorkspace::removeSet(const char* name)
1034 {
1035  // First check if set exists
1036  if (!set(name)) {
1037  coutE(InputArguments) << "RooWorkspace::removeSet(" << GetName() << ") ERROR a set with name " << name
1038  << " does not exist" << endl ;
1039  return kTRUE ;
1040  }
1041 
1042  // Remove set with given name
1043  _namedSets.erase(name) ;
1044 
1045  return kFALSE ;
1046 }
1047 
1048 
1049 
1050 
1051 ////////////////////////////////////////////////////////////////////////////////
1052 /// Open an import transaction operations. Returns kTRUE if successful, kFALSE
1053 /// if there is already an ongoing transaction
1054 
1055 Bool_t RooWorkspace::startTransaction()
1056 {
1057  // Check that there was no ongoing transaction
1058  if (_openTrans) {
1059  return kFALSE ;
1060  }
1061 
1062  // Open transaction
1063  _openTrans = kTRUE ;
1064  return kTRUE ;
1065 }
1066 
1067 
1068 
1069 
1070 ////////////////////////////////////////////////////////////////////////////////
1071 /// Cancel an ongoing import transaction. All objects imported since startTransaction()
1072 /// will be removed and the transaction will be terminated. Return kTRUE if cancel operation
1073 /// succeeds, return kFALSE if there was no open transaction
1074 
1075 Bool_t RooWorkspace::cancelTransaction()
1076 {
1077  // Check that there is an ongoing transaction
1078  if (!_openTrans) {
1079  return kFALSE ;
1080  }
1081 
1082  // Delete all objects in the sandbox
1083  TIterator* iter = _sandboxNodes.createIterator() ;
1084  RooAbsArg* tmpArg ;
1085  while((tmpArg=(RooAbsArg*)iter->Next())) {
1086  _allOwnedNodes.remove(*tmpArg) ;
1087  }
1088  delete iter ;
1089  _sandboxNodes.removeAll() ;
1090 
1091  // Mark transaction as finished
1092  _openTrans = kFALSE ;
1093 
1094  return kTRUE ;
1095 }
1096 
1097 Bool_t RooWorkspace::commitTransaction()
1098 {
1099  // Commit an ongoing import transaction. Returns kTRUE if commit succeeded,
1100  // return kFALSE if there was no ongoing transaction
1101 
1102  // Check that there is an ongoing transaction
1103  if (!_openTrans) {
1104  return kFALSE ;
1105  }
1106 
1107  // Publish sandbox nodes in directory and/or CINT if requested
1108  TIterator* iter = _sandboxNodes.createIterator() ;
1109  RooAbsArg* sarg ;
1110  while((sarg=(RooAbsArg*)iter->Next())) {
1111  if (_dir && sarg->IsA() != RooConstVar::Class()) {
1112  _dir->InternalAppend(sarg) ;
1113  }
1114  if (_doExport && sarg->IsA() != RooConstVar::Class()) {
1115  exportObj(sarg) ;
1116  }
1117  }
1118  delete iter ;
1119 
1120  // Remove all committed objects from the sandbox
1121  _sandboxNodes.removeAll() ;
1122 
1123  // Mark transaction as finished
1124  _openTrans = kFALSE ;
1125 
1126  return kTRUE ;
1127 }
1128 
1129 
1130 
1131 
1132 ////////////////////////////////////////////////////////////////////////////////
1133 
1134 Bool_t RooWorkspace::importClassCode(TClass* theClass, Bool_t doReplace)
1135 {
1136  return _classes.autoImportClass(theClass,doReplace) ;
1137 }
1138 
1139 
1140 
1141 ////////////////////////////////////////////////////////////////////////////////
1142 /// Inport code of all classes in the workspace that have a class name
1143 /// that matches pattern 'pat' and which are not found to be part of
1144 /// the standard ROOT distribution. If doReplace is true any existing
1145 /// class code saved in the workspace is replaced
1146 
1147 Bool_t RooWorkspace::importClassCode(const char* pat, Bool_t doReplace)
1148 {
1149  Bool_t ret(kTRUE) ;
1150 
1151  TRegexp re(pat,kTRUE) ;
1152  TIterator* iter = componentIterator() ;
1153  RooAbsArg* carg ;
1154  while((carg=(RooAbsArg*)iter->Next())) {
1155  TString className = carg->IsA()->GetName() ;
1156  if (className.Index(re)>=0 && !_classes.autoImportClass(carg->IsA(),doReplace)) {
1157  coutW(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") WARNING: problems import class code of object "
1158  << carg->IsA()->GetName() << "::" << carg->GetName() << ", reading of workspace will require external definition of class" << endl ;
1159  ret = kFALSE ;
1160  }
1161  }
1162  delete iter ;
1163 
1164  return ret ;
1165 }
1166 
1167 
1168 
1169 
1170 
1171 ////////////////////////////////////////////////////////////////////////////////
1172 /// Save snapshot of values and attributes (including "Constant") of given parameters.
1173 /// \param[in] name Name of the snapshot.
1174 /// \param[in] paramNames Comma-separated list of parameter names to be snapshot.
1175 Bool_t RooWorkspace::saveSnapshot(const char* name, const char* paramNames)
1176 {
1177  return saveSnapshot(name,argSet(paramNames),kFALSE) ;
1178 }
1179 
1180 
1181 
1182 
1183 
1184 ////////////////////////////////////////////////////////////////////////////////
1185 /// Save snapshot of values and attributes (including "Constant") of parameters 'params'.
1186 /// If importValues is FALSE, the present values from the object in the workspace are
1187 /// saved. If importValues is TRUE, the values of the objects passed in the 'params'
1188 /// argument are saved
1189 
1190 Bool_t RooWorkspace::saveSnapshot(const char* name, const RooArgSet& params, Bool_t importValues)
1191 {
1192  RooArgSet* actualParams = (RooArgSet*) _allOwnedNodes.selectCommon(params) ;
1193  RooArgSet* snapshot = (RooArgSet*) actualParams->snapshot() ;
1194  delete actualParams ;
1195 
1196  snapshot->setName(name) ;
1197 
1198  if (importValues) {
1199  *snapshot = params ;
1200  }
1201 
1202  RooArgSet* oldSnap = (RooArgSet*) _snapshots.FindObject(name) ;
1203  if (oldSnap) {
1204  coutI(ObjectHandling) << "RooWorkspace::saveSnaphot(" << GetName() << ") replacing previous snapshot with name " << name << endl ;
1205  _snapshots.Remove(oldSnap) ;
1206  delete oldSnap ;
1207  }
1208 
1209  _snapshots.Add(snapshot) ;
1210 
1211  return kTRUE ;
1212 }
1213 
1214 
1215 
1216 
1217 ////////////////////////////////////////////////////////////////////////////////
1218 /// Load the values and attributes of the parameters in the snapshot saved with
1219 /// the given name
1220 
1221 Bool_t RooWorkspace::loadSnapshot(const char* name)
1222 {
1223  RooArgSet* snap = (RooArgSet*) _snapshots.find(name) ;
1224  if (!snap) {
1225  coutE(ObjectHandling) << "RooWorkspace::loadSnapshot(" << GetName() << ") no snapshot with name " << name << " is available" << endl ;
1226  return kFALSE ;
1227  }
1228 
1229  RooArgSet* actualParams = (RooArgSet*) _allOwnedNodes.selectCommon(*snap) ;
1230  *actualParams = *snap ;
1231  delete actualParams ;
1232 
1233  return kTRUE ;
1234 }
1235 
1236 
1237 ////////////////////////////////////////////////////////////////////////////////
1238 /// Return the RooArgSet containing a snapshot of variables contained in the workspace
1239 ///
1240 /// Note that the variables of the objects in the snapshots are **copies** of the
1241 /// variables in the workspace. To load the values of a snapshot in the workspace
1242 /// variables, use loadSnapshot() instead.
1243 
1244 const RooArgSet* RooWorkspace::getSnapshot(const char* name) const
1245 {
1246  RooArgSet* snap = (RooArgSet*) _snapshots.find(name) ;
1247  if (!snap) {
1248  coutE(ObjectHandling) << "RooWorkspace::loadSnapshot(" << GetName() << ") no snapshot with name " << name << " is available" << endl ;
1249  return 0 ;
1250  }
1251 
1252  return snap ;
1253 }
1254 
1255 
1256 
1257 // //_____________________________________________________________________________
1258 // RooAbsPdf* RooWorkspace::joinPdf(const char* jointPdfName, const char* indexName, const char* inputMapping)
1259 // {
1260 // // Join given list of p.d.f.s into a simultaneous p.d.f with given name. If the named index category
1261 // // does not exist, it is created.
1262 // //
1263 // // Example : joinPdf("simPdf","expIndex","A=pdfA,B=pdfB") ;
1264 // //
1265 // // will return a RooSimultaneous named 'simPdf' with index category 'expIndex' with
1266 // // state names A and B. Pdf 'pdfA' will be associated with state A and pdf 'pdfB'
1267 // // will be associated with state B
1268 // //
1269 // return 0 ;
1270 // }
1271 
1272 // //_____________________________________________________________________________
1273 // RooAbsData* RooWorkspace::joinData(const char* jointDataName, const char* indexName, const char* inputMapping)
1274 // {
1275 // // Join given list of dataset into a joint dataset for use with a simultaneous pdf
1276 // // (as e.g. created by joingPdf"
1277 // //
1278 // // Example : joinData("simData","expIndex","A=dataA,B=dataB") ;
1279 // //
1280 // // will return a RooDataSet named 'simData' that consist of the entries of both
1281 // // dataA and dataB. An extra category column 'expIndex' is added that labels
1282 // // each entry with state 'A' and 'B' to indicate the originating dataset
1283 // return 0 ;
1284 // }
1285 
1286 
1287 ////////////////////////////////////////////////////////////////////////////////
1288 /// Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found
1289 
1290 RooAbsPdf* RooWorkspace::pdf(const char* name) const
1291 {
1292  return dynamic_cast<RooAbsPdf*>(_allOwnedNodes.find(name)) ;
1293 }
1294 
1295 
1296 ////////////////////////////////////////////////////////////////////////////////
1297 /// Retrieve function (RooAbsReal) with given name. Note that all RooAbsPdfs are also RooAbsReals. A null pointer is returned if not found.
1298 
1299 RooAbsReal* RooWorkspace::function(const char* name) const
1300 {
1301  return dynamic_cast<RooAbsReal*>(_allOwnedNodes.find(name)) ;
1302 }
1303 
1304 
1305 ////////////////////////////////////////////////////////////////////////////////
1306 /// Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found
1307 
1308 RooRealVar* RooWorkspace::var(const char* name) const
1309 {
1310  return dynamic_cast<RooRealVar*>(_allOwnedNodes.find(name)) ;
1311 }
1312 
1313 
1314 ////////////////////////////////////////////////////////////////////////////////
1315 /// Retrieve discrete variable (RooCategory) with given name. A null pointer is returned if not found
1316 
1317 RooCategory* RooWorkspace::cat(const char* name) const
1318 {
1319  return dynamic_cast<RooCategory*>(_allOwnedNodes.find(name)) ;
1320 }
1321 
1322 
1323 ////////////////////////////////////////////////////////////////////////////////
1324 /// Retrieve discrete function (RooAbsCategory) with given name. A null pointer is returned if not found
1325 
1326 RooAbsCategory* RooWorkspace::catfunc(const char* name) const
1327 {
1328  return dynamic_cast<RooAbsCategory*>(_allOwnedNodes.find(name)) ;
1329 }
1330 
1331 
1332 
1333 ////////////////////////////////////////////////////////////////////////////////
1334 /// Return RooAbsArg with given name. A null pointer is returned if none is found.
1335 
1336 RooAbsArg* RooWorkspace::arg(const char* name) const
1337 {
1338  return _allOwnedNodes.find(name) ;
1339 }
1340 
1341 
1342 
1343 ////////////////////////////////////////////////////////////////////////////////
1344 /// Return set of RooAbsArgs matching to given list of names
1345 
1346 RooArgSet RooWorkspace::argSet(const char* nameList) const
1347 {
1348  RooArgSet ret ;
1349 
1350  for (const std::string& token : RooHelpers::tokenise(nameList, ",")) {
1351  RooAbsArg* oneArg = arg(token.c_str()) ;
1352  if (oneArg) {
1353  ret.add(*oneArg) ;
1354  } else {
1355  coutE(InputArguments) << " RooWorkspace::argSet(" << GetName() << ") no RooAbsArg named \"" << token << "\" in workspace" << endl ;
1356  }
1357  }
1358  return ret ;
1359 }
1360 
1361 
1362 
1363 ////////////////////////////////////////////////////////////////////////////////
1364 /// Return fundamental (i.e. non-derived) RooAbsArg with given name. Fundamental types
1365 /// are e.g. RooRealVar, RooCategory. A null pointer is returned if none is found.
1366 
1367 RooAbsArg* RooWorkspace::fundArg(const char* name) const
1368 {
1369  RooAbsArg* tmp = arg(name) ;
1370  if (!tmp) {
1371  return 0 ;
1372  }
1373  return tmp->isFundamental() ? tmp : 0 ;
1374 }
1375 
1376 
1377 
1378 ////////////////////////////////////////////////////////////////////////////////
1379 /// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1380 
1381 RooAbsData* RooWorkspace::data(const char* name) const
1382 {
1383  return (RooAbsData*)_dataList.FindObject(name) ;
1384 }
1385 
1386 
1387 ////////////////////////////////////////////////////////////////////////////////
1388 /// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1389 
1390 RooAbsData* RooWorkspace::embeddedData(const char* name) const
1391 {
1392  return (RooAbsData*)_embeddedDataList.FindObject(name) ;
1393 }
1394 
1395 
1396 
1397 
1398 ////////////////////////////////////////////////////////////////////////////////
1399 /// Return set with all variable objects
1400 
1401 RooArgSet RooWorkspace::allVars() const
1402 {
1403  RooArgSet ret ;
1404 
1405  // Split list of components in pdfs, functions and variables
1406  TIterator* iter = _allOwnedNodes.createIterator() ;
1407  RooAbsArg* parg ;
1408  while((parg=(RooAbsArg*)iter->Next())) {
1409  if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1410  ret.add(*parg) ;
1411  }
1412  }
1413  delete iter ;
1414 
1415  return ret ;
1416 }
1417 
1418 
1419 ////////////////////////////////////////////////////////////////////////////////
1420 /// Return set with all category objects
1421 
1422 RooArgSet RooWorkspace::allCats() const
1423 {
1424  RooArgSet ret ;
1425 
1426  // Split list of components in pdfs, functions and variables
1427  TIterator* iter = _allOwnedNodes.createIterator() ;
1428  RooAbsArg* parg ;
1429  while((parg=(RooAbsArg*)iter->Next())) {
1430  if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
1431  ret.add(*parg) ;
1432  }
1433  }
1434  delete iter ;
1435 
1436  return ret ;
1437 }
1438 
1439 
1440 
1441 ////////////////////////////////////////////////////////////////////////////////
1442 /// Return set with all function objects
1443 
1444 RooArgSet RooWorkspace::allFunctions() const
1445 {
1446  RooArgSet ret ;
1447 
1448  // Split list of components in pdfs, functions and variables
1449  TIterator* iter = _allOwnedNodes.createIterator() ;
1450  RooAbsArg* parg ;
1451  while((parg=(RooAbsArg*)iter->Next())) {
1452  if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
1453  !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1454  !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
1455  !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1456  ret.add(*parg) ;
1457  }
1458  }
1459 
1460  return ret ;
1461 }
1462 
1463 
1464 ////////////////////////////////////////////////////////////////////////////////
1465 /// Return set with all category function objects
1466 
1467 RooArgSet RooWorkspace::allCatFunctions() const
1468 {
1469  RooArgSet ret ;
1470 
1471  // Split list of components in pdfs, functions and variables
1472  TIterator* iter = _allOwnedNodes.createIterator() ;
1473  RooAbsArg* parg ;
1474  while((parg=(RooAbsArg*)iter->Next())) {
1475  if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
1476  !parg->IsA()->InheritsFrom(RooCategory::Class())) {
1477  ret.add(*parg) ;
1478  }
1479  }
1480  return ret ;
1481 }
1482 
1483 
1484 
1485 ////////////////////////////////////////////////////////////////////////////////
1486 /// Return set with all resolution model objects
1487 
1488 RooArgSet RooWorkspace::allResolutionModels() const
1489 {
1490  RooArgSet ret ;
1491 
1492  // Split list of components in pdfs, functions and variables
1493  TIterator* iter = _allOwnedNodes.createIterator() ;
1494  RooAbsArg* parg ;
1495  while((parg=(RooAbsArg*)iter->Next())) {
1496  if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1497  if (!((RooResolutionModel*)parg)->isConvolved()) {
1498  ret.add(*parg) ;
1499  }
1500  }
1501  }
1502  return ret ;
1503 }
1504 
1505 
1506 ////////////////////////////////////////////////////////////////////////////////
1507 /// Return set with all probability density function objects
1508 
1509 RooArgSet RooWorkspace::allPdfs() const
1510 {
1511  RooArgSet ret ;
1512 
1513  // Split list of components in pdfs, functions and variables
1514  TIterator* iter = _allOwnedNodes.createIterator() ;
1515  RooAbsArg* parg ;
1516  while((parg=(RooAbsArg*)iter->Next())) {
1517  if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1518  !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1519  ret.add(*parg) ;
1520  }
1521  }
1522  return ret ;
1523 }
1524 
1525 
1526 
1527 ////////////////////////////////////////////////////////////////////////////////
1528 /// Return list of all dataset in the workspace
1529 
1530 list<RooAbsData*> RooWorkspace::allData() const
1531 {
1532  list<RooAbsData*> ret ;
1533  TIterator* iter = _dataList.MakeIterator() ;
1534  RooAbsData* dat ;
1535  while((dat=(RooAbsData*)iter->Next())) {
1536  ret.push_back(dat) ;
1537  }
1538  delete iter ;
1539  return ret ;
1540 }
1541 
1542 
1543 ////////////////////////////////////////////////////////////////////////////////
1544 /// Return list of all dataset in the workspace
1545 
1546 list<RooAbsData*> RooWorkspace::allEmbeddedData() const
1547 {
1548  list<RooAbsData*> ret ;
1549  TIterator* iter = _embeddedDataList.MakeIterator() ;
1550  RooAbsData* dat ;
1551  while((dat=(RooAbsData*)iter->Next())) {
1552  ret.push_back(dat) ;
1553  }
1554  delete iter ;
1555  return ret ;
1556 }
1557 
1558 
1559 
1560 ////////////////////////////////////////////////////////////////////////////////
1561 /// Return list of all generic objects in the workspace
1562 
1563 list<TObject*> RooWorkspace::allGenericObjects() const
1564 {
1565  list<TObject*> ret ;
1566  TIterator* iter = _genObjects.MakeIterator() ;
1567  TObject* gobj ;
1568  while((gobj=(RooAbsData*)iter->Next())) {
1569 
1570  // If found object is wrapper, return payload
1571  if (gobj->IsA()==RooTObjWrap::Class()) {
1572  ret.push_back(((RooTObjWrap*)gobj)->obj()) ;
1573  } else {
1574  ret.push_back(gobj) ;
1575  }
1576  }
1577  delete iter ;
1578  return ret ;
1579 }
1580 
1581 
1582 
1583 
1584 ////////////////////////////////////////////////////////////////////////////////
1585 /// Import code of class 'tc' into the repository. If code is already in repository it is only imported
1586 /// again if doReplace is false. The names and location of the source files is determined from the information
1587 /// in TClass. If no location is found in the TClass information, the files are searched in the workspace
1588 /// search path, defined by addClassDeclImportDir() and addClassImplImportDir() for declaration and implementation
1589 /// files respectively. If files cannot be found, abort with error status, otherwise update the internal
1590 /// class-to-file map and import the contents of the files, if they are not imported yet.
1591 
1592 Bool_t RooWorkspace::CodeRepo::autoImportClass(TClass* tc, Bool_t doReplace)
1593 {
1594 
1595  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") request to import code of class " << tc->GetName() << endl ;
1596 
1597  // *** PHASE 1 *** Check if file needs to be imported, or is in ROOT distribution, and check if it can be persisted
1598 
1599  // Check if we already have the class (i.e. it is in the classToFile map)
1600  if (!doReplace && _c2fmap.find(tc->GetName())!=_c2fmap.end()) {
1601  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " already imported, skipping" << endl ;
1602  return kTRUE ;
1603  }
1604 
1605  // Check if class is listed in a ROOTMAP file - if so we can skip it because it is in the root distribtion
1606  const char* mapEntry = gInterpreter->GetClassSharedLibs(tc->GetName()) ;
1607  if (mapEntry && strlen(mapEntry)>0) {
1608  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << endl ;
1609  return kTRUE ;
1610  }
1611 
1612  // Retrieve file names through ROOT TClass interface
1613  string implfile = tc->GetImplFileName() ;
1614  string declfile = tc->GetDeclFileName() ;
1615 
1616  // Check that file names are not empty
1617  if (implfile.empty() || declfile.empty()) {
1618  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") ERROR: cannot retrieve code file names for class "
1619  << tc->GetName() << " through ROOT TClass interface, unable to import code" << endl ;
1620  return kFALSE ;
1621  }
1622 
1623  // Check if header filename is found in ROOT distribution, if so, do not import class
1624  TString rootsys = gSystem->Getenv("ROOTSYS") ;
1625  if (TString(implfile.c_str()).Index(rootsys)>=0) {
1626  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << endl ;
1627  return kTRUE ;
1628  }
1629  const char* implpath=0 ;
1630 
1631  // Require that class meets technical criteria to be persistable (i.e it has a default ctor)
1632  // (We also need a default ctor of abstract classes, but cannot check that through is interface
1633  // as TClass::HasDefaultCtor only returns true for callable default ctors)
1634  if (!(tc->Property() & kIsAbstract) && !tc->HasDefaultConstructor()) {
1635  oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING cannot import class "
1636  << tc->GetName() << " : it cannot be persisted because it doesn't have a default constructor. Please fix " << endl ;
1637  return kFALSE ;
1638  }
1639 
1640 
1641  // *** PHASE 2 *** Check if declaration and implementation files can be located
1642 
1643  char* declpath = 0 ;
1644 
1645  // Check if header file can be found in specified location
1646  // If not, scan through list of 'class declaration' paths in RooWorkspace
1647  if (gSystem->AccessPathName(declfile.c_str())) {
1648 
1649  // Check list of additional declaration paths
1650  list<string>::iterator diter = RooWorkspace::_classDeclDirList.begin() ;
1651 
1652  while(diter!= RooWorkspace::_classDeclDirList.end()) {
1653 
1654  declpath = gSystem->ConcatFileName(diter->c_str(),declfile.c_str()) ;
1655  if (!gSystem->AccessPathName(declpath)) {
1656  // found declaration file
1657  break ;
1658  }
1659  // cleanup and continue ;
1660  delete[] declpath ;
1661  declpath=0 ;
1662 
1663  ++diter ;
1664  }
1665 
1666  // Header file cannot be found anywhere, warn user and abort operation
1667  if (!declpath) {
1668  oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1669  << tc->GetName() << " because header file " << declfile << " is not found in current directory nor in $ROOTSYS" ;
1670  if (_classDeclDirList.size()>0) {
1671  ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1672  diter = RooWorkspace::_classDeclDirList.begin() ;
1673 
1674  while(diter!= RooWorkspace::_classDeclDirList.end()) {
1675 
1676  if (diter!=RooWorkspace::_classDeclDirList.begin()) {
1677  ooccoutW(_wspace,ObjectHandling) << "," ;
1678  }
1679  ooccoutW(_wspace,ObjectHandling) << diter->c_str() ;
1680  ++diter ;
1681  }
1682  }
1683  ooccoutW(_wspace,ObjectHandling) << ". To fix this problem, add the required directory to the search "
1684  << "path using RooWorkspace::addClassDeclImportDir(const char* dir)" << endl ;
1685 
1686  return kFALSE ;
1687  }
1688  }
1689 
1690 
1691  // Check if implementation file can be found in specified location
1692  // If not, scan through list of 'class implementation' paths in RooWorkspace
1693  if (gSystem->AccessPathName(implfile.c_str())) {
1694 
1695  // Check list of additional declaration paths
1696  list<string>::iterator iiter = RooWorkspace::_classImplDirList.begin() ;
1697 
1698  while(iiter!= RooWorkspace::_classImplDirList.end()) {
1699 
1700  implpath = gSystem->ConcatFileName(iiter->c_str(),implfile.c_str()) ;
1701  if (!gSystem->AccessPathName(implpath)) {
1702  // found implementation file
1703  break ;
1704  }
1705  // cleanup and continue ;
1706  delete[] implpath ;
1707  implpath=0 ;
1708 
1709  ++iiter ;
1710  }
1711 
1712  // Implementation file cannot be found anywhere, warn user and abort operation
1713  if (!implpath) {
1714  oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1715  << tc->GetName() << " because implementation file " << implfile << " is not found in current directory nor in $ROOTSYS" ;
1716  if (_classDeclDirList.size()>0) {
1717  ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1718  iiter = RooWorkspace::_classImplDirList.begin() ;
1719 
1720  while(iiter!= RooWorkspace::_classImplDirList.end()) {
1721 
1722  if (iiter!=RooWorkspace::_classImplDirList.begin()) {
1723  ooccoutW(_wspace,ObjectHandling) << "," ;
1724  }
1725  ooccoutW(_wspace,ObjectHandling) << iiter->c_str() ;
1726  ++iiter ;
1727  }
1728  }
1729  ooccoutW(_wspace,ObjectHandling) << ". To fix this problem add the required directory to the search "
1730  << "path using RooWorkspace::addClassImplImportDir(const char* dir)" << endl ;
1731  return kFALSE ;
1732  }
1733  }
1734 
1735  char buf[64000];
1736 
1737  // *** Phase 3 *** Prepare to import code from files into STL string buffer
1738  //
1739  // Code storage is organized in two linked maps
1740  //
1741  // _fmap contains stl strings with code, indexed on declaration file name
1742  //
1743  // _c2fmap contains list of declaration file names and list of base classes
1744  // and is indexed on class name
1745  //
1746  // Phase 3 is skipped if fmap already contains an entry with given filebasename
1747 
1748  string declfilename = declpath?gSystem->BaseName(declpath):gSystem->BaseName(declfile.c_str()) ;
1749 
1750  // Split in base and extension
1751  int dotpos2 = strrchr(declfilename.c_str(),'.') - declfilename.c_str() ;
1752  string declfilebase = declfilename.substr(0,dotpos2) ;
1753  string declfileext = declfilename.substr(dotpos2+1) ;
1754 
1755  list<string> extraHeaders ;
1756 
1757  // If file has not beed stored yet, enter stl strings with implementation and declaration in file map
1758  if (_fmap.find(declfilebase) == _fmap.end()) {
1759 
1760  // Open declaration file
1761  fstream fdecl(declpath?declpath:declfile.c_str()) ;
1762 
1763  // Abort import if declaration file cannot be opened
1764  if (!fdecl) {
1765  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1766  << ") ERROR opening declaration file " << declfile << endl ;
1767  return kFALSE ;
1768  }
1769 
1770  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1771  << ") importing code of class " << tc->GetName()
1772  << " from " << (implpath?implpath:implfile.c_str())
1773  << " and " << (declpath?declpath:declfile.c_str()) << endl ;
1774 
1775 
1776  // Read entire file into an stl string
1777  string decl ;
1778  while(fdecl.getline(buf,1023)) {
1779 
1780  // Look for include state of self
1781  Bool_t processedInclude = kFALSE ;
1782  char* extincfile = 0 ;
1783 
1784  // Look for include of declaration file corresponding to this implementation file
1785  if (strstr(buf,"#include")) {
1786  // Process #include statements here
1787  char tmp[64000];
1788  strlcpy(tmp, buf, 64000);
1789  Bool_t stdinclude = strchr(buf, '<');
1790  strtok(tmp, " <\"");
1791  char *incfile = strtok(0, " <>\"");
1792 
1793  if (!stdinclude) {
1794  // check if it lives in $ROOTSYS/include
1795  TString hpath = gSystem->Getenv("ROOTSYS");
1796  hpath += "/include/";
1797  hpath += incfile;
1798  if (gSystem->AccessPathName(hpath.Data())) {
1799  oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1800  << ") scheduling include file " << incfile << " for import" << endl;
1801  extraHeaders.push_back(incfile);
1802  extincfile = incfile;
1803  processedInclude = kTRUE;
1804  }
1805  }
1806  }
1807 
1808  if (processedInclude) {
1809  decl += "// external include file below retrieved from workspace code storage\n" ;
1810  decl += Form("#include \"%s\"\n",extincfile) ;
1811  } else {
1812  decl += buf ;
1813  decl += '\n' ;
1814  }
1815  }
1816 
1817  // Open implementation file
1818  fstream fimpl(implpath?implpath:implfile.c_str()) ;
1819 
1820  // Abort import if implementation file cannot be opened
1821  if (!fimpl) {
1822  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1823  << ") ERROR opening implementation file " << implfile << endl ;
1824  return kFALSE ;
1825  }
1826 
1827 
1828  // Import entire implentation file into stl string
1829  string impl ;
1830  while(fimpl.getline(buf,1023)) {
1831  // Process #include statements here
1832 
1833  // Look for include state of self
1834  Bool_t foundSelfInclude=kFALSE ;
1835  Bool_t processedInclude = kFALSE ;
1836  char* extincfile = 0 ;
1837 
1838  // Look for include of declaration file corresponding to this implementation file
1839  if (strstr(buf,"#include")) {
1840  // Process #include statements here
1841  char tmp[64000];
1842  strlcpy(tmp, buf, 64000);
1843  Bool_t stdinclude = strchr(buf, '<');
1844  strtok(tmp, " <\"");
1845  char *incfile = strtok(0, " <>\"");
1846 
1847  if (strstr(incfile, declfilename.c_str())) {
1848  foundSelfInclude = kTRUE;
1849  }
1850 
1851  if (!stdinclude && !foundSelfInclude) {
1852  // check if it lives in $ROOTSYS/include
1853  TString hpath = gSystem->Getenv("ROOTSYS");
1854  hpath += "/include/";
1855  hpath += incfile;
1856 
1857  if (gSystem->AccessPathName(hpath.Data())) {
1858  oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1859  << ") scheduling include file " << incfile << " for import" << endl;
1860  extraHeaders.push_back(incfile);
1861  extincfile = incfile;
1862  processedInclude = kTRUE;
1863  }
1864  }
1865  }
1866 
1867  // Explicitly rewrite include of own declaration file to string
1868  // any directory prefixes, copy all other lines verbatim in stl string
1869  if (foundSelfInclude) {
1870  // If include of self is found, substitute original include
1871  // which may have directory structure with a plain include
1872  impl += "// class declaration include file below retrieved from workspace code storage\n" ;
1873  impl += Form("#include \"%s.%s\"\n",declfilebase.c_str(),declfileext.c_str()) ;
1874  } else if (processedInclude) {
1875  impl += "// external include file below retrieved from workspace code storage\n" ;
1876  impl += Form("#include \"%s\"\n",extincfile) ;
1877  } else {
1878  impl += buf ;
1879  impl += '\n' ;
1880  }
1881  }
1882 
1883  // Create entry in file map
1884  _fmap[declfilebase]._hfile = decl ;
1885  _fmap[declfilebase]._cxxfile = impl ;
1886  _fmap[declfilebase]._hext = declfileext ;
1887 
1888  // Process extra includes now
1889  for (list<string>::iterator ehiter = extraHeaders.begin() ; ehiter != extraHeaders.end() ; ++ehiter ) {
1890  if (_ehmap.find(*ehiter) == _ehmap.end()) {
1891 
1892  ExtraHeader eh ;
1893  eh._hname = ehiter->c_str() ;
1894  fstream fehdr(ehiter->c_str()) ;
1895  string ehimpl ;
1896  char buf2[1024] ;
1897  while(fehdr.getline(buf2,1023)) {
1898 
1899  // Look for include of declaration file corresponding to this implementation file
1900  if (strstr(buf2,"#include")) {
1901  // Process #include statements here
1902  char tmp[64000];
1903  strlcpy(tmp, buf2, 64000);
1904  Bool_t stdinclude = strchr(buf, '<');
1905  strtok(tmp, " <\"");
1906  char *incfile = strtok(0, " <>\"");
1907 
1908  if (!stdinclude) {
1909  // check if it lives in $ROOTSYS/include
1910  TString hpath = gSystem->Getenv("ROOTSYS");
1911  hpath += "/include/";
1912  hpath += incfile;
1913  if (gSystem->AccessPathName(hpath.Data())) {
1914  oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1915  << ") scheduling recursive include file " << incfile << " for import"
1916  << endl;
1917  extraHeaders.push_back(incfile);
1918  }
1919  }
1920  }
1921 
1922  ehimpl += buf2;
1923  ehimpl += '\n';
1924  }
1925  eh._hfile = ehimpl.c_str();
1926 
1927  _ehmap[ehiter->c_str()] = eh;
1928  }
1929  }
1930 
1931  } else {
1932 
1933  // Inform that existing file entry is being recycled because it already contained class code
1934  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1935  << ") code of class " << tc->GetName()
1936  << " was already imported from " << (implpath?implpath:implfile.c_str())
1937  << " and " << (declpath?declpath:declfile.c_str()) << endl ;
1938 
1939  }
1940 
1941 
1942  // *** PHASE 4 *** Import stl strings with code into workspace
1943  //
1944  // If multiple classes are declared in a single code unit, there will be
1945  // multiple _c2fmap entries all pointing to the same _fmap entry.
1946 
1947  // Make list of all immediate base classes of this class
1948  TString baseNameList ;
1949  TList* bl = tc->GetListOfBases() ;
1950  TIterator* iter = bl->MakeIterator() ;
1951  TBaseClass* base ;
1952  list<TClass*> bases ;
1953  while((base=(TBaseClass*)iter->Next())) {
1954  if (baseNameList.Length()>0) {
1955  baseNameList += "," ;
1956  }
1957  baseNameList += base->GetClassPointer()->GetName() ;
1958  bases.push_back(base->GetClassPointer()) ;
1959  }
1960 
1961  // Map class name to above _fmap entries, along with list of base classes
1962  // in _c2fmap
1963  _c2fmap[tc->GetName()]._baseName = baseNameList ;
1964  _c2fmap[tc->GetName()]._fileBase = declfilebase ;
1965 
1966  // Recursive store all base classes.
1967  list<TClass*>::iterator biter = bases.begin() ;
1968  while(biter!=bases.end()) {
1969  autoImportClass(*biter,doReplace) ;
1970  ++biter ;
1971  }
1972 
1973  // Cleanup
1974  if (implpath) {
1975  delete[] implpath ;
1976  }
1977  if (declpath) {
1978  delete[] declpath ;
1979  }
1980 
1981 
1982  return kTRUE ;
1983 }
1984 
1985 
1986 ////////////////////////////////////////////////////////////////////////////////
1987 /// Create transient TDirectory representation of this workspace. This directory
1988 /// will appear as a subdirectory of the directory that contains the workspace
1989 /// and will have the name of the workspace suffixed with "Dir". The TDirectory
1990 /// interface is read-only. Any attempt to insert objects into the workspace
1991 /// directory representation will result in an error message. Note that some
1992 /// ROOT object like TH1 automatically insert themselves into the current directory
1993 /// when constructed. This will give error messages when done in a workspace
1994 /// directory.
1995 
1996 Bool_t RooWorkspace::makeDir()
1997 {
1998  if (_dir) return kTRUE ;
1999 
2000  TString title= Form("TDirectory representation of RooWorkspace %s",GetName()) ;
2001  _dir = new WSDir(GetName(),title.Data(),this) ;
2002 
2003  TIterator* iter = componentIterator() ;
2004  RooAbsArg* darg ;
2005  while((darg=(RooAbsArg*)iter->Next())) {
2006  if (darg->IsA() != RooConstVar::Class()) {
2007  _dir->InternalAppend(darg) ;
2008  }
2009  }
2010 
2011  return kTRUE ;
2012 }
2013 
2014 
2015 
2016 ////////////////////////////////////////////////////////////////////////////////
2017 /// Import a clone of a generic TObject into workspace generic object container. Imported
2018 /// object can be retrieved by name through the obj() method. The object is cloned upon
2019 /// importation and the input argument does not need to live beyond the import call
2020 ///
2021 /// Returns kTRUE if an error has occurred.
2022 
2023 Bool_t RooWorkspace::import(TObject& object, Bool_t replaceExisting)
2024 {
2025  // First check if object with given name already exists
2026  TObject* oldObj = _genObjects.FindObject(object.GetName()) ;
2027  if (oldObj && !replaceExisting) {
2028  coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
2029  << object.GetName() << " is already in workspace and replaceExisting flag is set to false" << endl ;
2030  return kTRUE ;
2031  }
2032 
2033  // Grab the current state of the directory Auto-Add
2034  ROOT::DirAutoAdd_t func = object.IsA()->GetDirectoryAutoAdd();
2035  object.IsA()->SetDirectoryAutoAdd(0);
2036  Bool_t tmp = RooPlot::setAddDirectoryStatus(kFALSE) ;
2037 
2038  if (oldObj) {
2039  _genObjects.Replace(oldObj,object.Clone()) ;
2040  delete oldObj ;
2041  } else {
2042  _genObjects.Add(object.Clone()) ;
2043  }
2044 
2045  // Reset the state of the directory Auto-Add
2046  object.IsA()->SetDirectoryAutoAdd(func);
2047  RooPlot::setAddDirectoryStatus(tmp) ;
2048 
2049  return kFALSE ;
2050 }
2051 
2052 
2053 
2054 
2055 ////////////////////////////////////////////////////////////////////////////////
2056 /// Import a clone of a generic TObject into workspace generic object container.
2057 /// The imported object will be stored under the given alias name rather than its
2058 /// own name. Imported object can be retrieved its alias name through the obj() method.
2059 /// The object is cloned upon importation and the input argument does not need to live beyond the import call
2060 /// This method is mostly useful for importing objects that do not have a settable name such as TMatrix
2061 ///
2062 /// Returns kTRUE if an error has occurred.
2063 
2064 Bool_t RooWorkspace::import(TObject& object, const char* aliasName, Bool_t replaceExisting)
2065 {
2066  // First check if object with given name already exists
2067  TObject* oldObj = _genObjects.FindObject(object.GetName()) ;
2068  if (oldObj && !replaceExisting) {
2069  coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
2070  << object.GetName() << " is already in workspace and replaceExisting flag is set to false" << endl ;
2071  return kTRUE ;
2072  }
2073 
2074  TH1::AddDirectory(kFALSE) ;
2075  RooTObjWrap* wrapper = new RooTObjWrap(object.Clone()) ;
2076  TH1::AddDirectory(kTRUE) ;
2077  wrapper->setOwning(kTRUE) ;
2078  wrapper->SetName(aliasName) ;
2079  wrapper->SetTitle(aliasName) ;
2080 
2081  if (oldObj) {
2082  _genObjects.Replace(oldObj,wrapper) ;
2083  delete oldObj ;
2084  } else {
2085  _genObjects.Add(wrapper) ;
2086  }
2087  return kFALSE ;
2088 }
2089 
2090 
2091 
2092 
2093 ////////////////////////////////////////////////////////////////////////////////
2094 /// Insert RooStudyManager module
2095 
2096 Bool_t RooWorkspace::addStudy(RooAbsStudy& study)
2097 {
2098  RooAbsStudy* clone = (RooAbsStudy*) study.Clone() ;
2099  _studyMods.Add(clone) ;
2100  return kFALSE ;
2101 }
2102 
2103 
2104 
2105 
2106 ////////////////////////////////////////////////////////////////////////////////
2107 /// Remove all RooStudyManager modules
2108 
2109 void RooWorkspace::clearStudies()
2110 {
2111  _studyMods.Delete() ;
2112 }
2113 
2114 
2115 
2116 
2117 ////////////////////////////////////////////////////////////////////////////////
2118 /// Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
2119 
2120 TObject* RooWorkspace::obj(const char* name) const
2121 {
2122  // Try RooAbsArg first
2123  TObject* ret = arg(name) ;
2124  if (ret) return ret ;
2125 
2126  // Then try RooAbsData
2127  ret = data(name) ;
2128  if (ret) return ret ;
2129 
2130  // Finally try generic object store
2131  return genobj(name) ;
2132 }
2133 
2134 
2135 
2136 ////////////////////////////////////////////////////////////////////////////////
2137 /// Return generic object with given name
2138 
2139 TObject* RooWorkspace::genobj(const char* name) const
2140 {
2141  // Find object by name
2142  TObject* gobj = _genObjects.FindObject(name) ;
2143 
2144  // Exit here if not found
2145  if (!gobj) return 0 ;
2146 
2147  // If found object is wrapper, return payload
2148  if (gobj->IsA()==RooTObjWrap::Class()) return ((RooTObjWrap*)gobj)->obj() ;
2149 
2150  return gobj ;
2151 }
2152 
2153 
2154 
2155 ////////////////////////////////////////////////////////////////////////////////
2156 
2157 Bool_t RooWorkspace::cd(const char* path)
2158 {
2159  makeDir() ;
2160  return _dir->cd(path) ;
2161 }
2162 
2163 
2164 
2165 ////////////////////////////////////////////////////////////////////////////////
2166 /// Save this current workspace into given file
2167 
2168 Bool_t RooWorkspace::writeToFile(const char* fileName, Bool_t recreate)
2169 {
2170  TFile f(fileName,recreate?"RECREATE":"UPDATE") ;
2171  Write() ;
2172  return kFALSE ;
2173 }
2174 
2175 
2176 
2177 ////////////////////////////////////////////////////////////////////////////////
2178 /// Return instance to factory tool
2179 
2180 RooFactoryWSTool& RooWorkspace::factory()
2181 {
2182  if (_factory) {
2183  return *_factory;
2184  }
2185  cxcoutD(ObjectHandling) << "INFO: Creating RooFactoryWSTool associated with this workspace" << endl ;
2186  _factory = make_unique<RooFactoryWSTool>(*this);
2187  return *_factory;
2188 }
2189 
2190 
2191 
2192 
2193 ////////////////////////////////////////////////////////////////////////////////
2194 /// Short-hand function for factory()->process(expr) ;
2195 
2196 RooAbsArg* RooWorkspace::factory(const char* expr)
2197 {
2198  return factory().process(expr) ;
2199 }
2200 
2201 
2202 
2203 
2204 ////////////////////////////////////////////////////////////////////////////////
2205 /// Print contents of the workspace
2206 
2207 void RooWorkspace::Print(Option_t* opts) const
2208 {
2209  Bool_t treeMode(kFALSE) ;
2210  Bool_t verbose(kFALSE);
2211  if (TString(opts).Contains("t")) {
2212  treeMode=kTRUE ;
2213  }
2214  if (TString(opts).Contains("v")) {
2215  verbose = kTRUE;
2216  }
2217 
2218  cout << endl << "RooWorkspace(" << GetName() << ") " << GetTitle() << " contents" << endl << endl ;
2219 
2220  RooAbsArg* parg ;
2221 
2222  RooArgSet pdfSet ;
2223  RooArgSet funcSet ;
2224  RooArgSet varSet ;
2225  RooArgSet catfuncSet ;
2226  RooArgSet convResoSet ;
2227  RooArgSet resoSet ;
2228 
2229 
2230  // Split list of components in pdfs, functions and variables
2231  TIterator* iter = _allOwnedNodes.createIterator() ;
2232  while((parg=(RooAbsArg*)iter->Next())) {
2233 
2234  //---------------
2235 
2236  if (treeMode) {
2237 
2238  // In tree mode, only add nodes with no clients to the print lists
2239 
2240  if (parg->IsA()->InheritsFrom(RooAbsPdf::Class())) {
2241  if (!parg->hasClients()) {
2242  pdfSet.add(*parg) ;
2243  }
2244  }
2245 
2246  if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2247  !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2248  !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2249  !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2250  if (!parg->hasClients()) {
2251  funcSet.add(*parg) ;
2252  }
2253  }
2254 
2255 
2256  if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2257  !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2258  if (!parg->hasClients()) {
2259  catfuncSet.add(*parg) ;
2260  }
2261  }
2262 
2263  } else {
2264 
2265  if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2266  if (((RooResolutionModel*)parg)->isConvolved()) {
2267  convResoSet.add(*parg) ;
2268  } else {
2269  resoSet.add(*parg) ;
2270  }
2271  }
2272 
2273  if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2274  !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2275  pdfSet.add(*parg) ;
2276  }
2277 
2278  if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2279  !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2280  !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2281  !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2282  funcSet.add(*parg) ;
2283  }
2284 
2285  if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2286  !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2287  catfuncSet.add(*parg) ;
2288  }
2289  }
2290 
2291  if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2292  varSet.add(*parg) ;
2293  }
2294 
2295  if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
2296  varSet.add(*parg) ;
2297  }
2298 
2299  }
2300  delete iter ;
2301 
2302 
2303  RooFit::MsgLevel oldLevel = RooMsgService::instance().globalKillBelow() ;
2304  RooMsgService::instance().setGlobalKillBelow(RooFit::WARNING) ;
2305 
2306  if (varSet.getSize()>0) {
2307  varSet.sort() ;
2308  cout << "variables" << endl ;
2309  cout << "---------" << endl ;
2310  cout << varSet << endl ;
2311  cout << endl ;
2312  }
2313 
2314  if (pdfSet.getSize()>0) {
2315  cout << "p.d.f.s" << endl ;
2316  cout << "-------" << endl ;
2317  pdfSet.sort() ;
2318  iter = pdfSet.createIterator() ;
2319  while((parg=(RooAbsArg*)iter->Next())) {
2320  if (treeMode) {
2321  parg->printComponentTree() ;
2322  } else {
2323  parg->Print() ;
2324  }
2325  }
2326  delete iter ;
2327  cout << endl ;
2328  }
2329 
2330  if (!treeMode) {
2331  if (resoSet.getSize()>0) {
2332  cout << "analytical resolution models" << endl ;
2333  cout << "----------------------------" << endl ;
2334  resoSet.sort() ;
2335  iter = resoSet.createIterator() ;
2336  while((parg=(RooAbsArg*)iter->Next())) {
2337  parg->Print() ;
2338  }
2339  delete iter ;
2340  // iter = convResoSet.createIterator() ;
2341  // while((parg=(RooAbsArg*)iter->Next())) {
2342  // parg->Print() ;
2343  // }
2344  // delete iter ;
2345  cout << endl ;
2346  }
2347  }
2348 
2349  if (funcSet.getSize()>0) {
2350  cout << "functions" << endl ;
2351  cout << "--------" << endl ;
2352  funcSet.sort() ;
2353  iter = funcSet.createIterator() ;
2354  while((parg=(RooAbsArg*)iter->Next())) {
2355  if (treeMode) {
2356  parg->printComponentTree() ;
2357  } else {
2358  parg->Print() ;
2359  }
2360  }
2361  delete iter ;
2362  cout << endl ;
2363  }
2364 
2365  if (catfuncSet.getSize()>0) {
2366  cout << "category functions" << endl ;
2367  cout << "------------------" << endl ;
2368  catfuncSet.sort() ;
2369  iter = catfuncSet.createIterator() ;
2370  while((parg=(RooAbsArg*)iter->Next())) {
2371  if (treeMode) {
2372  parg->printComponentTree() ;
2373  } else {
2374  parg->Print() ;
2375  }
2376  }
2377  delete iter ;
2378  cout << endl ;
2379  }
2380 
2381  if (_dataList.GetSize()>0) {
2382  cout << "datasets" << endl ;
2383  cout << "--------" << endl ;
2384  iter = _dataList.MakeIterator() ;
2385  RooAbsData* data2 ;
2386  while((data2=(RooAbsData*)iter->Next())) {
2387  cout << data2->IsA()->GetName() << "::" << data2->GetName() << *data2->get() << endl ;
2388  }
2389  delete iter ;
2390  cout << endl ;
2391  }
2392 
2393  if (_embeddedDataList.GetSize()>0) {
2394  cout << "embedded datasets (in pdfs and functions)" << endl ;
2395  cout << "-----------------------------------------" << endl ;
2396  iter = _embeddedDataList.MakeIterator() ;
2397  RooAbsData* data2 ;
2398  while((data2=(RooAbsData*)iter->Next())) {
2399  cout << data2->IsA()->GetName() << "::" << data2->GetName() << *data2->get() << endl ;
2400  }
2401  delete iter ;
2402  cout << endl ;
2403  }
2404 
2405  if (_snapshots.GetSize()>0) {
2406  cout << "parameter snapshots" << endl ;
2407  cout << "-------------------" << endl ;
2408  iter = _snapshots.MakeIterator() ;
2409  RooArgSet* snap ;
2410  while((snap=(RooArgSet*)iter->Next())) {
2411  cout << snap->GetName() << " = (" ;
2412  TIterator* aiter = snap->createIterator() ;
2413  RooAbsArg* a ;
2414  Bool_t first(kTRUE) ;
2415  while((a=(RooAbsArg*)aiter->Next())) {
2416  if (first) { first=kFALSE ; } else { cout << "," ; }
2417  cout << a->GetName() << "=" ;
2418  a->printValue(cout) ;
2419  if (a->isConstant()) {
2420  cout << "[C]" ;
2421  }
2422  }
2423  cout << ")" << endl ;
2424  delete aiter ;
2425  }
2426  delete iter ;
2427  cout << endl ;
2428  }
2429 
2430 
2431  if (_namedSets.size()>0) {
2432  cout << "named sets" << endl ;
2433  cout << "----------" << endl ;
2434  for (map<string,RooArgSet>::const_iterator it = _namedSets.begin() ; it != _namedSets.end() ; ++it) {
2435  if (verbose || !TString(it->first.c_str()).BeginsWith("CACHE_")) {
2436  cout << it->first << ":" << it->second << endl;
2437  }
2438  }
2439 
2440  cout << endl ;
2441  }
2442 
2443 
2444  if (_genObjects.GetSize()>0) {
2445  cout << "generic objects" << endl ;
2446  cout << "---------------" << endl ;
2447  iter = _genObjects.MakeIterator() ;
2448  TObject* gobj ;
2449  while((gobj=(TObject*)iter->Next())) {
2450  if (gobj->IsA()==RooTObjWrap::Class()) {
2451  cout << ((RooTObjWrap*)gobj)->obj()->IsA()->GetName() << "::" << gobj->GetName() << endl ;
2452  } else {
2453  cout << gobj->IsA()->GetName() << "::" << gobj->GetName() << endl ;
2454  }
2455  }
2456  delete iter ;
2457  cout << endl ;
2458 
2459  }
2460 
2461  if (_studyMods.GetSize()>0) {
2462  cout << "study modules" << endl ;
2463  cout << "-------------" << endl ;
2464  iter = _studyMods.MakeIterator() ;
2465  TObject* smobj ;
2466  while((smobj=(TObject*)iter->Next())) {
2467  cout << smobj->IsA()->GetName() << "::" << smobj->GetName() << endl ;
2468  }
2469  delete iter ;
2470  cout << endl ;
2471 
2472  }
2473 
2474  if (_classes.listOfClassNames().size()>0) {
2475  cout << "embedded class code" << endl ;
2476  cout << "-------------------" << endl ;
2477  cout << _classes.listOfClassNames() << endl ;
2478  cout << endl ;
2479  }
2480 
2481  if (_eocache.size()>0) {
2482  cout << "embedded precalculated expensive components" << endl ;
2483  cout << "-------------------------------------------" << endl ;
2484  _eocache.print() ;
2485  }
2486 
2487  RooMsgService::instance().setGlobalKillBelow(oldLevel) ;
2488 
2489  return ;
2490 }
2491 
2492 
2493 ////////////////////////////////////////////////////////////////////////////////
2494 /// Custom streamer for the workspace. Stream contents of workspace
2495 /// and code repository. When reading, read code repository first
2496 /// and compile missing classes before proceeding with streaming
2497 /// of workspace contents to avoid errors.
2498 
2499 void RooWorkspace::CodeRepo::Streamer(TBuffer &R__b)
2500 {
2501  typedef ::RooWorkspace::CodeRepo thisClass;
2502 
2503  // Stream an object of class RooWorkspace::CodeRepo.
2504  if (R__b.IsReading()) {
2505 
2506  UInt_t R__s, R__c;
2507  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2508 
2509  // Stream contents of ClassFiles map
2510  Int_t count(0) ;
2511  R__b >> count ;
2512  while(count--) {
2513  TString name ;
2514  name.Streamer(R__b) ;
2515  _fmap[name]._hext.Streamer(R__b) ;
2516  _fmap[name]._hfile.Streamer(R__b) ;
2517  _fmap[name]._cxxfile.Streamer(R__b) ;
2518  }
2519 
2520  // Stream contents of ClassRelInfo map
2521  count=0 ;
2522  R__b >> count ;
2523  while(count--) {
2524  TString name ;
2525  name.Streamer(R__b) ;
2526  _c2fmap[name]._baseName.Streamer(R__b) ;
2527  _c2fmap[name]._fileBase.Streamer(R__b) ;
2528  }
2529 
2530  if (R__v==2) {
2531 
2532  count=0 ;
2533  R__b >> count ;
2534  while(count--) {
2535  TString name ;
2536  name.Streamer(R__b) ;
2537  _ehmap[name]._hname.Streamer(R__b) ;
2538  _ehmap[name]._hfile.Streamer(R__b) ;
2539  }
2540  }
2541 
2542  R__b.CheckByteCount(R__s, R__c, thisClass::IsA());
2543 
2544  // Instantiate any classes that are not defined in current session
2545  _compiledOK = !compileClasses() ;
2546 
2547  } else {
2548 
2549  UInt_t R__c;
2550  R__c = R__b.WriteVersion(thisClass::IsA(), kTRUE);
2551 
2552  // Stream contents of ClassFiles map
2553  UInt_t count = _fmap.size() ;
2554  R__b << count ;
2555  map<TString,ClassFiles>::iterator iter = _fmap.begin() ;
2556  while(iter!=_fmap.end()) {
2557  TString key_copy(iter->first) ;
2558  key_copy.Streamer(R__b) ;
2559  iter->second._hext.Streamer(R__b) ;
2560  iter->second._hfile.Streamer(R__b);
2561  iter->second._cxxfile.Streamer(R__b);
2562 
2563  ++iter ;
2564  }
2565 
2566  // Stream contents of ClassRelInfo map
2567  count = _c2fmap.size() ;
2568  R__b << count ;
2569  map<TString,ClassRelInfo>::iterator iter2 = _c2fmap.begin() ;
2570  while(iter2!=_c2fmap.end()) {
2571  TString key_copy(iter2->first) ;
2572  key_copy.Streamer(R__b) ;
2573  iter2->second._baseName.Streamer(R__b) ;
2574  iter2->second._fileBase.Streamer(R__b);
2575  ++iter2 ;
2576  }
2577 
2578  // Stream contents of ExtraHeader map
2579  count = _ehmap.size() ;
2580  R__b << count ;
2581  map<TString,ExtraHeader>::iterator iter3 = _ehmap.begin() ;
2582  while(iter3!=_ehmap.end()) {
2583  TString key_copy(iter3->first) ;
2584  key_copy.Streamer(R__b) ;
2585  iter3->second._hname.Streamer(R__b) ;
2586  iter3->second._hfile.Streamer(R__b);
2587  ++iter3 ;
2588  }
2589 
2590  R__b.SetByteCount(R__c, kTRUE);
2591 
2592  }
2593 }
2594 
2595 
2596 ////////////////////////////////////////////////////////////////////////////////
2597 /// Stream an object of class RooWorkspace. This is a standard ROOT streamer for the
2598 /// I/O part. This custom function exists to detach all external client links
2599 /// from the payload prior to writing the payload so that these client links
2600 /// are not persisted. (Client links occur if external function objects use
2601 /// objects contained in the workspace as input)
2602 /// After the actual writing, these client links are restored.
2603 
2604 void RooWorkspace::Streamer(TBuffer &R__b)
2605 {
2606  if (R__b.IsReading()) {
2607 
2608  R__b.ReadClassBuffer(RooWorkspace::Class(),this);
2609 
2610  // Perform any pass-2 schema evolution here
2611  RooFIter fiter = _allOwnedNodes.fwdIterator() ;
2612  RooAbsArg* node ;
2613  while((node=fiter.next())) {
2614  node->ioStreamerPass2() ;
2615  }
2616  RooAbsArg::ioStreamerPass2Finalize() ;
2617 
2618  // Make expensive object cache of all objects point to intermal copy.
2619  // Somehow this doesn't work OK automatically
2620  TIterator* iter = _allOwnedNodes.createIterator() ;
2621  while((node=(RooAbsArg*)iter->Next())) {
2622  node->setExpensiveObjectCache(_eocache) ;
2623  node->setWorkspace(*this);
2624  if (node->IsA()->InheritsFrom(RooAbsOptTestStatistic::Class())) {
2625  RooAbsOptTestStatistic *tmp = (RooAbsOptTestStatistic *)node;
2626  if (tmp->isSealed() && tmp->sealNotice() && strlen(tmp->sealNotice()) > 0) {
2627  cout << "RooWorkspace::Streamer(" << GetName() << ") " << node->IsA()->GetName() << "::" << node->GetName()
2628  << " : " << tmp->sealNotice() << endl;
2629  }
2630  }
2631  }
2632  delete iter ;
2633 
2634 
2635  } else {
2636 
2637  // Make lists of external clients of WS objects, and remove those links temporarily
2638 
2639  map<RooAbsArg*,vector<RooAbsArg *> > extClients, extValueClients, extShapeClients ;
2640 
2641  TIterator* iter = _allOwnedNodes.createIterator() ;
2642  RooAbsArg* tmparg ;
2643  while((tmparg=(RooAbsArg*)iter->Next())) {
2644 
2645  // Loop over client list of this arg
2646  std::vector<RooAbsArg *> clientsTmp{tmparg->_clientList.begin(), tmparg->_clientList.end()};
2647  for (auto client : clientsTmp) {
2648  if (!_allOwnedNodes.containsInstance(*client)) {
2649 
2650  const auto refCount = tmparg->_clientList.refCount(client);
2651  auto& bufferVec = extClients[tmparg];
2652 
2653  bufferVec.insert(bufferVec.end(), refCount, client);
2654  tmparg->_clientList.Remove(client, true);
2655  }
2656  }
2657 
2658  // Loop over value client list of this arg
2659  clientsTmp.assign(tmparg->_clientListValue.begin(), tmparg->_clientListValue.end());
2660  for (auto vclient : clientsTmp) {
2661  if (!_allOwnedNodes.containsInstance(*vclient)) {
2662  cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2663  << " has external value client link to " << vclient << " (" << vclient->GetName() << ") with ref count " << tmparg->_clientListValue.refCount(vclient) << endl ;
2664 
2665  const auto refCount = tmparg->_clientListValue.refCount(vclient);
2666  auto& bufferVec = extValueClients[tmparg];
2667 
2668  bufferVec.insert(bufferVec.end(), refCount, vclient);
2669  tmparg->_clientListValue.Remove(vclient, true);
2670  }
2671  }
2672 
2673  // Loop over shape client list of this arg
2674  clientsTmp.assign(tmparg->_clientListShape.begin(), tmparg->_clientListShape.end());
2675  for (auto sclient : clientsTmp) {
2676  if (!_allOwnedNodes.containsInstance(*sclient)) {
2677  cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2678  << " has external shape client link to " << sclient << " (" << sclient->GetName() << ") with ref count " << tmparg->_clientListShape.refCount(sclient) << endl ;
2679 
2680  const auto refCount = tmparg->_clientListShape.refCount(sclient);
2681  auto& bufferVec = extShapeClients[tmparg];
2682 
2683  bufferVec.insert(bufferVec.end(), refCount, sclient);
2684  tmparg->_clientListShape.Remove(sclient, true);
2685  }
2686  }
2687 
2688  }
2689  delete iter ;
2690 
2691  R__b.WriteClassBuffer(RooWorkspace::Class(),this);
2692 
2693  // Reinstate clients here
2694 
2695 
2696  for (auto iterx : extClients) {
2697  for (auto client : iterx.second) {
2698  iterx.first->_clientList.Add(client);
2699  }
2700  }
2701 
2702  for (auto iterx : extValueClients) {
2703  for (auto client : iterx.second) {
2704  iterx.first->_clientListValue.Add(client);
2705  }
2706  }
2707 
2708  for (auto iterx : extShapeClients) {
2709  for (auto client : iterx.second) {
2710  iterx.first->_clientListShape.Add(client);
2711  }
2712  }
2713 
2714  }
2715 }
2716 
2717 
2718 
2719 
2720 ////////////////////////////////////////////////////////////////////////////////
2721 /// Return STL string with last of class names contained in the code repository
2722 
2723 std::string RooWorkspace::CodeRepo::listOfClassNames() const
2724 {
2725  string ret ;
2726  map<TString,ClassRelInfo>::const_iterator iter = _c2fmap.begin() ;
2727  while(iter!=_c2fmap.end()) {
2728  if (ret.size()>0) {
2729  ret += ", " ;
2730  }
2731  ret += iter->first ;
2732  ++iter ;
2733  }
2734 
2735  return ret ;
2736 }
2737 
2738 namespace {
2739 UInt_t crc32(const char* data, ULong_t sz, UInt_t crc)
2740 {
2741  // update CRC32 with new data
2742 
2743  // use precomputed table, rather than computing it on the fly
2744  static const UInt_t crctab[256] = { 0x00000000,
2745  0x04c11db7, 0x09823b6e, 0x0d4326d9, 0x130476dc, 0x17c56b6b,
2746  0x1a864db2, 0x1e475005, 0x2608edb8, 0x22c9f00f, 0x2f8ad6d6,
2747  0x2b4bcb61, 0x350c9b64, 0x31cd86d3, 0x3c8ea00a, 0x384fbdbd,
2748  0x4c11db70, 0x48d0c6c7, 0x4593e01e, 0x4152fda9, 0x5f15adac,
2749  0x5bd4b01b, 0x569796c2, 0x52568b75, 0x6a1936c8, 0x6ed82b7f,
2750  0x639b0da6, 0x675a1011, 0x791d4014, 0x7ddc5da3, 0x709f7b7a,
2751  0x745e66cd, 0x9823b6e0, 0x9ce2ab57, 0x91a18d8e, 0x95609039,
2752  0x8b27c03c, 0x8fe6dd8b, 0x82a5fb52, 0x8664e6e5, 0xbe2b5b58,
2753  0xbaea46ef, 0xb7a96036, 0xb3687d81, 0xad2f2d84, 0xa9ee3033,
2754  0xa4ad16ea, 0xa06c0b5d, 0xd4326d90, 0xd0f37027, 0xddb056fe,
2755  0xd9714b49, 0xc7361b4c, 0xc3f706fb, 0xceb42022, 0xca753d95,
2756  0xf23a8028, 0xf6fb9d9f, 0xfbb8bb46, 0xff79a6f1, 0xe13ef6f4,
2757  0xe5ffeb43, 0xe8bccd9a, 0xec7dd02d, 0x34867077, 0x30476dc0,
2758  0x3d044b19, 0x39c556ae, 0x278206ab, 0x23431b1c, 0x2e003dc5,
2759  0x2ac12072, 0x128e9dcf, 0x164f8078, 0x1b0ca6a1, 0x1fcdbb16,
2760  0x018aeb13, 0x054bf6a4, 0x0808d07d, 0x0cc9cdca, 0x7897ab07,
2761  0x7c56b6b0, 0x71159069, 0x75d48dde, 0x6b93dddb, 0x6f52c06c,
2762  0x6211e6b5, 0x66d0fb02, 0x5e9f46bf, 0x5a5e5b08, 0x571d7dd1,
2763  0x53dc6066, 0x4d9b3063, 0x495a2dd4, 0x44190b0d, 0x40d816ba,
2764  0xaca5c697, 0xa864db20, 0xa527fdf9, 0xa1e6e04e, 0xbfa1b04b,
2765  0xbb60adfc, 0xb6238b25, 0xb2e29692, 0x8aad2b2f, 0x8e6c3698,
2766  0x832f1041, 0x87ee0df6, 0x99a95df3, 0x9d684044, 0x902b669d,
2767  0x94ea7b2a, 0xe0b41de7, 0xe4750050, 0xe9362689, 0xedf73b3e,
2768  0xf3b06b3b, 0xf771768c, 0xfa325055, 0xfef34de2, 0xc6bcf05f,
2769  0xc27dede8, 0xcf3ecb31, 0xcbffd686, 0xd5b88683, 0xd1799b34,
2770  0xdc3abded, 0xd8fba05a, 0x690ce0ee, 0x6dcdfd59, 0x608edb80,
2771  0x644fc637, 0x7a089632, 0x7ec98b85, 0x738aad5c, 0x774bb0eb,
2772  0x4f040d56, 0x4bc510e1, 0x46863638, 0x42472b8f, 0x5c007b8a,
2773  0x58c1663d, 0x558240e4, 0x51435d53, 0x251d3b9e, 0x21dc2629,
2774  0x2c9f00f0, 0x285e1d47, 0x36194d42, 0x32d850f5, 0x3f9b762c,
2775  0x3b5a6b9b, 0x0315d626, 0x07d4cb91, 0x0a97ed48, 0x0e56f0ff,
2776  0x1011a0fa, 0x14d0bd4d, 0x19939b94, 0x1d528623, 0xf12f560e,
2777  0xf5ee4bb9, 0xf8ad6d60, 0xfc6c70d7, 0xe22b20d2, 0xe6ea3d65,
2778  0xeba91bbc, 0xef68060b, 0xd727bbb6, 0xd3e6a601, 0xdea580d8,
2779  0xda649d6f, 0xc423cd6a, 0xc0e2d0dd, 0xcda1f604, 0xc960ebb3,
2780  0xbd3e8d7e, 0xb9ff90c9, 0xb4bcb610, 0xb07daba7, 0xae3afba2,
2781  0xaafbe615, 0xa7b8c0cc, 0xa379dd7b, 0x9b3660c6, 0x9ff77d71,
2782  0x92b45ba8, 0x9675461f, 0x8832161a, 0x8cf30bad, 0x81b02d74,
2783  0x857130c3, 0x5d8a9099, 0x594b8d2e, 0x5408abf7, 0x50c9b640,
2784  0x4e8ee645, 0x4a4ffbf2, 0x470cdd2b, 0x43cdc09c, 0x7b827d21,
2785  0x7f436096, 0x7200464f, 0x76c15bf8, 0x68860bfd, 0x6c47164a,
2786  0x61043093, 0x65c52d24, 0x119b4be9, 0x155a565e, 0x18197087,
2787  0x1cd86d30, 0x029f3d35, 0x065e2082, 0x0b1d065b, 0x0fdc1bec,
2788  0x3793a651, 0x3352bbe6, 0x3e119d3f, 0x3ad08088, 0x2497d08d,
2789  0x2056cd3a, 0x2d15ebe3, 0x29d4f654, 0xc5a92679, 0xc1683bce,
2790  0xcc2b1d17, 0xc8ea00a0, 0xd6ad50a5, 0xd26c4d12, 0xdf2f6bcb,
2791  0xdbee767c, 0xe3a1cbc1, 0xe760d676, 0xea23f0af, 0xeee2ed18,
2792  0xf0a5bd1d, 0xf464a0aa, 0xf9278673, 0xfde69bc4, 0x89b8fd09,
2793  0x8d79e0be, 0x803ac667, 0x84fbdbd0, 0x9abc8bd5, 0x9e7d9662,
2794  0x933eb0bb, 0x97ffad0c, 0xafb010b1, 0xab710d06, 0xa6322bdf,
2795  0xa2f33668, 0xbcb4666d, 0xb8757bda, 0xb5365d03, 0xb1f740b4
2796  };
2797 
2798  crc = ~crc;
2799  while (sz--) crc = (crc << 8) ^ UInt_t(*data++) ^ crctab[crc >> 24];
2800 
2801  return ~crc;
2802 }
2803 
2804 UInt_t crc32(const char* data)
2805 {
2806  // Calculate crc32 checksum on given string
2807  unsigned long sz = strlen(data);
2808  switch (strlen(data)) {
2809  case 0:
2810  return 0;
2811  case 1:
2812  return data[0];
2813  case 2:
2814  return (data[0] << 8) | data[1];
2815  case 3:
2816  return (data[0] << 16) | (data[1] << 8) | data[2];
2817  case 4:
2818  return (data[0] << 24) | (data[1] << 16) | (data[2] << 8) | data[3];
2819  default:
2820  return crc32(data + 4, sz - 4, (data[0] << 24) | (data[1] << 16) |
2821  (data[2] << 8) | data[3]);
2822  }
2823 }
2824 
2825 }
2826 
2827 ////////////////////////////////////////////////////////////////////////////////
2828 /// For all classes in the workspace for which no class definition is
2829 /// found in the ROOT class table extract source code stored in code
2830 /// repository into temporary directory set by
2831 /// setClassFileExportDir(), compile classes and link them with
2832 /// current ROOT session. If a compilation error occurs print
2833 /// instructions for user how to fix errors and recover workspace and
2834 /// abort import procedure.
2835 
2836 Bool_t RooWorkspace::CodeRepo::compileClasses()
2837 {
2838  Bool_t haveDir=kFALSE ;
2839 
2840  // Retrieve name of directory in which to export code files
2841  string dirName = Form(_classFileExportDir.c_str(),_wspace->uuid().AsString(),_wspace->GetName()) ;
2842 
2843  Bool_t writeExtraHeaders(kFALSE) ;
2844 
2845  // Process all class entries in repository
2846  map<TString,ClassRelInfo>::iterator iter = _c2fmap.begin() ;
2847  while(iter!=_c2fmap.end()) {
2848 
2849  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing class " << iter->first.Data() << endl ;
2850 
2851  // If class is already known, don't load
2852  if (gClassTable->GetDict(iter->first.Data())) {
2853  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Embedded class "
2854  << iter->first << " already in ROOT class table, skipping" << endl ;
2855  ++iter ;
2856  continue ;
2857  }
2858 
2859  // Check that export directory exists
2860  if (!haveDir) {
2861 
2862  // If not, make local directory to extract files
2863  if (!gSystem->AccessPathName(dirName.c_str())) {
2864  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() reusing code export directory " << dirName.c_str()
2865  << " to extract coded embedded in workspace" << endl ;
2866  } else {
2867  if (gSystem->MakeDirectory(dirName.c_str())==0) {
2868  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() creating code export directory " << dirName.c_str()
2869  << " to extract coded embedded in workspace" << endl ;
2870  } else {
2871  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR creating code export directory " << dirName.c_str()
2872  << " to extract coded embedded in workspace" << endl ;
2873  return kFALSE ;
2874  }
2875  }
2876  haveDir=kTRUE ;
2877 
2878  }
2879 
2880  // First write any extra header files
2881  if (!writeExtraHeaders) {
2882  writeExtraHeaders = kTRUE ;
2883 
2884  map<TString,ExtraHeader>::iterator eiter = _ehmap.begin() ;
2885  while(eiter!=_ehmap.end()) {
2886 
2887  // Check if identical declaration file (header) is already written
2888  Bool_t needEHWrite=kTRUE ;
2889  string fdname = Form("%s/%s",dirName.c_str(),eiter->second._hname.Data()) ;
2890  ifstream ifdecl(fdname.c_str()) ;
2891  if (ifdecl) {
2892  TString contents ;
2893  char buf[64000];
2894  while (ifdecl.getline(buf, 64000)) {
2895  contents += buf;
2896  contents += "\n";
2897  }
2898  UInt_t crcFile = crc32(contents.Data());
2899  UInt_t crcWS = crc32(eiter->second._hfile.Data());
2900  needEHWrite = (crcFile != crcWS);
2901  }
2902 
2903  // Write declaration file if required
2904  if (needEHWrite) {
2905  oocoutI(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting extra header file "
2906  << fdname << endl;
2907 
2908  // Extra headers may contain non-existing path - create first to be sure
2909  gSystem->MakeDirectory(gSystem->DirName(fdname.c_str()));
2910 
2911  ofstream fdecl(fdname.c_str());
2912  if (!fdecl) {
2913  oocoutE(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file " << fdname
2914  << " for writing" << endl;
2915  return kFALSE;
2916  }
2917  fdecl << eiter->second._hfile.Data();
2918  fdecl.close();
2919  }
2920  ++eiter;
2921  }
2922  }
2923 
2924 
2925  // Navigate from class to file
2926  ClassFiles& cfinfo = _fmap[iter->second._fileBase] ;
2927 
2928  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing file with base " << iter->second._fileBase << endl ;
2929 
2930  // If file is already processed, skip to next class
2931  if (cfinfo._extracted) {
2932  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() file with base name " << iter->second._fileBase
2933  << " has already been extracted, skipping to next class" << endl ;
2934  continue ;
2935  }
2936 
2937  // Check if identical declaration file (header) is already written
2938  Bool_t needDeclWrite=kTRUE ;
2939  string fdname = Form("%s/%s.%s",dirName.c_str(),iter->second._fileBase.Data(),cfinfo._hext.Data()) ;
2940  ifstream ifdecl(fdname.c_str()) ;
2941  if (ifdecl) {
2942  TString contents ;
2943  char buf[64000];
2944  while (ifdecl.getline(buf, 64000)) {
2945  contents += buf;
2946  contents += "\n";
2947  }
2948  UInt_t crcFile = crc32(contents.Data()) ;
2949  UInt_t crcWS = crc32(cfinfo._hfile.Data()) ;
2950  needDeclWrite = (crcFile!=crcWS) ;
2951  }
2952 
2953  // Write declaration file if required
2954  if (needDeclWrite) {
2955  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting declaration code of class " << iter->first << ", file " << fdname << endl ;
2956  ofstream fdecl(fdname.c_str()) ;
2957  if (!fdecl) {
2958  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file "
2959  << fdname << " for writing" << endl ;
2960  return kFALSE ;
2961  }
2962  fdecl << cfinfo._hfile ;
2963  fdecl.close() ;
2964  }
2965 
2966  // Check if identical implementation file is already written
2967  Bool_t needImplWrite=kTRUE ;
2968  string finame = Form("%s/%s.cxx",dirName.c_str(),iter->second._fileBase.Data()) ;
2969  ifstream ifimpl(finame.c_str()) ;
2970  if (ifimpl) {
2971  TString contents ;
2972  char buf[64000];
2973  while (ifimpl.getline(buf, 64000)) {
2974  contents += buf;
2975  contents += "\n";
2976  }
2977  UInt_t crcFile = crc32(contents.Data()) ;
2978  UInt_t crcWS = crc32(cfinfo._cxxfile.Data()) ;
2979  needImplWrite = (crcFile!=crcWS) ;
2980  }
2981 
2982  // Write implementation file if required
2983  if (needImplWrite) {
2984  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting implementation code of class " << iter->first << ", file " << finame << endl ;
2985  ofstream fimpl(finame.c_str()) ;
2986  if (!fimpl) {
2987  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file"
2988  << finame << " for writing" << endl ;
2989  return kFALSE ;
2990  }
2991  fimpl << cfinfo._cxxfile ;
2992  fimpl.close() ;
2993  }
2994 
2995  // Mark this file as extracted
2996  cfinfo._extracted = kTRUE ;
2997  oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() marking code unit " << iter->second._fileBase << " as extracted" << endl ;
2998 
2999  // Compile class
3000  oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Compiling code unit " << iter->second._fileBase.Data() << " to define class " << iter->first << endl ;
3001  Bool_t ok = gSystem->CompileMacro(finame.c_str(),"k") ;
3002 
3003  if (!ok) {
3004  oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR compiling class " << iter->first.Data() << ", to fix this you can do the following: " << endl
3005  << " 1) Fix extracted source code files in directory " << dirName.c_str() << "/" << endl
3006  << " 2) In clean ROOT session compiled fixed classes by hand using '.x " << dirName.c_str() << "/ClassName.cxx+'" << endl
3007  << " 3) Reopen file with RooWorkspace with broken source code in UPDATE mode. Access RooWorkspace to force loading of class" << endl
3008  << " Broken instances in workspace will _not_ be compiled, instead precompiled fixed instances will be used." << endl
3009  << " 4) Reimport fixed code in workspace using 'RooWorkspace::importClassCode(\"*\",kTRUE)' method, Write() updated workspace to file and close file" << endl
3010  << " 5) Reopen file in clean ROOT session to confirm that problems are fixed" << endl ;
3011  return kFALSE ;
3012  }
3013 
3014  ++iter ;
3015  }
3016 
3017  return kTRUE ;
3018 }
3019 
3020 
3021 
3022 ////////////////////////////////////////////////////////////////////////////////
3023 /// Internal access to TDirectory append method
3024 
3025 void RooWorkspace::WSDir::InternalAppend(TObject* obj)
3026 {
3027  TDirectory::Append(obj,kFALSE) ;
3028 }
3029 
3030 
3031 ////////////////////////////////////////////////////////////////////////////////
3032 /// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
3033 
3034 void RooWorkspace::WSDir::Add(TObject* obj,Bool_t)
3035 {
3036  if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
3037  coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << endl ;
3038  } else {
3039  InternalAppend(obj) ;
3040  }
3041 }
3042 
3043 
3044 ////////////////////////////////////////////////////////////////////////////////
3045 /// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
3046 
3047 void RooWorkspace::WSDir::Append(TObject* obj,Bool_t)
3048 {
3049  if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
3050  coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << endl ;
3051  } else {
3052  InternalAppend(obj) ;
3053  }
3054 }
3055 
3056 
3057 
3058 ////////////////////////////////////////////////////////////////////////////////
3059 /// Activate export of workspace symbols to CINT in a namespace with given name. If no name
3060 /// is given the namespace will have the same name as the workspace
3061 
3062 void RooWorkspace::exportToCint(const char* nsname)
3063 {
3064  // If export is already active, do nothing
3065  if (_doExport) {
3066  coutE(ObjectHandling) << "RooWorkspace::exportToCint(" << GetName() << ") WARNING: repeated calls to exportToCint() have no effect" << endl ;
3067  return ;
3068  }
3069 
3070  // Set flag so that future import to workspace are automatically exported to CINT
3071  _doExport = kTRUE ;
3072 
3073  // If no name is provided choose name of workspace
3074  if (!nsname) nsname = GetName() ;
3075  _exportNSName = nsname ;
3076 
3077  coutI(ObjectHandling) << "RooWorkspace::exportToCint(" << GetName()
3078  << ") INFO: references to all objects in this workspace will be created in CINT in 'namespace " << _exportNSName << "'" << endl ;
3079 
3080  // Export present contents of workspace to CINT
3081  TIterator* iter = _allOwnedNodes.createIterator() ;
3082  TObject* wobj ;
3083  while((wobj=iter->Next())) {
3084  exportObj(wobj) ;
3085  }
3086  delete iter ;
3087  iter = _dataList.MakeIterator() ;
3088  while((wobj=iter->Next())) {
3089  exportObj(wobj) ;
3090  }
3091  delete iter ;
3092 }
3093 
3094 
3095 ////////////////////////////////////////////////////////////////////////////////
3096 /// Export reference to given workspace object to CINT
3097 
3098 void RooWorkspace::exportObj(TObject* wobj)
3099 {
3100  // Do nothing if export flag is not set
3101  if (!_doExport) return ;
3102 
3103  // Do not export RooConstVars
3104  if (wobj->IsA() == RooConstVar::Class()) {
3105  return ;
3106  }
3107 
3108 
3109  // Determine if object name is a valid C++ identifier name
3110 
3111  // Do not export objects that have names that are not valid C++ identifiers
3112  if (!isValidCPPID(wobj->GetName())) {
3113  cxcoutD(ObjectHandling) << "RooWorkspace::exportObj(" << GetName() << ") INFO: Workspace object name " << wobj->GetName() << " is not a valid C++ identifier and is not exported to CINT" << endl ;
3114  return ;
3115  }
3116 
3117  // Declare correctly typed reference to object in CINT in the namespace associated with this workspace
3118  string cintExpr = Form("namespace %s { %s& %s = *(%s *)0x%lx ; }",_exportNSName.c_str(),wobj->IsA()->GetName(),wobj->GetName(),wobj->IsA()->GetName(),(ULong_t)wobj) ;
3119  gROOT->ProcessLine(cintExpr.c_str()) ;
3120 }
3121 
3122 
3123 
3124 ////////////////////////////////////////////////////////////////////////////////
3125 /// Return true if given name is a valid C++ identifier name
3126 
3127 Bool_t RooWorkspace::isValidCPPID(const char* name)
3128 {
3129  string oname(name) ;
3130  if (isdigit(oname[0])) {
3131  return kFALSE ;
3132  } else {
3133  for (UInt_t i=0 ; i<oname.size() ; i++) {
3134  char c = oname[i] ;
3135  if (!isalnum(c) && (c!='_')) {
3136  return kFALSE ;
3137  }
3138  }
3139  }
3140  return kTRUE ;
3141 }
3142 
3143 ////////////////////////////////////////////////////////////////////////////////
3144 /// Delete exported reference in CINT namespace
3145 
3146 void RooWorkspace::unExport()
3147 {
3148  char buf[64000];
3149  TIterator *iter = _allOwnedNodes.createIterator();
3150  TObject *wobj;
3151  while ((wobj = iter->Next())) {
3152  if (isValidCPPID(wobj->GetName())) {
3153  strlcpy(buf, Form("%s::%s", _exportNSName.c_str(), wobj->GetName()), 64000);
3154  gInterpreter->DeleteVariable(buf);
3155  }
3156  }
3157  delete iter ;
3158 }
3159 
3160 ////////////////////////////////////////////////////////////////////////////////
3161 /// If one of the TObject we have a referenced to is deleted, remove the
3162 /// reference.
3163 
3164 void RooWorkspace::RecursiveRemove(TObject *removedObj)
3165 {
3166  _dataList.RecursiveRemove(removedObj);
3167  if (removedObj == _dir) _dir = nullptr;
3168 
3169  _allOwnedNodes.RecursiveRemove(removedObj); // RooArgSet
3170 
3171  _dataList.RecursiveRemove(removedObj);
3172  _embeddedDataList.RecursiveRemove(removedObj);
3173  _views.RecursiveRemove(removedObj);
3174  _snapshots.RecursiveRemove(removedObj);
3175  _genObjects.RecursiveRemove(removedObj);
3176  _studyMods.RecursiveRemove(removedObj);
3177 
3178  for(auto c : _namedSets) {
3179  c.second.RecursiveRemove(removedObj);
3180  }
3181 
3182  _eocache.RecursiveRemove(removedObj); // RooExpensiveObjectCache
3183 }