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