Logo ROOT   6.16/01
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
755 const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
756 const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
757 const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
758
759{
760
761 coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") importing dataset " << inData.GetName() << endl ;
762
763 RooLinkedList args ;
764 args.Add((TObject*)&arg1) ;
765 args.Add((TObject*)&arg2) ;
766 args.Add((TObject*)&arg3) ;
767 args.Add((TObject*)&arg4) ;
768 args.Add((TObject*)&arg5) ;
769 args.Add((TObject*)&arg6) ;
770 args.Add((TObject*)&arg7) ;
771 args.Add((TObject*)&arg8) ;
772 args.Add((TObject*)&arg9) ;
773
774 // Select the pdf-specific commands
775 RooCmdConfig pc(Form("RooWorkspace::import(%s)",GetName())) ;
776
777 pc.defineString("dsetName","Rename",0,"") ;
778 pc.defineString("varChangeIn","RenameVar",0,"",kTRUE) ;
779 pc.defineString("varChangeOut","RenameVar",1,"",kTRUE) ;
780 pc.defineInt("embedded","Embedded",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
794 // Transform emtpy string into null pointer
795 if (dsetName && strlen(dsetName)==0) {
796 dsetName=0 ;
797 }
798
799 RooLinkedList& dataList = embedded ? _embeddedDataList : _dataList ;
800
801 // Check that no dataset with target name already exists
802 if (dsetName && dataList.FindObject(dsetName)) {
803 coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << dsetName << " already exists in workspace, import aborted" << endl ;
804 return kTRUE ;
805 }
806 if (!dsetName && dataList.FindObject(inData.GetName())) {
807 coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << inData.GetName() << " already exists in workspace, import aborted" << endl ;
808 return kTRUE ;
809 }
810
811 // Rename dataset if required
812 RooAbsData* clone ;
813 if (dsetName) {
814 coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset from " << inData.GetName() << " to " << dsetName << endl ;
815 clone = (RooAbsData*) inData.Clone(dsetName) ;
816 } else {
817 clone = (RooAbsData*) inData.Clone(inData.GetName()) ;
818 }
819
820
821 // Process any change in variable names
822 if (strlen(varChangeIn)>0) {
823 // Parse comma separated lists of variable name changes
824 const std::vector<std::string> tokIn = RooHelpers::tokenise(varChangeIn, ",");
825 const std::vector<std::string> tokOut = RooHelpers::tokenise(varChangeOut, ",");
826 for (unsigned int i=0; i < tokIn.size(); ++i) {
827 coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset observable " << tokIn[i] << " to " << tokOut[i] << endl ;
828 clone->changeObservableName(tokIn[i].c_str(), tokOut[i].c_str());
829 }
830 }
831
832 // Now import the dataset observables, unless dataset is embedded
833 RooAbsArg* carg ;
834 if (!embedded) {
835 TIterator* iter = clone->get()->createIterator() ;
836 while((carg=(RooAbsArg*)iter->Next())) {
837 if (!arg(carg->GetName())) {
838 import(*carg) ;
839 }
840 }
841 delete iter ;
842 }
843
844 dataList.Add(clone) ;
845 if (_dir) {
846 _dir->InternalAppend(clone) ;
847 }
848 if (_doExport) {
849 exportObj(clone) ;
850 }
851
852 // Set expensive object cache of dataset internal buffers to that of workspace
853 RooFIter iter2 = clone->get()->fwdIterator() ;
854 while ((carg=iter2.next())) {
856 }
857
858
859 return kFALSE ;
860}
861
862
863
864
865////////////////////////////////////////////////////////////////////////////////
866/// Define a named RooArgSet with given constituents. If importMissing is true, any constituents
867/// of aset that are not in the workspace will be imported, otherwise an error is returned
868/// for missing components
869
870Bool_t RooWorkspace::defineSet(const char* name, const RooArgSet& aset, Bool_t importMissing)
871{
872 // Check if set was previously defined, if so print warning
873 map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
874 if (i!=_namedSets.end()) {
875 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << endl ;
876 }
877
878 RooArgSet wsargs ;
879
880 // Check all constituents of provided set
881 TIterator* iter = aset.createIterator() ;
882 RooAbsArg* sarg ;
883 while((sarg=(RooAbsArg*)iter->Next())) {
884 // If missing, either import or report error
885 if (!arg(sarg->GetName())) {
886 if (importMissing) {
887 import(*sarg) ;
888 } else {
889 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR set constituent \"" << sarg->GetName()
890 << "\" is not in workspace and importMissing option is disabled" << endl ;
891 return kTRUE ;
892 }
893 }
894 wsargs.add(*arg(sarg->GetName())) ;
895 }
896 delete iter ;
897
898 // Install named set
899 _namedSets[name].removeAll() ;
900 _namedSets[name].add(wsargs) ;
901
902 return kFALSE ;
903}
904
905//_____________________________________________________________________________
907{
908 // Define a named RooArgSet with given constituents. If importMissing is true, any constituents
909 // of aset that are not in the workspace will be imported, otherwise an error is returned
910 // for missing components
911
912 // Check if set was previously defined, if so print warning
913 map<string, RooArgSet>::iterator i = _namedSets.find(name);
914 if (i != _namedSets.end()) {
915 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName()
916 << ") WARNING redefining previously defined named set " << name << endl;
917 }
918
919 // Install named set
920 _namedSets[name].removeAll();
921 _namedSets[name].add(aset);
922
923 return kFALSE;
924}
925
926////////////////////////////////////////////////////////////////////////////////
927/// Define a named set in the work space through a comma separated list of
928/// names of objects already in the workspace
929
930Bool_t RooWorkspace::defineSet(const char* name, const char* contentList)
931{
932 // Check if set was previously defined, if so print warning
933 map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
934 if (i!=_namedSets.end()) {
935 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << endl ;
936 }
937
938 RooArgSet wsargs ;
939
940 // Check all constituents of provided set
941 for (const std::string& token : RooHelpers::tokenise(contentList, ",")) {
942 // If missing, either import or report error
943 if (!arg(token.c_str())) {
944 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
945 << "\" is not in workspace" << endl ;
946 return kTRUE ;
947 }
948 wsargs.add(*arg(token.c_str())) ;
949 }
950
951 // Install named set
952 _namedSets[name].removeAll() ;
953 _namedSets[name].add(wsargs) ;
954
955 return kFALSE ;
956}
957
958
959
960
961////////////////////////////////////////////////////////////////////////////////
962/// Define a named set in the work space through a comma separated list of
963/// names of objects already in the workspace
964
965Bool_t RooWorkspace::extendSet(const char* name, const char* newContents)
966{
967 RooArgSet wsargs ;
968
969 // Check all constituents of provided set
970 for (const std::string& token : RooHelpers::tokenise(newContents, ",")) {
971 // If missing, either import or report error
972 if (!arg(token.c_str())) {
973 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
974 << "\" is not in workspace" << endl ;
975 return kTRUE ;
976 }
977 wsargs.add(*arg(token.c_str())) ;
978 }
979
980 // Extend named set
981 _namedSets[name].add(wsargs,kTRUE) ;
982
983 return kFALSE ;
984}
985
986
987
988////////////////////////////////////////////////////////////////////////////////
989/// Return pointer to previously defined named set with given nmame
990/// If no such set is found a null pointer is returned
991
992const RooArgSet* RooWorkspace::set(const char* name)
993{
994 map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
995 return (i!=_namedSets.end()) ? &(i->second) : 0 ;
996}
997
998
999
1000
1001////////////////////////////////////////////////////////////////////////////////
1002/// Rename set to a new name
1003
1004Bool_t RooWorkspace::renameSet(const char* name, const char* newName)
1005{
1006 // First check if set exists
1007 if (!set(name)) {
1008 coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << name
1009 << " does not exist" << endl ;
1010 return kTRUE ;
1011 }
1012
1013 // Check if no set exists with new name
1014 if (set(newName)) {
1015 coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << newName
1016 << " already exists" << endl ;
1017 return kTRUE ;
1018 }
1019
1020 // Copy entry under 'name' to 'newName'
1021 _namedSets[newName].add(_namedSets[name]) ;
1022
1023 // Remove entry under old name
1024 _namedSets.erase(name) ;
1025
1026 return kFALSE ;
1027}
1028
1029
1030
1031
1032////////////////////////////////////////////////////////////////////////////////
1033/// Remove a named set from the workspace
1034
1036{
1037 // First check if set exists
1038 if (!set(name)) {
1039 coutE(InputArguments) << "RooWorkspace::removeSet(" << GetName() << ") ERROR a set with name " << name
1040 << " does not exist" << endl ;
1041 return kTRUE ;
1042 }
1043
1044 // Remove set with given name
1045 _namedSets.erase(name) ;
1046
1047 return kFALSE ;
1048}
1049
1050
1051
1052
1053////////////////////////////////////////////////////////////////////////////////
1054/// Open an import transaction operations. Returns kTRUE if successful, kFALSE
1055/// if there is already an ongoing transaction
1056
1058{
1059 // Check that there was no ongoing transaction
1060 if (_openTrans) {
1061 return kFALSE ;
1062 }
1063
1064 // Open transaction
1065 _openTrans = kTRUE ;
1066 return kTRUE ;
1067}
1068
1069
1070
1071
1072////////////////////////////////////////////////////////////////////////////////
1073/// Cancel an ongoing import transaction. All objects imported since startTransaction()
1074/// will be removed and the transaction will be terminated. Return kTRUE if cancel operation
1075/// succeeds, return kFALSE if there was no open transaction
1076
1078{
1079 // Check that there is an ongoing transaction
1080 if (!_openTrans) {
1081 return kFALSE ;
1082 }
1083
1084 // Delete all objects in the sandbox
1086 RooAbsArg* tmpArg ;
1087 while((tmpArg=(RooAbsArg*)iter->Next())) {
1088 _allOwnedNodes.remove(*tmpArg) ;
1089 }
1090 delete iter ;
1092
1093 // Mark transaction as finished
1094 _openTrans = kFALSE ;
1095
1096 return kTRUE ;
1097}
1098
1100{
1101 // Commit an ongoing import transaction. Returns kTRUE if commit succeeded,
1102 // return kFALSE if there was no ongoing transaction
1103
1104 // Check that there is an ongoing transaction
1105 if (!_openTrans) {
1106 return kFALSE ;
1107 }
1108
1109 // Publish sandbox nodes in directory and/or CINT if requested
1111 RooAbsArg* sarg ;
1112 while((sarg=(RooAbsArg*)iter->Next())) {
1113 if (_dir && sarg->IsA() != RooConstVar::Class()) {
1114 _dir->InternalAppend(sarg) ;
1115 }
1116 if (_doExport && sarg->IsA() != RooConstVar::Class()) {
1117 exportObj(sarg) ;
1118 }
1119 }
1120 delete iter ;
1121
1122 // Remove all committed objects from the sandbox
1124
1125 // Mark transaction as finished
1126 _openTrans = kFALSE ;
1127
1128 return kTRUE ;
1129}
1130
1131
1132
1133
1134////////////////////////////////////////////////////////////////////////////////
1135
1137{
1138 return _classes.autoImportClass(theClass,doReplace) ;
1139}
1140
1141
1142
1143////////////////////////////////////////////////////////////////////////////////
1144/// Inport code of all classes in the workspace that have a class name
1145/// that matches pattern 'pat' and which are not found to be part of
1146/// the standard ROOT distribution. If doReplace is true any existing
1147/// class code saved in the workspace is replaced
1148
1149Bool_t RooWorkspace::importClassCode(const char* pat, Bool_t doReplace)
1150{
1151 Bool_t ret(kTRUE) ;
1152
1153 TRegexp re(pat,kTRUE) ;
1154 TIterator* iter = componentIterator() ;
1155 RooAbsArg* carg ;
1156 while((carg=(RooAbsArg*)iter->Next())) {
1157 TString className = carg->IsA()->GetName() ;
1158 if (className.Index(re)>=0 && !_classes.autoImportClass(carg->IsA(),doReplace)) {
1159 coutW(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") WARNING: problems import class code of object "
1160 << carg->IsA()->GetName() << "::" << carg->GetName() << ", reading of workspace will require external definition of class" << endl ;
1161 ret = kFALSE ;
1162 }
1163 }
1164 delete iter ;
1165
1166 return ret ;
1167}
1168
1169
1170
1171
1172
1173////////////////////////////////////////////////////////////////////////////////
1174/// Save snapshot of values and attributes (including "Constant") of parameters 'params'
1175/// If importValues is FALSE, the present values from the object in the workspace are
1176/// saved. If importValues is TRUE, the values of the objects passed in the 'params'
1177/// argument are saved
1178
1179Bool_t RooWorkspace::saveSnapshot(const char* name, const char* paramNames)
1180{
1181 return saveSnapshot(name,argSet(paramNames),kFALSE) ;
1182}
1183
1184
1185
1186
1187
1188////////////////////////////////////////////////////////////////////////////////
1189/// Save snapshot of values and attributes (including "Constant") of parameters 'params'
1190/// If importValues is FALSE, the present values from the object in the workspace are
1191/// saved. If importValues is TRUE, the values of the objects passed in the 'params'
1192/// argument are saved
1193
1194Bool_t RooWorkspace::saveSnapshot(const char* name, const RooArgSet& params, Bool_t importValues)
1195{
1196 RooArgSet* actualParams = (RooArgSet*) _allOwnedNodes.selectCommon(params) ;
1197 RooArgSet* snapshot = (RooArgSet*) actualParams->snapshot() ;
1198 delete actualParams ;
1199
1200 snapshot->setName(name) ;
1201
1202 if (importValues) {
1203 *snapshot = params ;
1204 }
1205
1207 if (oldSnap) {
1208 coutI(ObjectHandling) << "RooWorkspace::saveSnaphot(" << GetName() << ") replacing previous snapshot with name " << name << endl ;
1209 _snapshots.Remove(oldSnap) ;
1210 delete oldSnap ;
1211 }
1212
1213 _snapshots.Add(snapshot) ;
1214
1215 return kTRUE ;
1216}
1217
1218
1219
1220
1221////////////////////////////////////////////////////////////////////////////////
1222/// Load the values and attributes of the parameters in the snapshot saved with
1223/// the given name
1224
1226{
1227 RooArgSet* snap = (RooArgSet*) _snapshots.find(name) ;
1228 if (!snap) {
1229 coutE(ObjectHandling) << "RooWorkspace::loadSnapshot(" << GetName() << ") no snapshot with name " << name << " is available" << endl ;
1230 return kFALSE ;
1231 }
1232
1233 RooArgSet* actualParams = (RooArgSet*) _allOwnedNodes.selectCommon(*snap) ;
1234 *actualParams = *snap ;
1235 delete actualParams ;
1236
1237 return kTRUE ;
1238}
1239
1240
1241////////////////////////////////////////////////////////////////////////////////
1242/// Return the RooArgSet containgin a snapshot of variables contained in the workspace
1243///
1244/// Note that the variables of the objects in the snapshots are _copies_ of the
1245/// variables in the workspace. To load the values of a snapshot in the workspace
1246/// variables use loadSnapshot() instead
1247
1249{
1250 RooArgSet* snap = (RooArgSet*) _snapshots.find(name) ;
1251 if (!snap) {
1252 coutE(ObjectHandling) << "RooWorkspace::loadSnapshot(" << GetName() << ") no snapshot with name " << name << " is available" << endl ;
1253 return 0 ;
1254 }
1255
1256 return snap ;
1257}
1258
1259
1260
1261// //_____________________________________________________________________________
1262// RooAbsPdf* RooWorkspace::joinPdf(const char* jointPdfName, const char* indexName, const char* inputMapping)
1263// {
1264// // Join given list of p.d.f.s into a simultaneous p.d.f with given name. If the named index category
1265// // does not exist, it is created.
1266// //
1267// // Example : joinPdf("simPdf","expIndex","A=pdfA,B=pdfB") ;
1268// //
1269// // will return a RooSimultaneous named 'simPdf' with index category 'expIndex' with
1270// // state names A and B. Pdf 'pdfA' will be associated with state A and pdf 'pdfB'
1271// // will be associated with state B
1272// //
1273// return 0 ;
1274// }
1275
1276// //_____________________________________________________________________________
1277// RooAbsData* RooWorkspace::joinData(const char* jointDataName, const char* indexName, const char* inputMapping)
1278// {
1279// // Join given list of dataset into a joint dataset for use with a simultaneous pdf
1280// // (as e.g. created by joingPdf"
1281// //
1282// // Example : joinData("simData","expIndex","A=dataA,B=dataB") ;
1283// //
1284// // will return a RooDataSet named 'simData' that consist of the entries of both
1285// // dataA and dataB. An extra category column 'expIndex' is added that labels
1286// // each entry with state 'A' and 'B' to indicate the originating dataset
1287// return 0 ;
1288// }
1289
1290
1291////////////////////////////////////////////////////////////////////////////////
1292/// Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found
1293
1295{
1296 return dynamic_cast<RooAbsPdf*>(_allOwnedNodes.find(name)) ;
1297}
1298
1299
1300////////////////////////////////////////////////////////////////////////////////
1301/// Retrieve function (RooAbsReal) with given name. Note that all RooAbsPdfs are also RooAbsReals. A null pointer is returned if not found.
1302
1304{
1305 return dynamic_cast<RooAbsReal*>(_allOwnedNodes.find(name)) ;
1306}
1307
1308
1309////////////////////////////////////////////////////////////////////////////////
1310/// Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found
1311
1313{
1314 return dynamic_cast<RooRealVar*>(_allOwnedNodes.find(name)) ;
1315}
1316
1317
1318////////////////////////////////////////////////////////////////////////////////
1319/// Retrieve discrete variable (RooCategory) with given name. A null pointer is returned if not found
1320
1322{
1323 return dynamic_cast<RooCategory*>(_allOwnedNodes.find(name)) ;
1324}
1325
1326
1327////////////////////////////////////////////////////////////////////////////////
1328/// Retrieve discrete function (RooAbsCategory) with given name. A null pointer is returned if not found
1329
1331{
1332 return dynamic_cast<RooAbsCategory*>(_allOwnedNodes.find(name)) ;
1333}
1334
1335
1336
1337////////////////////////////////////////////////////////////////////////////////
1338/// Return RooAbsArg with given name. A null pointer is returned if none is found.
1339
1341{
1342 return _allOwnedNodes.find(name) ;
1343}
1344
1345
1346
1347////////////////////////////////////////////////////////////////////////////////
1348/// Return set of RooAbsArgs matching to given list of names
1349
1350RooArgSet RooWorkspace::argSet(const char* nameList) const
1351{
1352 RooArgSet ret ;
1353
1354 for (const std::string& token : RooHelpers::tokenise(nameList, ",")) {
1355 RooAbsArg* oneArg = arg(token.c_str()) ;
1356 if (oneArg) {
1357 ret.add(*oneArg) ;
1358 } else {
1359 coutE(InputArguments) << " RooWorkspace::argSet(" << GetName() << ") no RooAbsArg named \"" << token << "\" in workspace" << endl ;
1360 }
1361 }
1362 return ret ;
1363}
1364
1365
1366
1367////////////////////////////////////////////////////////////////////////////////
1368/// Return fundamental (i.e. non-derived) RooAbsArg with given name. Fundamental types
1369/// are e.g. RooRealVar, RooCategory. A null pointer is returned if none is found.
1370
1372{
1373 RooAbsArg* tmp = arg(name) ;
1374 if (!tmp) {
1375 return 0 ;
1376 }
1377 return tmp->isFundamental() ? tmp : 0 ;
1378}
1379
1380
1381
1382////////////////////////////////////////////////////////////////////////////////
1383/// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1384
1386{
1388}
1389
1390
1391////////////////////////////////////////////////////////////////////////////////
1392/// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1393
1395{
1397}
1398
1399
1400
1401
1402////////////////////////////////////////////////////////////////////////////////
1403/// Return set with all variable objects
1404
1406{
1407 RooArgSet ret ;
1408
1409 // Split list of components in pdfs, functions and variables
1411 RooAbsArg* parg ;
1412 while((parg=(RooAbsArg*)iter->Next())) {
1413 if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1414 ret.add(*parg) ;
1415 }
1416 }
1417 delete iter ;
1418
1419 return ret ;
1420}
1421
1422
1423////////////////////////////////////////////////////////////////////////////////
1424/// Return set with all category objects
1425
1427{
1428 RooArgSet ret ;
1429
1430 // Split list of components in pdfs, functions and variables
1432 RooAbsArg* parg ;
1433 while((parg=(RooAbsArg*)iter->Next())) {
1434 if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
1435 ret.add(*parg) ;
1436 }
1437 }
1438 delete iter ;
1439
1440 return ret ;
1441}
1442
1443
1444
1445////////////////////////////////////////////////////////////////////////////////
1446/// Return set with all function objects
1447
1449{
1450 RooArgSet ret ;
1451
1452 // Split list of components in pdfs, functions and variables
1454 RooAbsArg* parg ;
1455 while((parg=(RooAbsArg*)iter->Next())) {
1456 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
1457 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1458 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
1459 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1460 ret.add(*parg) ;
1461 }
1462 }
1463
1464 return ret ;
1465}
1466
1467
1468////////////////////////////////////////////////////////////////////////////////
1469/// Return set with all category function objects
1470
1472{
1473 RooArgSet ret ;
1474
1475 // Split list of components in pdfs, functions and variables
1477 RooAbsArg* parg ;
1478 while((parg=(RooAbsArg*)iter->Next())) {
1479 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
1480 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
1481 ret.add(*parg) ;
1482 }
1483 }
1484 return ret ;
1485}
1486
1487
1488
1489////////////////////////////////////////////////////////////////////////////////
1490/// Return set with all resolution model objects
1491
1493{
1494 RooArgSet ret ;
1495
1496 // Split list of components in pdfs, functions and variables
1498 RooAbsArg* parg ;
1499 while((parg=(RooAbsArg*)iter->Next())) {
1500 if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1501 if (!((RooResolutionModel*)parg)->isConvolved()) {
1502 ret.add(*parg) ;
1503 }
1504 }
1505 }
1506 return ret ;
1507}
1508
1509
1510////////////////////////////////////////////////////////////////////////////////
1511/// Return set with all probability density function objects
1512
1514{
1515 RooArgSet ret ;
1516
1517 // Split list of components in pdfs, functions and variables
1519 RooAbsArg* parg ;
1520 while((parg=(RooAbsArg*)iter->Next())) {
1521 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1522 !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1523 ret.add(*parg) ;
1524 }
1525 }
1526 return ret ;
1527}
1528
1529
1530
1531////////////////////////////////////////////////////////////////////////////////
1532/// Return list of all dataset in the workspace
1533
1534list<RooAbsData*> RooWorkspace::allData() const
1535{
1536 list<RooAbsData*> ret ;
1537 TIterator* iter = _dataList.MakeIterator() ;
1538 RooAbsData* dat ;
1539 while((dat=(RooAbsData*)iter->Next())) {
1540 ret.push_back(dat) ;
1541 }
1542 delete iter ;
1543 return ret ;
1544}
1545
1546
1547////////////////////////////////////////////////////////////////////////////////
1548/// Return list of all dataset in the workspace
1549
1550list<RooAbsData*> RooWorkspace::allEmbeddedData() const
1551{
1552 list<RooAbsData*> ret ;
1554 RooAbsData* dat ;
1555 while((dat=(RooAbsData*)iter->Next())) {
1556 ret.push_back(dat) ;
1557 }
1558 delete iter ;
1559 return ret ;
1560}
1561
1562
1563
1564////////////////////////////////////////////////////////////////////////////////
1565/// Return list of all generic objects in the workspace
1566
1567list<TObject*> RooWorkspace::allGenericObjects() const
1568{
1569 list<TObject*> ret ;
1571 TObject* gobj ;
1572 while((gobj=(RooAbsData*)iter->Next())) {
1573
1574 // If found object is wrapper, return payload
1575 if (gobj->IsA()==RooTObjWrap::Class()) {
1576 ret.push_back(((RooTObjWrap*)gobj)->obj()) ;
1577 } else {
1578 ret.push_back(gobj) ;
1579 }
1580 }
1581 delete iter ;
1582 return ret ;
1583}
1584
1585
1586
1587
1588////////////////////////////////////////////////////////////////////////////////
1589/// Import code of class 'tc' into the repository. If code is already in repository it is only imported
1590/// again if doReplace is false. The names and location of the source files is determined from the information
1591/// in TClass. If no location is found in the TClass information, the files are searched in the workspace
1592/// search path, defined by addClassDeclImportDir() and addClassImplImportDir() for declaration and implementation
1593/// files respectively. If files cannot be found, abort with error status, otherwise update the internal
1594/// class-to-file map and import the contents of the files, if they are not imported yet.
1595
1597{
1598
1599 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") request to import code of class " << tc->GetName() << endl ;
1600
1601 // *** PHASE 1 *** Check if file needs to be imported, or is in ROOT distribution, and check if it can be persisted
1602
1603 // Check if we already have the class (i.e. it is in the classToFile map)
1604 if (!doReplace && _c2fmap.find(tc->GetName())!=_c2fmap.end()) {
1605 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " already imported, skipping" << endl ;
1606 return kTRUE ;
1607 }
1608
1609 // Check if class is listed in a ROOTMAP file - if so we can skip it because it is in the root distribtion
1610 const char* mapEntry = gInterpreter->GetClassSharedLibs(tc->GetName()) ;
1611 if (mapEntry && strlen(mapEntry)>0) {
1612 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << endl ;
1613 return kTRUE ;
1614 }
1615
1616 // Retrieve file names through ROOT TClass interface
1617 string implfile = tc->GetImplFileName() ;
1618 string declfile = tc->GetDeclFileName() ;
1619
1620 // Check that file names are not empty
1621 if (implfile.empty() || declfile.empty()) {
1622 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") ERROR: cannot retrieve code file names for class "
1623 << tc->GetName() << " through ROOT TClass interface, unable to import code" << endl ;
1624 return kFALSE ;
1625 }
1626
1627 // Check if header filename is found in ROOT distribution, if so, do not import class
1628 TString rootsys = gSystem->Getenv("ROOTSYS") ;
1629 if (TString(implfile.c_str()).Index(rootsys)>=0) {
1630 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << endl ;
1631 return kTRUE ;
1632 }
1633 const char* implpath=0 ;
1634
1635 // Require that class meets technical criteria to be persistable (i.e it has a default ctor)
1636 // (We also need a default ctor of abstract classes, but cannot check that through is interface
1637 // as TClass::HasDefaultCtor only returns true for callable default ctors)
1638 if (!(tc->Property() & kIsAbstract) && !tc->HasDefaultConstructor()) {
1639 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING cannot import class "
1640 << tc->GetName() << " : it cannot be persisted because it doesn't have a default constructor. Please fix " << endl ;
1641 return kFALSE ;
1642 }
1643
1644
1645 // *** PHASE 2 *** Check if declaration and implementation files can be located
1646
1647 char* declpath = 0 ;
1648
1649 // Check if header file can be found in specified location
1650 // If not, scan through list of 'class declaration' paths in RooWorkspace
1651 if (gSystem->AccessPathName(declfile.c_str())) {
1652
1653 // Check list of additional declaration paths
1654 list<string>::iterator diter = RooWorkspace::_classDeclDirList.begin() ;
1655
1656 while(diter!= RooWorkspace::_classDeclDirList.end()) {
1657
1658 declpath = gSystem->ConcatFileName(diter->c_str(),declfile.c_str()) ;
1659 if (!gSystem->AccessPathName(declpath)) {
1660 // found declaration file
1661 break ;
1662 }
1663 // cleanup and continue ;
1664 delete[] declpath ;
1665 declpath=0 ;
1666
1667 ++diter ;
1668 }
1669
1670 // Header file cannot be found anywhere, warn user and abort operation
1671 if (!declpath) {
1672 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1673 << tc->GetName() << " because header file " << declfile << " is not found in current directory nor in $ROOTSYS" ;
1674 if (_classDeclDirList.size()>0) {
1675 ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1676 diter = RooWorkspace::_classDeclDirList.begin() ;
1677
1678 while(diter!= RooWorkspace::_classDeclDirList.end()) {
1679
1680 if (diter!=RooWorkspace::_classDeclDirList.begin()) {
1682 }
1683 ooccoutW(_wspace,ObjectHandling) << diter->c_str() ;
1684 ++diter ;
1685 }
1686 }
1687 ooccoutW(_wspace,ObjectHandling) << ". To fix this problem add the required directory to the search "
1688 << "path using RooWorkspace::addClassDeclDir(const char* dir)" << endl ;
1689
1690 return kFALSE ;
1691 }
1692 }
1693
1694
1695 // Check if implementation file can be found in specified location
1696 // If not, scan through list of 'class implementation' paths in RooWorkspace
1697 if (gSystem->AccessPathName(implfile.c_str())) {
1698
1699 // Check list of additional declaration paths
1700 list<string>::iterator iiter = RooWorkspace::_classImplDirList.begin() ;
1701
1702 while(iiter!= RooWorkspace::_classImplDirList.end()) {
1703
1704 implpath = gSystem->ConcatFileName(iiter->c_str(),implfile.c_str()) ;
1705 if (!gSystem->AccessPathName(implpath)) {
1706 // found implementation file
1707 break ;
1708 }
1709 // cleanup and continue ;
1710 delete[] implpath ;
1711 implpath=0 ;
1712
1713 ++iiter ;
1714 }
1715
1716 // Implementation file cannot be found anywhere, warn user and abort operation
1717 if (!implpath) {
1718 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1719 << tc->GetName() << " because implementation file " << implfile << " is not found in current directory nor in $ROOTSYS" ;
1720 if (_classDeclDirList.size()>0) {
1721 ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1722 iiter = RooWorkspace::_classImplDirList.begin() ;
1723
1724 while(iiter!= RooWorkspace::_classImplDirList.end()) {
1725
1726 if (iiter!=RooWorkspace::_classImplDirList.begin()) {
1728 }
1729 ooccoutW(_wspace,ObjectHandling) << iiter->c_str() ;
1730 ++iiter ;
1731 }
1732 }
1733 ooccoutW(_wspace,ObjectHandling) << ". To fix this problem add the required directory to the search "
1734 << "path using RooWorkspace::addClassImplDir(const char* dir)" << endl ;
1735 return kFALSE ;
1736 }
1737 }
1738
1739 char buf[64000];
1740
1741 // *** Phase 3 *** Prepare to import code from files into STL string buffer
1742 //
1743 // Code storage is organized in two linked maps
1744 //
1745 // _fmap contains stl strings with code, indexed on declaration file name
1746 //
1747 // _c2fmap contains list of declaration file names and list of base classes
1748 // and is indexed on class name
1749 //
1750 // Phase 3 is skipped if fmap already contains an entry with given filebasename
1751
1752 string declfilename = declpath?gSystem->BaseName(declpath):gSystem->BaseName(declfile.c_str()) ;
1753
1754 // Split in base and extension
1755 int dotpos2 = strrchr(declfilename.c_str(),'.') - declfilename.c_str() ;
1756 string declfilebase = declfilename.substr(0,dotpos2) ;
1757 string declfileext = declfilename.substr(dotpos2+1) ;
1758
1759 list<string> extraHeaders ;
1760
1761 // If file has not beed stored yet, enter stl strings with implementation and declaration in file map
1762 if (_fmap.find(declfilebase) == _fmap.end()) {
1763
1764 // Open declaration file
1765 fstream fdecl(declpath?declpath:declfile.c_str()) ;
1766
1767 // Abort import if declaration file cannot be opened
1768 if (!fdecl) {
1769 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1770 << ") ERROR opening declaration file " << declfile << endl ;
1771 return kFALSE ;
1772 }
1773
1774 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1775 << ") importing code of class " << tc->GetName()
1776 << " from " << (implpath?implpath:implfile.c_str())
1777 << " and " << (declpath?declpath:declfile.c_str()) << endl ;
1778
1779
1780 // Read entire file into an stl string
1781 string decl ;
1782 while(fdecl.getline(buf,1023)) {
1783
1784 // Look for include state of self
1785 Bool_t processedInclude = kFALSE ;
1786 char* extincfile = 0 ;
1787
1788 // Look for include of declaration file corresponding to this implementation file
1789 if (strstr(buf,"#include")) {
1790 // Process #include statements here
1791 char tmp[64000];
1792 strlcpy(tmp, buf, 64000);
1793 Bool_t stdinclude = strchr(buf, '<');
1794 strtok(tmp, " <\"");
1795 char *incfile = strtok(0, " <>\"");
1796
1797 if (!stdinclude) {
1798 // check if it lives in $ROOTSYS/include
1799 TString hpath = gSystem->Getenv("ROOTSYS");
1800 hpath += "/include/";
1801 hpath += incfile;
1802 if (gSystem->AccessPathName(hpath.Data())) {
1803 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1804 << ") scheduling include file " << incfile << " for import" << endl;
1805 extraHeaders.push_back(incfile);
1806 extincfile = incfile;
1807 processedInclude = kTRUE;
1808 }
1809 }
1810 }
1811
1812 if (processedInclude) {
1813 decl += "// external include file below retrieved from workspace code storage\n" ;
1814 decl += Form("#include \"%s\"\n",extincfile) ;
1815 } else {
1816 decl += buf ;
1817 decl += '\n' ;
1818 }
1819 }
1820
1821 // Open implementation file
1822 fstream fimpl(implpath?implpath:implfile.c_str()) ;
1823
1824 // Abort import if implementation file cannot be opened
1825 if (!fimpl) {
1826 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1827 << ") ERROR opening implementation file " << implfile << endl ;
1828 return kFALSE ;
1829 }
1830
1831
1832 // Import entire implentation file into stl string
1833 string impl ;
1834 while(fimpl.getline(buf,1023)) {
1835 // Process #include statements here
1836
1837 // Look for include state of self
1838 Bool_t foundSelfInclude=kFALSE ;
1839 Bool_t processedInclude = kFALSE ;
1840 char* extincfile = 0 ;
1841
1842 // Look for include of declaration file corresponding to this implementation file
1843 if (strstr(buf,"#include")) {
1844 // Process #include statements here
1845 char tmp[64000];
1846 strlcpy(tmp, buf, 64000);
1847 Bool_t stdinclude = strchr(buf, '<');
1848 strtok(tmp, " <\"");
1849 char *incfile = strtok(0, " <>\"");
1850
1851 if (strstr(incfile, declfilename.c_str())) {
1852 foundSelfInclude = kTRUE;
1853 }
1854
1855 if (!stdinclude && !foundSelfInclude) {
1856 // check if it lives in $ROOTSYS/include
1857 TString hpath = gSystem->Getenv("ROOTSYS");
1858 hpath += "/include/";
1859 hpath += incfile;
1860
1861 if (gSystem->AccessPathName(hpath.Data())) {
1862 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1863 << ") scheduling include file " << incfile << " for import" << endl;
1864 extraHeaders.push_back(incfile);
1865 extincfile = incfile;
1866 processedInclude = kTRUE;
1867 }
1868 }
1869 }
1870
1871 // Explicitly rewrite include of own declaration file to string
1872 // any directory prefixes, copy all other lines verbatim in stl string
1873 if (foundSelfInclude) {
1874 // If include of self is found, substitute original include
1875 // which may have directory structure with a plain include
1876 impl += "// class declaration include file below retrieved from workspace code storage\n" ;
1877 impl += Form("#include \"%s.%s\"\n",declfilebase.c_str(),declfileext.c_str()) ;
1878 } else if (processedInclude) {
1879 impl += "// external include file below retrieved from workspace code storage\n" ;
1880 impl += Form("#include \"%s\"\n",extincfile) ;
1881 } else {
1882 impl += buf ;
1883 impl += '\n' ;
1884 }
1885 }
1886
1887 // Create entry in file map
1888 _fmap[declfilebase]._hfile = decl ;
1889 _fmap[declfilebase]._cxxfile = impl ;
1890 _fmap[declfilebase]._hext = declfileext ;
1891
1892 // Process extra includes now
1893 for (list<string>::iterator ehiter = extraHeaders.begin() ; ehiter != extraHeaders.end() ; ++ehiter ) {
1894 if (_ehmap.find(*ehiter) == _ehmap.end()) {
1895
1896 ExtraHeader eh ;
1897 eh._hname = ehiter->c_str() ;
1898 fstream fehdr(ehiter->c_str()) ;
1899 string ehimpl ;
1900 char buf2[1024] ;
1901 while(fehdr.getline(buf2,1023)) {
1902
1903 // Look for include of declaration file corresponding to this implementation file
1904 if (strstr(buf2,"#include")) {
1905 // Process #include statements here
1906 char tmp[64000];
1907 strlcpy(tmp, buf2, 64000);
1908 Bool_t stdinclude = strchr(buf, '<');
1909 strtok(tmp, " <\"");
1910 char *incfile = strtok(0, " <>\"");
1911
1912 if (!stdinclude) {
1913 // check if it lives in $ROOTSYS/include
1914 TString hpath = gSystem->Getenv("ROOTSYS");
1915 hpath += "/include/";
1916 hpath += incfile;
1917 if (gSystem->AccessPathName(hpath.Data())) {
1918 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1919 << ") scheduling recursive include file " << incfile << " for import"
1920 << endl;
1921 extraHeaders.push_back(incfile);
1922 }
1923 }
1924 }
1925
1926 ehimpl += buf2;
1927 ehimpl += '\n';
1928 }
1929 eh._hfile = ehimpl.c_str();
1930
1931 _ehmap[ehiter->c_str()] = eh;
1932 }
1933 }
1934
1935 } else {
1936
1937 // Inform that existing file entry is being recycled because it already contained class code
1938 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1939 << ") code of class " << tc->GetName()
1940 << " was already imported from " << (implpath?implpath:implfile.c_str())
1941 << " and " << (declpath?declpath:declfile.c_str()) << endl ;
1942
1943 }
1944
1945
1946 // *** PHASE 4 *** Import stl strings with code into workspace
1947 //
1948 // If multiple classes are declared in a single code unit, there will be
1949 // multiple _c2fmap entries all pointing to the same _fmap entry.
1950
1951 // Make list of all immediate base classes of this class
1952 TString baseNameList ;
1953 TList* bl = tc->GetListOfBases() ;
1954 TIterator* iter = bl->MakeIterator() ;
1955 TBaseClass* base ;
1956 list<TClass*> bases ;
1957 while((base=(TBaseClass*)iter->Next())) {
1958 if (baseNameList.Length()>0) {
1959 baseNameList += "," ;
1960 }
1961 baseNameList += base->GetClassPointer()->GetName() ;
1962 bases.push_back(base->GetClassPointer()) ;
1963 }
1964
1965 // Map class name to above _fmap entries, along with list of base classes
1966 // in _c2fmap
1967 _c2fmap[tc->GetName()]._baseName = baseNameList ;
1968 _c2fmap[tc->GetName()]._fileBase = declfilebase ;
1969
1970 // Recursive store all base classes.
1971 list<TClass*>::iterator biter = bases.begin() ;
1972 while(biter!=bases.end()) {
1973 autoImportClass(*biter,doReplace) ;
1974 ++biter ;
1975 }
1976
1977 // Cleanup
1978 if (implpath) {
1979 delete[] implpath ;
1980 }
1981 if (declpath) {
1982 delete[] declpath ;
1983 }
1984
1985
1986 return kTRUE ;
1987}
1988
1989
1990////////////////////////////////////////////////////////////////////////////////
1991/// Create transient TDirectory representation of this workspace. This directory
1992/// will appear as a subdirectory of the directory that contains the workspace
1993/// and will have the name of the workspace suffixed with "Dir". The TDirectory
1994/// interface is read-only. Any attempt to insert objects into the workspace
1995/// directory representation will result in an error message. Note that some
1996/// ROOT object like TH1 automatically insert themselves into the current directory
1997/// when constructed. This will give error messages when done in a workspace
1998/// directory.
1999
2001{
2002 if (_dir) return kTRUE ;
2003
2004 TString title= Form("TDirectory representation of RooWorkspace %s",GetName()) ;
2005 _dir = new WSDir(GetName(),title.Data(),this) ;
2006
2007 TIterator* iter = componentIterator() ;
2008 RooAbsArg* darg ;
2009 while((darg=(RooAbsArg*)iter->Next())) {
2010 if (darg->IsA() != RooConstVar::Class()) {
2011 _dir->InternalAppend(darg) ;
2012 }
2013 }
2014
2015 return kTRUE ;
2016}
2017
2018
2019
2020////////////////////////////////////////////////////////////////////////////////
2021/// Import a clone of a generic TObject into workspace generic object container. Imported
2022/// object can be retrieved by name through the obj() method. The object is cloned upon
2023/// importation and the input argument does not need to live beyond the import call
2024///
2025/// Returns kTRUE if an error has occurred.
2026
2027Bool_t RooWorkspace::import(TObject& object, Bool_t replaceExisting)
2028{
2029 // First check if object with given name already exists
2030 TObject* oldObj = _genObjects.FindObject(object.GetName()) ;
2031 if (oldObj && !replaceExisting) {
2032 coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
2033 << object.GetName() << " is already in workspace and replaceExisting flag is set to false" << endl ;
2034 return kTRUE ;
2035 }
2036
2037 // Grab the current state of the directory Auto-Add
2038 ROOT::DirAutoAdd_t func = object.IsA()->GetDirectoryAutoAdd();
2039 object.IsA()->SetDirectoryAutoAdd(0);
2041
2042 if (oldObj) {
2043 _genObjects.Replace(oldObj,object.Clone()) ;
2044 delete oldObj ;
2045 } else {
2046 _genObjects.Add(object.Clone()) ;
2047 }
2048
2049 // Reset the state of the directory Auto-Add
2050 object.IsA()->SetDirectoryAutoAdd(func);
2052
2053 return kFALSE ;
2054}
2055
2056
2057
2058
2059////////////////////////////////////////////////////////////////////////////////
2060/// Import a clone of a generic TObject into workspace generic object container.
2061/// The imported object will be stored under the given alias name rather than its
2062/// own name. Imported object can be retrieved its alias name through the obj() method.
2063/// The object is cloned upon importation and the input argument does not need to live beyond the import call
2064/// This method is mostly useful for importing objects that do not have a settable name such as TMatrix
2065///
2066/// Returns kTRUE if an error has occurred.
2067
2068Bool_t RooWorkspace::import(TObject& object, const char* aliasName, Bool_t replaceExisting)
2069{
2070 // First check if object with given name already exists
2071 TObject* oldObj = _genObjects.FindObject(object.GetName()) ;
2072 if (oldObj && !replaceExisting) {
2073 coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
2074 << object.GetName() << " is already in workspace and replaceExisting flag is set to false" << endl ;
2075 return kTRUE ;
2076 }
2077
2079 RooTObjWrap* wrapper = new RooTObjWrap(object.Clone()) ;
2081 wrapper->setOwning(kTRUE) ;
2082 wrapper->SetName(aliasName) ;
2083 wrapper->SetTitle(aliasName) ;
2084
2085 if (oldObj) {
2086 _genObjects.Replace(oldObj,wrapper) ;
2087 delete oldObj ;
2088 } else {
2089 _genObjects.Add(wrapper) ;
2090 }
2091 return kFALSE ;
2092}
2093
2094
2095
2096
2097////////////////////////////////////////////////////////////////////////////////
2098/// Insert RooStudyManager module
2099
2101{
2102 RooAbsStudy* clone = (RooAbsStudy*) study.Clone() ;
2103 _studyMods.Add(clone) ;
2104 return kFALSE ;
2105}
2106
2107
2108
2109
2110////////////////////////////////////////////////////////////////////////////////
2111/// Remove all RooStudyManager modules
2112
2114{
2115 _studyMods.Delete() ;
2116}
2117
2118
2119
2120
2121////////////////////////////////////////////////////////////////////////////////
2122/// Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
2123
2124TObject* RooWorkspace::obj(const char* name) const
2125{
2126 // Try RooAbsArg first
2127 TObject* ret = arg(name) ;
2128 if (ret) return ret ;
2129
2130 // Then try RooAbsData
2131 ret = data(name) ;
2132 if (ret) return ret ;
2133
2134 // Finally try generic object store
2135 return genobj(name) ;
2136}
2137
2138
2139
2140////////////////////////////////////////////////////////////////////////////////
2141/// Return generic object with given name
2142
2144{
2145 // Find object by name
2147
2148 // Exit here if not found
2149 if (!gobj) return 0 ;
2150
2151 // If found object is wrapper, return payload
2152 if (gobj->IsA()==RooTObjWrap::Class()) return ((RooTObjWrap*)gobj)->obj() ;
2153
2154 return gobj ;
2155}
2156
2157
2158
2159////////////////////////////////////////////////////////////////////////////////
2160
2161Bool_t RooWorkspace::cd(const char* path)
2162{
2163 makeDir() ;
2164 return _dir->cd(path) ;
2165}
2166
2167
2168
2169////////////////////////////////////////////////////////////////////////////////
2170/// Save this current workspace into given file
2171
2172Bool_t RooWorkspace::writeToFile(const char* fileName, Bool_t recreate)
2173{
2174 TFile f(fileName,recreate?"RECREATE":"UPDATE") ;
2175 Write() ;
2176 return kFALSE ;
2177}
2178
2179
2180
2181////////////////////////////////////////////////////////////////////////////////
2182/// Return instance to factory tool
2183
2185{
2186 if (_factory) {
2187 return *_factory;
2188 }
2189 cxcoutD(ObjectHandling) << "INFO: Creating RooFactoryWSTool associated with this workspace" << endl ;
2190 _factory = make_unique<RooFactoryWSTool>(*this);
2191 return *_factory;
2192}
2193
2194
2195
2196
2197////////////////////////////////////////////////////////////////////////////////
2198/// Short-hand function for factory()->process(expr) ;
2199
2201{
2202 return factory().process(expr) ;
2203}
2204
2205
2206
2207
2208////////////////////////////////////////////////////////////////////////////////
2209/// Print contents of the workspace
2210
2212{
2213 Bool_t treeMode(kFALSE) ;
2214 Bool_t verbose(kFALSE);
2215 if (TString(opts).Contains("t")) {
2216 treeMode=kTRUE ;
2217 }
2218 if (TString(opts).Contains("v")) {
2219 verbose = kTRUE;
2220 }
2221
2222 cout << endl << "RooWorkspace(" << GetName() << ") " << GetTitle() << " contents" << endl << endl ;
2223
2224 RooAbsArg* parg ;
2225
2226 RooArgSet pdfSet ;
2227 RooArgSet funcSet ;
2228 RooArgSet varSet ;
2229 RooArgSet catfuncSet ;
2230 RooArgSet convResoSet ;
2231 RooArgSet resoSet ;
2232
2233
2234 // Split list of components in pdfs, functions and variables
2236 while((parg=(RooAbsArg*)iter->Next())) {
2237
2238 //---------------
2239
2240 if (treeMode) {
2241
2242 // In tree mode, only add nodes with no clients to the print lists
2243
2244 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class())) {
2245 if (!parg->hasClients()) {
2246 pdfSet.add(*parg) ;
2247 }
2248 }
2249
2250 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2251 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2252 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2253 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2254 if (!parg->hasClients()) {
2255 funcSet.add(*parg) ;
2256 }
2257 }
2258
2259
2260 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2261 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2262 if (!parg->hasClients()) {
2263 catfuncSet.add(*parg) ;
2264 }
2265 }
2266
2267 } else {
2268
2269 if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2270 if (((RooResolutionModel*)parg)->isConvolved()) {
2271 convResoSet.add(*parg) ;
2272 } else {
2273 resoSet.add(*parg) ;
2274 }
2275 }
2276
2277 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2278 !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2279 pdfSet.add(*parg) ;
2280 }
2281
2282 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2283 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2284 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2285 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2286 funcSet.add(*parg) ;
2287 }
2288
2289 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2290 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2291 catfuncSet.add(*parg) ;
2292 }
2293 }
2294
2295 if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2296 varSet.add(*parg) ;
2297 }
2298
2299 if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
2300 varSet.add(*parg) ;
2301 }
2302
2303 }
2304 delete iter ;
2305
2306
2309
2310 if (varSet.getSize()>0) {
2311 varSet.sort() ;
2312 cout << "variables" << endl ;
2313 cout << "---------" << endl ;
2314 cout << varSet << endl ;
2315 cout << endl ;
2316 }
2317
2318 if (pdfSet.getSize()>0) {
2319 cout << "p.d.f.s" << endl ;
2320 cout << "-------" << endl ;
2321 pdfSet.sort() ;
2322 iter = pdfSet.createIterator() ;
2323 while((parg=(RooAbsArg*)iter->Next())) {
2324 if (treeMode) {
2325 parg->printComponentTree() ;
2326 } else {
2327 parg->Print() ;
2328 }
2329 }
2330 delete iter ;
2331 cout << endl ;
2332 }
2333
2334 if (!treeMode) {
2335 if (resoSet.getSize()>0) {
2336 cout << "analytical resolution models" << endl ;
2337 cout << "----------------------------" << endl ;
2338 resoSet.sort() ;
2339 iter = resoSet.createIterator() ;
2340 while((parg=(RooAbsArg*)iter->Next())) {
2341 parg->Print() ;
2342 }
2343 delete iter ;
2344 // iter = convResoSet.createIterator() ;
2345 // while((parg=(RooAbsArg*)iter->Next())) {
2346 // parg->Print() ;
2347 // }
2348 // delete iter ;
2349 cout << endl ;
2350 }
2351 }
2352
2353 if (funcSet.getSize()>0) {
2354 cout << "functions" << endl ;
2355 cout << "--------" << endl ;
2356 funcSet.sort() ;
2357 iter = funcSet.createIterator() ;
2358 while((parg=(RooAbsArg*)iter->Next())) {
2359 if (treeMode) {
2360 parg->printComponentTree() ;
2361 } else {
2362 parg->Print() ;
2363 }
2364 }
2365 delete iter ;
2366 cout << endl ;
2367 }
2368
2369 if (catfuncSet.getSize()>0) {
2370 cout << "category functions" << endl ;
2371 cout << "------------------" << endl ;
2372 catfuncSet.sort() ;
2373 iter = catfuncSet.createIterator() ;
2374 while((parg=(RooAbsArg*)iter->Next())) {
2375 if (treeMode) {
2376 parg->printComponentTree() ;
2377 } else {
2378 parg->Print() ;
2379 }
2380 }
2381 delete iter ;
2382 cout << endl ;
2383 }
2384
2385 if (_dataList.GetSize()>0) {
2386 cout << "datasets" << endl ;
2387 cout << "--------" << endl ;
2388 iter = _dataList.MakeIterator() ;
2389 RooAbsData* data2 ;
2390 while((data2=(RooAbsData*)iter->Next())) {
2391 cout << data2->IsA()->GetName() << "::" << data2->GetName() << *data2->get() << endl ;
2392 }
2393 delete iter ;
2394 cout << endl ;
2395 }
2396
2397 if (_embeddedDataList.GetSize()>0) {
2398 cout << "embedded datasets (in pdfs and functions)" << endl ;
2399 cout << "-----------------------------------------" << endl ;
2401 RooAbsData* data2 ;
2402 while((data2=(RooAbsData*)iter->Next())) {
2403 cout << data2->IsA()->GetName() << "::" << data2->GetName() << *data2->get() << endl ;
2404 }
2405 delete iter ;
2406 cout << endl ;
2407 }
2408
2409 if (_snapshots.GetSize()>0) {
2410 cout << "parameter snapshots" << endl ;
2411 cout << "-------------------" << endl ;
2412 iter = _snapshots.MakeIterator() ;
2413 RooArgSet* snap ;
2414 while((snap=(RooArgSet*)iter->Next())) {
2415 cout << snap->GetName() << " = (" ;
2416 TIterator* aiter = snap->createIterator() ;
2417 RooAbsArg* a ;
2418 Bool_t first(kTRUE) ;
2419 while((a=(RooAbsArg*)aiter->Next())) {
2420 if (first) { first=kFALSE ; } else { cout << "," ; }
2421 cout << a->GetName() << "=" ;
2422 a->printValue(cout) ;
2423 if (a->isConstant()) {
2424 cout << "[C]" ;
2425 }
2426 }
2427 cout << ")" << endl ;
2428 delete aiter ;
2429 }
2430 delete iter ;
2431 cout << endl ;
2432 }
2433
2434
2435 if (_namedSets.size()>0) {
2436 cout << "named sets" << endl ;
2437 cout << "----------" << endl ;
2438 for (map<string,RooArgSet>::const_iterator it = _namedSets.begin() ; it != _namedSets.end() ; ++it) {
2439 if (verbose || !TString(it->first.c_str()).BeginsWith("CACHE_")) {
2440 cout << it->first << ":" << it->second << endl;
2441 }
2442 }
2443
2444 cout << endl ;
2445 }
2446
2447
2448 if (_genObjects.GetSize()>0) {
2449 cout << "generic objects" << endl ;
2450 cout << "---------------" << endl ;
2451 iter = _genObjects.MakeIterator() ;
2452 TObject* gobj ;
2453 while((gobj=(TObject*)iter->Next())) {
2454 if (gobj->IsA()==RooTObjWrap::Class()) {
2455 cout << ((RooTObjWrap*)gobj)->obj()->IsA()->GetName() << "::" << gobj->GetName() << endl ;
2456 } else {
2457 cout << gobj->IsA()->GetName() << "::" << gobj->GetName() << endl ;
2458 }
2459 }
2460 delete iter ;
2461 cout << endl ;
2462
2463 }
2464
2465 if (_studyMods.GetSize()>0) {
2466 cout << "study modules" << endl ;
2467 cout << "-------------" << endl ;
2468 iter = _studyMods.MakeIterator() ;
2469 TObject* smobj ;
2470 while((smobj=(TObject*)iter->Next())) {
2471 cout << smobj->IsA()->GetName() << "::" << smobj->GetName() << endl ;
2472 }
2473 delete iter ;
2474 cout << endl ;
2475
2476 }
2477
2478 if (_classes.listOfClassNames().size()>0) {
2479 cout << "embedded class code" << endl ;
2480 cout << "-------------------" << endl ;
2481 cout << _classes.listOfClassNames() << endl ;
2482 cout << endl ;
2483 }
2484
2485 if (_eocache.size()>0) {
2486 cout << "embedded precalculated expensive components" << endl ;
2487 cout << "-------------------------------------------" << endl ;
2488 _eocache.print() ;
2489 }
2490
2492
2493 return ;
2494}
2495
2496
2497////////////////////////////////////////////////////////////////////////////////
2498/// Custom streamer for the workspace. Stream contents of workspace
2499/// and code repository. When reading, read code repository first
2500/// and compile missing classes before proceeding with streaming
2501/// of workspace contents to avoid errors.
2502
2503void RooWorkspace::CodeRepo::Streamer(TBuffer &R__b)
2504{
2505 typedef ::RooWorkspace::CodeRepo thisClass;
2506
2507 // Stream an object of class RooWorkspace::CodeRepo.
2508 if (R__b.IsReading()) {
2509
2510 UInt_t R__s, R__c;
2511 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2512
2513 // Stream contents of ClassFiles map
2514 Int_t count(0) ;
2515 R__b >> count ;
2516 while(count--) {
2517 TString name ;
2518 name.Streamer(R__b) ;
2519 _fmap[name]._hext.Streamer(R__b) ;
2520 _fmap[name]._hfile.Streamer(R__b) ;
2521 _fmap[name]._cxxfile.Streamer(R__b) ;
2522 }
2523
2524 // Stream contents of ClassRelInfo map
2525 count=0 ;
2526 R__b >> count ;
2527 while(count--) {
2528 TString name ;
2529 name.Streamer(R__b) ;
2530 _c2fmap[name]._baseName.Streamer(R__b) ;
2531 _c2fmap[name]._fileBase.Streamer(R__b) ;
2532 }
2533
2534 if (R__v==2) {
2535
2536 count=0 ;
2537 R__b >> count ;
2538 while(count--) {
2539 TString name ;
2540 name.Streamer(R__b) ;
2541 _ehmap[name]._hname.Streamer(R__b) ;
2542 _ehmap[name]._hfile.Streamer(R__b) ;
2543 }
2544 }
2545
2546 R__b.CheckByteCount(R__s, R__c, thisClass::IsA());
2547
2548 // Instantiate any classes that are not defined in current session
2549 _compiledOK = !compileClasses() ;
2550
2551 } else {
2552
2553 UInt_t R__c;
2554 R__c = R__b.WriteVersion(thisClass::IsA(), kTRUE);
2555
2556 // Stream contents of ClassFiles map
2557 UInt_t count = _fmap.size() ;
2558 R__b << count ;
2559 map<TString,ClassFiles>::iterator iter = _fmap.begin() ;
2560 while(iter!=_fmap.end()) {
2561 TString key_copy(iter->first) ;
2562 key_copy.Streamer(R__b) ;
2563 iter->second._hext.Streamer(R__b) ;
2564 iter->second._hfile.Streamer(R__b);
2565 iter->second._cxxfile.Streamer(R__b);
2566
2567 ++iter ;
2568 }
2569
2570 // Stream contents of ClassRelInfo map
2571 count = _c2fmap.size() ;
2572 R__b << count ;
2573 map<TString,ClassRelInfo>::iterator iter2 = _c2fmap.begin() ;
2574 while(iter2!=_c2fmap.end()) {
2575 TString key_copy(iter2->first) ;
2576 key_copy.Streamer(R__b) ;
2577 iter2->second._baseName.Streamer(R__b) ;
2578 iter2->second._fileBase.Streamer(R__b);
2579 ++iter2 ;
2580 }
2581
2582 // Stream contents of ExtraHeader map
2583 count = _ehmap.size() ;
2584 R__b << count ;
2585 map<TString,ExtraHeader>::iterator iter3 = _ehmap.begin() ;
2586 while(iter3!=_ehmap.end()) {
2587 TString key_copy(iter3->first) ;
2588 key_copy.Streamer(R__b) ;
2589 iter3->second._hname.Streamer(R__b) ;
2590 iter3->second._hfile.Streamer(R__b);
2591 ++iter3 ;
2592 }
2593
2594 R__b.SetByteCount(R__c, kTRUE);
2595
2596 }
2597}
2598
2599
2600////////////////////////////////////////////////////////////////////////////////
2601/// Stream an object of class RooWorkspace. This is a standard ROOT streamer for the
2602/// I/O part. This custom function exists to detach all external client links
2603/// from the payload prior to writing the payload so that these client links
2604/// are not persisted. (Client links occur if external function objects use
2605/// objects contained in the workspace as input)
2606/// After the actual writing, these client links are restored.
2607
2608void RooWorkspace::Streamer(TBuffer &R__b)
2609{
2610 if (R__b.IsReading()) {
2611
2613
2614 // Perform any pass-2 schema evolution here
2616 RooAbsArg* node ;
2617 while((node=fiter.next())) {
2618 node->ioStreamerPass2() ;
2619 }
2621
2622 // Make expensive object cache of all objects point to intermal copy.
2623 // Somehow this doesn't work OK automatically
2625 while((node=(RooAbsArg*)iter->Next())) {
2627 node->setWorkspace(*this);
2628 if (node->IsA()->InheritsFrom(RooAbsOptTestStatistic::Class())) {
2630 if (tmp->isSealed() && tmp->sealNotice() && strlen(tmp->sealNotice()) > 0) {
2631 cout << "RooWorkspace::Streamer(" << GetName() << ") " << node->IsA()->GetName() << "::" << node->GetName()
2632 << " : " << tmp->sealNotice() << endl;
2633 }
2634 }
2635 }
2636 delete iter ;
2637
2638
2639 } else {
2640
2641 // Make lists of external clients of WS objects, and remove those links temporarily
2642
2643 map<RooAbsArg*,list<RooAbsArg*> > extClients, extValueClients, extShapeClients ;
2644
2646 RooAbsArg* tmparg ;
2647 while((tmparg=(RooAbsArg*)iter->Next())) {
2648
2649 // Loop over client list of this arg
2650 TIterator* clientIter = tmparg->_clientList.MakeIterator() ;
2651 RooAbsArg* client ;
2652 while((client=(RooAbsArg*)clientIter->Next())) {
2653 if (!_allOwnedNodes.containsInstance(*client)) {
2654 while(tmparg->_clientList.refCount(client)>0) {
2655 tmparg->_clientList.Remove(client) ;
2656 extClients[tmparg].push_back(client) ;
2657 }
2658 }
2659 }
2660 delete clientIter ;
2661
2662 // Loop over value client list of this arg
2663 TIterator* vclientIter = tmparg->_clientListValue.MakeIterator() ;
2664 RooAbsArg* vclient ;
2665 while((vclient=(RooAbsArg*)vclientIter->Next())) {
2666 if (!_allOwnedNodes.containsInstance(*vclient)) {
2667 cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2668 << " has external value client link to " << vclient << " (" << vclient->GetName() << ") with ref count " << tmparg->_clientListValue.refCount(vclient) << endl ;
2669 while(tmparg->_clientListValue.refCount(vclient)>0) {
2670 tmparg->_clientListValue.Remove(vclient) ;
2671 extValueClients[tmparg].push_back(vclient) ;
2672 }
2673 }
2674 }
2675 delete vclientIter ;
2676
2677 // Loop over shape client list of this arg
2678 TIterator* sclientIter = tmparg->_clientListShape.MakeIterator() ;
2679 RooAbsArg* sclient ;
2680 while((sclient=(RooAbsArg*)sclientIter->Next())) {
2681 if (!_allOwnedNodes.containsInstance(*sclient)) {
2682 cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2683 << " has external shape client link to " << sclient << " (" << sclient->GetName() << ") with ref count " << tmparg->_clientListShape.refCount(sclient) << endl ;
2684 while(tmparg->_clientListShape.refCount(sclient)>0) {
2685 tmparg->_clientListShape.Remove(sclient) ;
2686 extShapeClients[tmparg].push_back(sclient) ;
2687 }
2688 }
2689 }
2690 delete sclientIter ;
2691
2692 }
2693 delete iter ;
2694
2696
2697 // Reinstate clients here
2698
2699
2700 for (map<RooAbsArg*,list<RooAbsArg*> >::iterator iterx = extClients.begin() ; iterx!=extClients.end() ; ++iterx) {
2701 for (list<RooAbsArg*>::iterator citer = iterx->second.begin() ; citer!=iterx->second.end() ; ++citer) {
2702 iterx->first->_clientList.Add(*citer) ;
2703 }
2704 }
2705
2706 for (map<RooAbsArg*,list<RooAbsArg*> >::iterator iterx = extValueClients.begin() ; iterx!=extValueClients.end() ; ++iterx) {
2707 for (list<RooAbsArg*>::iterator citer = iterx->second.begin() ; citer!=iterx->second.end() ; ++citer) {
2708 iterx->first->_clientListValue.Add(*citer) ;
2709 }
2710 }
2711
2712 for (map<RooAbsArg*,list<RooAbsArg*> >::iterator iterx = extShapeClients.begin() ; iterx!=extShapeClients.end() ; ++iterx) {
2713 for (list<RooAbsArg*>::iterator citer = iterx->second.begin() ; citer!=iterx->second.end() ; ++citer) {
2714 iterx->first->_clientListShape.Add(*citer) ;
2715 }
2716 }
2717
2718 }
2719}
2720
2721
2722
2723
2724////////////////////////////////////////////////////////////////////////////////
2725/// Return STL string with last of class names contained in the code repository
2726
2728{
2729 string ret ;
2730 map<TString,ClassRelInfo>::const_iterator iter = _c2fmap.begin() ;
2731 while(iter!=_c2fmap.end()) {
2732 if (ret.size()>0) {
2733 ret += ", " ;
2734 }
2735 ret += iter->first ;
2736 ++iter ;
2737 }
2738
2739 return ret ;
2740}
2741
2742
2743
2744////////////////////////////////////////////////////////////////////////////////
2745/// For all classes in the workspace for which no class definition is
2746/// found in the ROOT class table extract source code stored in code
2747/// repository into temporary directory set by
2748/// setClassFileExportDir(), compile classes and link them with
2749/// current ROOT session. If a compilation error occurs print
2750/// instructions for user how to fix errors and recover workspace and
2751/// abort import procedure.
2752
2754{
2755 Bool_t haveDir=kFALSE ;
2756
2757 // Retrieve name of directory in which to export code files
2758 string dirName = Form(_classFileExportDir.c_str(),_wspace->uuid().AsString(),_wspace->GetName()) ;
2759
2760 Bool_t writeExtraHeaders(kFALSE) ;
2761
2762 // Process all class entries in repository
2763 map<TString,ClassRelInfo>::iterator iter = _c2fmap.begin() ;
2764 while(iter!=_c2fmap.end()) {
2765
2766 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing class " << iter->first.Data() << endl ;
2767
2768 // If class is already known, don't load
2769 if (gClassTable->GetDict(iter->first.Data())) {
2770 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Embedded class "
2771 << iter->first << " already in ROOT class table, skipping" << endl ;
2772 ++iter ;
2773 continue ;
2774 }
2775
2776 // Check that export directory exists
2777 if (!haveDir) {
2778
2779 // If not, make local directory to extract files
2780 if (!gSystem->AccessPathName(dirName.c_str())) {
2781 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() reusing code export directory " << dirName.c_str()
2782 << " to extract coded embedded in workspace" << endl ;
2783 } else {
2784 if (gSystem->MakeDirectory(dirName.c_str())==0) {
2785 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() creating code export directory " << dirName.c_str()
2786 << " to extract coded embedded in workspace" << endl ;
2787 } else {
2788 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR creating code export directory " << dirName.c_str()
2789 << " to extract coded embedded in workspace" << endl ;
2790 return kFALSE ;
2791 }
2792 }
2793 haveDir=kTRUE ;
2794
2795 }
2796
2797 // First write any extra header files
2798 if (!writeExtraHeaders) {
2799 writeExtraHeaders = kTRUE ;
2800
2801 map<TString,ExtraHeader>::iterator eiter = _ehmap.begin() ;
2802 while(eiter!=_ehmap.end()) {
2803
2804 // Check if identical declaration file (header) is already written
2805 Bool_t needEHWrite=kTRUE ;
2806 string fdname = Form("%s/%s",dirName.c_str(),eiter->second._hname.Data()) ;
2807 ifstream ifdecl(fdname.c_str()) ;
2808 if (ifdecl) {
2809 TString contents ;
2810 char buf[64000];
2811 while (ifdecl.getline(buf, 64000)) {
2812 contents += buf;
2813 contents += "\n";
2814 }
2815 UInt_t crcFile = RooAbsArg::crc32(contents.Data());
2816 UInt_t crcWS = RooAbsArg::crc32(eiter->second._hfile.Data());
2817 needEHWrite = (crcFile != crcWS);
2818 }
2819
2820 // Write declaration file if required
2821 if (needEHWrite) {
2822 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting extra header file "
2823 << fdname << endl;
2824
2825 // Extra headers may contain non-existing path - create first to be sure
2826 gSystem->MakeDirectory(gSystem->DirName(fdname.c_str()));
2827
2828 ofstream fdecl(fdname.c_str());
2829 if (!fdecl) {
2830 oocoutE(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file " << fdname
2831 << " for writing" << endl;
2832 return kFALSE;
2833 }
2834 fdecl << eiter->second._hfile.Data();
2835 fdecl.close();
2836 }
2837 ++eiter;
2838 }
2839 }
2840
2841
2842 // Navigate from class to file
2843 ClassFiles& cfinfo = _fmap[iter->second._fileBase] ;
2844
2845 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing file with base " << iter->second._fileBase << endl ;
2846
2847 // If file is already processed, skip to next class
2848 if (cfinfo._extracted) {
2849 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() file with base name " << iter->second._fileBase
2850 << " has already been extracted, skipping to next class" << endl ;
2851 continue ;
2852 }
2853
2854 // Check if identical declaration file (header) is already written
2855 Bool_t needDeclWrite=kTRUE ;
2856 string fdname = Form("%s/%s.%s",dirName.c_str(),iter->second._fileBase.Data(),cfinfo._hext.Data()) ;
2857 ifstream ifdecl(fdname.c_str()) ;
2858 if (ifdecl) {
2859 TString contents ;
2860 char buf[64000];
2861 while (ifdecl.getline(buf, 64000)) {
2862 contents += buf;
2863 contents += "\n";
2864 }
2865 UInt_t crcFile = RooAbsArg::crc32(contents.Data()) ;
2866 UInt_t crcWS = RooAbsArg::crc32(cfinfo._hfile.Data()) ;
2867 needDeclWrite = (crcFile!=crcWS) ;
2868 }
2869
2870 // Write declaration file if required
2871 if (needDeclWrite) {
2872 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting declaration code of class " << iter->first << ", file " << fdname << endl ;
2873 ofstream fdecl(fdname.c_str()) ;
2874 if (!fdecl) {
2875 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file "
2876 << fdname << " for writing" << endl ;
2877 return kFALSE ;
2878 }
2879 fdecl << cfinfo._hfile ;
2880 fdecl.close() ;
2881 }
2882
2883 // Check if identical implementation file is already written
2884 Bool_t needImplWrite=kTRUE ;
2885 string finame = Form("%s/%s.cxx",dirName.c_str(),iter->second._fileBase.Data()) ;
2886 ifstream ifimpl(finame.c_str()) ;
2887 if (ifimpl) {
2888 TString contents ;
2889 char buf[64000];
2890 while (ifimpl.getline(buf, 64000)) {
2891 contents += buf;
2892 contents += "\n";
2893 }
2894 UInt_t crcFile = RooAbsArg::crc32(contents.Data()) ;
2895 UInt_t crcWS = RooAbsArg::crc32(cfinfo._cxxfile.Data()) ;
2896 needImplWrite = (crcFile!=crcWS) ;
2897 }
2898
2899 // Write implementation file if required
2900 if (needImplWrite) {
2901 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting implementation code of class " << iter->first << ", file " << finame << endl ;
2902 ofstream fimpl(finame.c_str()) ;
2903 if (!fimpl) {
2904 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file"
2905 << finame << " for writing" << endl ;
2906 return kFALSE ;
2907 }
2908 fimpl << cfinfo._cxxfile ;
2909 fimpl.close() ;
2910 }
2911
2912 // Mark this file as extracted
2913 cfinfo._extracted = kTRUE ;
2914 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() marking code unit " << iter->second._fileBase << " as extracted" << endl ;
2915
2916 // Compile class
2917 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Compiling code unit " << iter->second._fileBase.Data() << " to define class " << iter->first << endl ;
2918 Bool_t ok = gSystem->CompileMacro(finame.c_str(),"k") ;
2919
2920 if (!ok) {
2921 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR compiling class " << iter->first.Data() << ", to fix this you can do the following: " << endl
2922 << " 1) Fix extracted source code files in directory " << dirName.c_str() << "/" << endl
2923 << " 2) In clean ROOT session compiled fixed classes by hand using '.x " << dirName.c_str() << "/ClassName.cxx+'" << endl
2924 << " 3) Reopen file with RooWorkspace with broken source code in UPDATE mode. Access RooWorkspace to force loading of class" << endl
2925 << " Broken instances in workspace will _not_ be compiled, instead precompiled fixed instances will be used." << endl
2926 << " 4) Reimport fixed code in workspace using 'RooWorkspace::importClassCode(\"*\",kTRUE)' method, Write() updated workspace to file and close file" << endl
2927 << " 5) Reopen file in clean ROOT session to confirm that problems are fixed" << endl ;
2928 return kFALSE ;
2929 }
2930
2931 ++iter ;
2932 }
2933
2934 return kTRUE ;
2935}
2936
2937
2938
2939////////////////////////////////////////////////////////////////////////////////
2940/// Internal access to TDirectory append method
2941
2943{
2944#if ROOT_VERSION_CODE <= ROOT_VERSION(5,19,02)
2946#else
2948#endif
2949
2950}
2951
2952
2953////////////////////////////////////////////////////////////////////////////////
2954/// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
2955
2956#if ROOT_VERSION_CODE <= ROOT_VERSION(5,19,02)
2958#else
2960#endif
2961{
2962 if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
2963 coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << endl ;
2964 } else {
2965 InternalAppend(obj) ;
2966 }
2967}
2968
2969
2970////////////////////////////////////////////////////////////////////////////////
2971/// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
2972
2973#if ROOT_VERSION_CODE <= ROOT_VERSION(5,19,02)
2975#else
2977#endif
2978{
2979 if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
2980 coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << endl ;
2981 } else {
2982 InternalAppend(obj) ;
2983 }
2984}
2985
2986
2987
2988////////////////////////////////////////////////////////////////////////////////
2989/// Activate export of workspace symbols to CINT in a namespace with given name. If no name
2990/// is given the namespace will have the same name as the workspace
2991
2992void RooWorkspace::exportToCint(const char* nsname)
2993{
2994 // If export is already active, do nothing
2995 if (_doExport) {
2996 coutE(ObjectHandling) << "RooWorkspace::exportToCint(" << GetName() << ") WARNING: repeated calls to exportToCint() have no effect" << endl ;
2997 return ;
2998 }
2999
3000 // Set flag so that future import to workspace are automatically exported to CINT
3001 _doExport = kTRUE ;
3002
3003 // If no name is provided choose name of workspace
3004 if (!nsname) nsname = GetName() ;
3005 _exportNSName = nsname ;
3006
3007 coutI(ObjectHandling) << "RooWorkspace::exportToCint(" << GetName()
3008 << ") INFO: references to all objects in this workspace will be created in CINT in 'namespace " << _exportNSName << "'" << endl ;
3009
3010 // Export present contents of workspace to CINT
3012 TObject* wobj ;
3013 while((wobj=iter->Next())) {
3014 exportObj(wobj) ;
3015 }
3016 delete iter ;
3017 iter = _dataList.MakeIterator() ;
3018 while((wobj=iter->Next())) {
3019 exportObj(wobj) ;
3020 }
3021 delete iter ;
3022}
3023
3024
3025////////////////////////////////////////////////////////////////////////////////
3026/// Export reference to given workspace object to CINT
3027
3029{
3030 // Do nothing if export flag is not set
3031 if (!_doExport) return ;
3032
3033 // Do not export RooConstVars
3034 if (wobj->IsA() == RooConstVar::Class()) {
3035 return ;
3036 }
3037
3038
3039 // Determine if object name is a valid C++ identifier name
3040
3041 // Do not export objects that have names that are not valid C++ identifiers
3042 if (!isValidCPPID(wobj->GetName())) {
3043 cxcoutD(ObjectHandling) << "RooWorkspace::exportObj(" << GetName() << ") INFO: Workspace object name " << wobj->GetName() << " is not a valid C++ identifier and is not exported to CINT" << endl ;
3044 return ;
3045 }
3046
3047 // Declare correctly typed reference to object in CINT in the namespace associated with this workspace
3048 string cintExpr = Form("namespace %s { %s& %s = *(%s *)0x%lx ; }",_exportNSName.c_str(),wobj->IsA()->GetName(),wobj->GetName(),wobj->IsA()->GetName(),(ULong_t)wobj) ;
3049 gROOT->ProcessLine(cintExpr.c_str()) ;
3050}
3051
3052
3053
3054////////////////////////////////////////////////////////////////////////////////
3055/// Return true if given name is a valid C++ identifier name
3056
3058{
3059 string oname(name) ;
3060 if (isdigit(oname[0])) {
3061 return kFALSE ;
3062 } else {
3063 for (UInt_t i=0 ; i<oname.size() ; i++) {
3064 char c = oname[i] ;
3065 if (!isalnum(c) && (c!='_')) {
3066 return kFALSE ;
3067 }
3068 }
3069 }
3070 return kTRUE ;
3071}
3072
3073////////////////////////////////////////////////////////////////////////////////
3074/// Delete exported reference in CINT namespace
3075
3077{
3078 char buf[64000];
3080 TObject *wobj;
3081 while ((wobj = iter->Next())) {
3082 if (isValidCPPID(wobj->GetName())) {
3083 strlcpy(buf, Form("%s::%s", _exportNSName.c_str(), wobj->GetName()), 64000);
3084 gInterpreter->DeleteVariable(buf);
3085 }
3086 }
3087 delete iter ;
3088}
3089
3090////////////////////////////////////////////////////////////////////////////////
3091/// If one of the TObject we have a referenced to is deleted, remove the
3092/// reference.
3093
3095{
3096 _dataList.RecursiveRemove(removedObj);
3097 if (removedObj == _dir) _dir = nullptr;
3098
3099 _allOwnedNodes.RecursiveRemove(removedObj); // RooArgSet
3100
3101 _dataList.RecursiveRemove(removedObj);
3103 _views.RecursiveRemove(removedObj);
3104 _snapshots.RecursiveRemove(removedObj);
3105 _genObjects.RecursiveRemove(removedObj);
3106 _studyMods.RecursiveRemove(removedObj);
3107
3108 for(auto c : _namedSets) {
3109 c.second.RecursiveRemove(removedObj);
3110 }
3111
3112 _eocache.RecursiveRemove(removedObj); // RooExpensiveObjectCache
3113}
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:363
R__EXTERN TClassTable * gClassTable
Definition: TClassTable.h:95
@ kIsAbstract
Definition: TDictionary.h:71
#define gInterpreter
Definition: TInterpreter.h:538
#define gROOT
Definition: TROOT.h:410
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
RooExpensiveObjectCache & expensiveObjectCache() const
Definition: RooAbsArg.cxx:2334
virtual Bool_t importWorkspaceHook(RooWorkspace &ws)
Definition: RooAbsArg.h:501
static UInt_t crc32(const char *data)
Definition: RooAbsArg.cxx:1898
Bool_t redirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t isRecursionStep=kFALSE)
Iterator over _clientListValue.
Definition: RooAbsArg.cxx:921
RooRefCountList _clientListValue
Definition: RooAbsArg.h:478
RooRefCountList _clientListShape
Definition: RooAbsArg.h:477
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
Definition: RooAbsArg.cxx:273
static void ioStreamerPass2Finalize()
Method called by workspace container to finalize schema evolution issues that cannot be handled in a ...
Definition: RooAbsArg.cxx:2467
void setWorkspace(RooWorkspace &ws)
Definition: RooAbsArg.h:418
void printComponentTree(const char *indent="", const char *namePat=0, Int_t nLevel=999)
Print tree structure of expression tree on given ostream, only branch nodes are printed.
Definition: RooAbsArg.cxx:1837
virtual void setExpensiveObjectCache(RooExpensiveObjectCache &cache)
Definition: RooAbsArg.h:499
virtual Bool_t isFundamental() const
Definition: RooAbsArg.h:157
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
Definition: RooAbsArg.cxx:286
void setAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
Definition: RooAbsArg.cxx:241
Bool_t hasClients() const
Definition: RooAbsArg.h:99
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2074
virtual void ioStreamerPass2()
In which workspace do I live, if any.
Definition: RooAbsArg.cxx:2437
Bool_t getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
Definition: RooAbsArg.cxx:264
void branchNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all branch nodes of the arg tree starting with ourself as top node.
Definition: RooAbsArg.cxx:477
void SetName(const char *name)
Set the name of the TNamed.
Definition: RooAbsArg.cxx:2382
RooRefCountList _clientList
Definition: RooAbsArg.h:476
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Bool_t containsInstance(const RooAbsArg &var) const
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
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooFIter fwdIterator() const
void sort(Bool_t ascend=kTRUE)
const char * GetName() const
Returns name of object.
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.
void setHashTableSize(Int_t i)
TIterator * createIterator(Bool_t dir=kIterForward) const
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:79
virtual Bool_t changeObservableName(const char *from, const char *to)
Definition: RooAbsData.cxx:265
RooAbsOptTestStatistic is the abstract base class for test statistics objects that evaluate a functio...
const char * sealNotice() const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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
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)
RooAbsArg * next()
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:35
Int_t GetSize() const
Definition: RooLinkedList.h:60
TObject * FindObject(const char *name) const
Return pointer to obejct with given name.
TIterator * MakeIterator(Bool_t dir=kTRUE) const
Return an iterator over 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:62
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
Int_t refCount(TObject *obj)
Return reference count associated with 'obj'.
virtual Bool_t Remove(TObject *obj)
Remove object from list and if reference count reaches zero delete object itself as well.
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 containgin 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:40
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition: TBuffer.h:83
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
static DictFuncPtr_t GetDict(const char *cname)
Given the class name returns the Dictionary() function of a class (uses hash of name).
The ROOT global object gROOT contains a list of all defined classes.
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:3505
Long_t Property() const
Set TObject::fBits and fStreamerType to cache information about the class.
Definition: TClass.cxx:5800
const char * GetDeclFileName() const
Definition: TClass.h:399
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
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:3975
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:718
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
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:114
@ 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
STL namespace.
auto * a
Definition: textangle.C:12