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