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