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 CINT 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 // Make final check list of conflicting nodes
629 for (const auto branch2 : branchSet2) {
630 if (_allOwnedNodes.find(branch2->GetName())) {
631 conflictNodes2.add(*branch2) ;
632 }
633 }
634
635 // Terminate here if there are conflicts and no resolution protocol
636 if (!conflictNodes2.empty()) {
637 coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR object named " << inArg.GetName() << ": component(s) "
638 << conflictNodes2 << " cause naming conflict after conflict resolution protocol was executed" << std::endl ;
639 return true ;
640 }
641
642 // Perform any auxiliary imports at this point
643 for (const auto node : cloneSet2) {
644 if (node->importWorkspaceHook(*this)) {
645 coutE(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") ERROR object named " << node->GetName()
646 << " has an error in importing in one or more of its auxiliary objects, aborting" << std::endl ;
647 return true ;
648 }
649 }
650
653 for (const auto node : cloneSet2) {
654 if (_autoClass) {
655 if (!_classes.autoImportClass(node->IsA())) {
656 coutW(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") WARNING: problems import class code of object "
657 << node->ClassName() << "::" << node->GetName() << ", reading of workspace will require external definition of class" << std::endl ;
658 }
659 }
660
661 // Point expensiveObjectCache to copy in this workspace
662 RooExpensiveObjectCache& oldCache = node->expensiveObjectCache() ;
663 node->setExpensiveObjectCache(_eocache) ;
664 _eocache.importCacheObjects(oldCache,node->GetName(),true) ;
665
666 // Check if node is already in workspace (can only happen for variables or identical instances, unless RecycleConflictNodes is specified)
667 RooAbsArg* wsnode = _allOwnedNodes.find(node->GetName()) ;
668
669 if (wsnode) {
670 // Do not import node, add not to list of nodes that require reconnection
671 if (!silence && useExistingNodes) {
672 coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") using existing copy of " << node->ClassName()
673 << "::" << node->GetName() << " for import of " << cloneTop2->ClassName() << "::"
674 << cloneTop2->GetName() << std::endl ;
675 }
676 recycledNodes.add(*_allOwnedNodes.find(node->GetName())) ;
677
678 // Delete clone of incoming node
679 nodesToBeDeleted.addOwned(std::unique_ptr<RooAbsArg>{node});
680
681 //cout << "WV: recycling existing node " << existingNode << " = " << existingNode->GetName() << " for imported node " << node << std::endl ;
682
683 } else {
684 // Import node
685 if (!silence) {
686 coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") importing " << node->ClassName() << "::"
687 << node->GetName() << std::endl ;
688 }
689 _allOwnedNodes.addOwned(std::unique_ptr<RooAbsArg>{node});
690 node->setWorkspace(*this);
691 if (_openTrans) {
692 _sandboxNodes.add(*node) ;
693 } else {
694 if (_dir && node->IsA() != RooConstVar::Class()) {
695 _dir->InternalAppend(node) ;
696 }
697 }
698 }
699 }
700
701 // Reconnect any nodes that need to be
702 if (!recycledNodes.empty()) {
703 for (const auto node : cloneSet2) {
704 node->redirectServers(recycledNodes) ;
705 }
706 }
707
708 cloneSet2.releaseOwnership() ;
709
710 return false ;
711}
712
713
714
715////////////////////////////////////////////////////////////////////////////////
716/// Import a dataset (RooDataSet or RooDataHist) into the workspace. The workspace will contain a copy of the data.
717/// The dataset and its variables can be renamed upon insertion with the options below
718///
719/// <table>
720/// <tr><th> Accepted arguments
721/// <tr><td> `Rename(const char* suffix)` <td> Rename dataset upon insertion
722/// <tr><td> `RenameVariable(const char* inputName, const char* outputName)` <td> Change names of observables in dataset upon insertion
723/// <tr><td> `Silence` <td> Be quiet, except in case of errors
724/// \note From python, use `Import()`, since `import` is a reserved keyword.
725/// \return due to historical reasons: false (0) on success and true (1) on failure
727 const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,
728 const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6,
729 const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9)
730
731{
732
733 RooLinkedList args ;
734 args.Add((TObject*)&arg1) ;
735 args.Add((TObject*)&arg2) ;
736 args.Add((TObject*)&arg3) ;
737 args.Add((TObject*)&arg4) ;
738 args.Add((TObject*)&arg5) ;
739 args.Add((TObject*)&arg6) ;
740 args.Add((TObject*)&arg7) ;
741 args.Add((TObject*)&arg8) ;
742 args.Add((TObject*)&arg9) ;
743
744 // Select the pdf-specific commands
745 RooCmdConfig pc(Form("RooWorkspace::import(%s)",GetName())) ;
746
747 pc.defineString("dsetName","Rename",0,"") ;
748 pc.defineString("varChangeIn","RenameVar",0,"",true) ;
749 pc.defineString("varChangeOut","RenameVar",1,"",true) ;
750 pc.defineInt("embedded","Embedded",0,0) ;
751 pc.defineInt("silence","Silence",0,0) ;
752
753 // Process and check varargs
754 pc.process(args) ;
755 if (!pc.ok(true)) {
756 return true ;
757 }
758
759 // Decode renaming logic into suffix string and boolean for conflictOnly mode
760 const char* dsetName = pc.getString("dsetName") ;
761 const char* varChangeIn = pc.getString("varChangeIn") ;
762 const char* varChangeOut = pc.getString("varChangeOut") ;
763 bool embedded = pc.getInt("embedded") ;
764 Int_t silence = pc.getInt("silence") ;
765
766 if (!silence)
767 coutI(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") importing dataset " << inData.GetName() << std::endl ;
768
769 // Transform empty string into null pointer
770 if (dsetName && strlen(dsetName)==0) {
771 dsetName=nullptr ;
772 }
773
775 if (dataList.size() > 50 && dataList.getHashTableSize() == 0) {
776 // When the workspaces get larger, traversing the linked list becomes a bottleneck:
777 dataList.setHashTableSize(200);
778 }
779
780 // Check that no dataset with target name already exists
781 if (dsetName && dataList.FindObject(dsetName)) {
782 coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << dsetName << " already exists in workspace, import aborted" << std::endl ;
783 return true ;
784 }
785 if (!dsetName && dataList.FindObject(inData.GetName())) {
786 coutE(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") ERROR dataset with name " << inData.GetName() << " already exists in workspace, import aborted" << std::endl ;
787 return true ;
788 }
789
790 // Rename dataset if required
791 RooAbsData* clone ;
792 if (dsetName) {
793 if (!silence)
794 coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset from " << inData.GetName() << " to " << dsetName << std::endl ;
795 clone = static_cast<RooAbsData*>(inData.Clone(dsetName)) ;
796 } else {
797 clone = static_cast<RooAbsData*>(inData.Clone(inData.GetName())) ;
798 }
799
800
801 // Process any change in variable names
802 if (strlen(varChangeIn)>0) {
803 // Parse comma separated lists of variable name changes
804 const std::vector<std::string> tokIn = ROOT::Split(varChangeIn, ",");
805 const std::vector<std::string> tokOut = ROOT::Split(varChangeOut, ",");
806 for (unsigned int i=0; i < tokIn.size(); ++i) {
807 if (!silence)
808 coutI(ObjectHandling) << "RooWorkSpace::import(" << GetName() << ") changing name of dataset observable " << tokIn[i] << " to " << tokOut[i] << std::endl ;
809 clone->changeObservableName(tokIn[i].c_str(), tokOut[i].c_str());
810 }
811 }
812
813 // Now import the dataset observables, unless dataset is embedded
814 if (!embedded) {
815 for(RooAbsArg* carg : *clone->get()) {
816 if (!arg(carg->GetName())) {
817 import(*carg) ;
818 }
819 }
820 }
821
822 dataList.Add(clone) ;
823 if (_dir) {
824 _dir->InternalAppend(clone) ;
825 }
826
827 // Set expensive object cache of dataset internal buffers to that of workspace
828 for(RooAbsArg* carg : *clone->get()) {
829 carg->setExpensiveObjectCache(expensiveObjectCache()) ;
830 }
831
832
833 return false ;
834}
835
836
837
838
839////////////////////////////////////////////////////////////////////////////////
840/// Define a named RooArgSet with given constituents. If importMissing is true, any constituents
841/// of aset that are not in the workspace will be imported, otherwise an error is returned
842/// for missing components
843/// \return due to historical reasons: false (0) on success and true (1) on failure
845{
846 // Check if set was previously defined, if so print warning
847 map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
848 if (i!=_namedSets.end()) {
849 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << std::endl ;
850 }
851
853
854 // Check all constituents of provided set
855 for (RooAbsArg* sarg : aset) {
856 // If missing, either import or report error
857 if (!arg(sarg->GetName())) {
858 if (importMissing) {
859 import(*sarg) ;
860 } else {
861 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR set constituent \"" << sarg->GetName()
862 << "\" is not in workspace and importMissing option is disabled" << std::endl ;
863 return true ;
864 }
865 }
866 wsargs.add(*arg(sarg->GetName())) ;
867 }
868
869
870 // Install named set
871 _namedSets[name].removeAll() ;
872 _namedSets[name].add(wsargs) ;
873
874 return false ;
875}
876
877//_____________________________________________________________________________
878// \return due to historical reasons: false (0) on success (always)
880{
881 // Define a named RooArgSet with given constituents. If importMissing is true, any constituents
882 // of aset that are not in the workspace will be imported, otherwise an error is returned
883 // for missing components
884
885 // Check if set was previously defined, if so print warning
886 map<string, RooArgSet>::iterator i = _namedSets.find(name);
887 if (i != _namedSets.end()) {
888 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName()
889 << ") WARNING redefining previously defined named set " << name << std::endl;
890 }
891
892 // Install named set
893 _namedSets[name].removeAll();
894 _namedSets[name].add(aset);
895
896 return false;
897}
898
899////////////////////////////////////////////////////////////////////////////////
900/// Define a named set in the workspace through a comma separated list of
901/// names of objects already in the workspace
902/// \return due to historical reasons: false (0) on success and true (1) on failure
903bool RooWorkspace::defineSet(const char* name, const char* contentList)
904{
905 // Check if set was previously defined, if so print warning
906 map<string,RooArgSet>::iterator i = _namedSets.find(name) ;
907 if (i!=_namedSets.end()) {
908 coutW(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") WARNING redefining previously defined named set " << name << std::endl ;
909 }
910
912
913 // Check all constituents of provided set
914 for (const std::string& token : ROOT::Split(contentList, ",")) {
915 // If missing, either import or report error
916 if (!arg(token.c_str())) {
917 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
918 << "\" is not in workspace" << std::endl ;
919 return true ;
920 }
921 wsargs.add(*arg(token.c_str())) ;
922 }
923
924 // Install named set
925 _namedSets[name].removeAll() ;
926 _namedSets[name].add(wsargs) ;
927
928 return false ;
929}
930
931
932
933
934////////////////////////////////////////////////////////////////////////////////
935/// Define a named set in the workspace through a comma separated list of
936/// names of objects already in the workspace
937/// \return due to historical reasons: false (0) on success and true (1) on failure
938bool RooWorkspace::extendSet(const char* name, const char* newContents)
939{
941
942 // Check all constituents of provided set
943 for (const std::string& token : ROOT::Split(newContents, ",")) {
944 // If missing, either import or report error
945 if (!arg(token.c_str())) {
946 coutE(InputArguments) << "RooWorkspace::defineSet(" << GetName() << ") ERROR proposed set constituent \"" << token
947 << "\" is not in workspace" << std::endl ;
948 return true ;
949 }
950 wsargs.add(*arg(token.c_str())) ;
951 }
952
953 // Extend named set
954 _namedSets[name].add(wsargs,true) ;
955
956 return false ;
957}
958
959
960
961////////////////////////////////////////////////////////////////////////////////
962/// Return pointer to previously defined named set with given nmame
963/// If no such set is found a null pointer is returned
964
966{
967 std::map<string,RooArgSet>::iterator i = _namedSets.find(name.c_str());
968 return (i!=_namedSets.end()) ? &(i->second) : nullptr;
969}
970
971
972
973
974////////////////////////////////////////////////////////////////////////////////
975/// Rename set to a new name
976/// \return due to historical reasons: false (0) on success and true (1) on failure
977bool RooWorkspace::renameSet(const char* name, const char* newName)
978{
979 // First check if set exists
980 if (!set(name)) {
981 coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << name
982 << " does not exist" << std::endl ;
983 return true ;
984 }
985
986 // Check if no set exists with new name
987 if (set(newName)) {
988 coutE(InputArguments) << "RooWorkspace::renameSet(" << GetName() << ") ERROR a set with name " << newName
989 << " already exists" << std::endl ;
990 return true ;
991 }
992
993 // Copy entry under 'name' to 'newName'
995
996 // Remove entry under old name
997 _namedSets.erase(name) ;
998
999 return false ;
1000}
1001
1002
1003
1004
1005////////////////////////////////////////////////////////////////////////////////
1006/// Remove a named set from the workspace
1007/// \return due to historical reasons: false (0) on success and true (1) on failure
1009{
1010 // First check if set exists
1011 if (!set(name)) {
1012 coutE(InputArguments) << "RooWorkspace::removeSet(" << GetName() << ") ERROR a set with name " << name
1013 << " does not exist" << std::endl ;
1014 return true ;
1015 }
1016
1017 // Remove set with given name
1018 _namedSets.erase(name) ;
1019
1020 return false ;
1021}
1022
1023
1024
1025
1026////////////////////////////////////////////////////////////////////////////////
1027/// Open an import transaction operations.
1028/// \return true if successful, false if there is already an ongoing transaction
1029
1031{
1032 // Check that there was no ongoing transaction
1033 if (_openTrans) {
1034 return false ;
1035 }
1036
1037 // Open transaction
1038 _openTrans = true ;
1039 return true ;
1040}
1041
1042
1043
1044
1045////////////////////////////////////////////////////////////////////////////////
1046/// Cancel an ongoing import transaction. All objects imported since startTransaction()
1047/// will be removed and the transaction will be terminated.
1048/// \return true if cancel operation succeeds, return false if there was no open transaction
1049
1051{
1052 // Check that there is an ongoing transaction
1053 if (!_openTrans) {
1054 return false ;
1055 }
1056
1057 // Delete all objects in the sandbox
1058 for(RooAbsArg * tmpArg : _sandboxNodes) {
1060 }
1062
1063 // Mark transaction as finished
1064 _openTrans = false ;
1065
1066 return true ;
1067}
1068
1069/// Commit an ongoing import transaction.
1070/// \return true if commit succeeded, return false if there was no ongoing transaction
1072{
1073 // Check that there is an ongoing transaction
1074 if (!_openTrans) {
1075 return false ;
1076 }
1077
1078 // Publish sandbox nodes in directory and/or CINT if requested
1079 for(RooAbsArg* sarg : _sandboxNodes) {
1080 if (_dir && sarg->IsA() != RooConstVar::Class()) {
1082 }
1083 }
1084
1085 // Remove all committed objects from the sandbox
1087
1088 // Mark transaction as finished
1089 _openTrans = false ;
1090
1091 return true ;
1092}
1093
1094
1095
1096
1097////////////////////////////////////////////////////////////////////////////////
1098/// \return true on success, false on failure
1099/// \see RooWorkspace::CodeRepo::autoImportClass
1104
1105
1106
1107////////////////////////////////////////////////////////////////////////////////
1108/// Import code of all classes in the workspace that have a class name
1109/// that matches pattern 'pat' and which are not found to be part of
1110/// the standard ROOT distribution. If doReplace is true any existing
1111/// class code saved in the workspace is replaced
1112/// \return true on success, false on failure
1114{
1115 bool ret(true) ;
1116
1117 TRegexp re(pat,true) ;
1118 for (RooAbsArg * carg : _allOwnedNodes) {
1119 TString className = carg->ClassName() ;
1120 if (className.Index(re)>=0 && !_classes.autoImportClass(carg->IsA(),doReplace)) {
1121 coutW(ObjectHandling) << "RooWorkspace::import(" << GetName() << ") WARNING: problems import class code of object "
1122 << carg->ClassName() << "::" << carg->GetName() << ", reading of workspace will require external definition of class" << std::endl ;
1123 ret = false ;
1124 }
1125 }
1126
1127 return ret ;
1128}
1129
1130
1131
1132
1133
1134////////////////////////////////////////////////////////////////////////////////
1135/// Save snapshot of values and attributes (including "Constant") of given parameters.
1136/// \param[in] name Name of the snapshot.
1137/// \param[in] paramNames Comma-separated list of parameter names to be snapshot.
1138/// \return true always (success)
1140{
1141 return saveSnapshot(name,argSet(paramNames),false) ;
1142}
1143
1144
1145
1146
1147
1148////////////////////////////////////////////////////////////////////////////////
1149/// Save snapshot of values and attributes (including "Constant") of parameters 'params'.
1150/// If importValues is FALSE, the present values from the object in the workspace are
1151/// saved. If importValues is TRUE, the values of the objects passed in the 'params'
1152/// argument are saved
1153/// \return true always (success)
1155{
1158 auto snapshot = new RooArgSet;
1159 actualParams.snapshot(*snapshot);
1160
1161 snapshot->setName(name.c_str()) ;
1162
1163 if (importValues) {
1164 snapshot->assign(params) ;
1165 }
1166
1167 if (std::unique_ptr<RooArgSet> oldSnap{static_cast<RooArgSet*>(_snapshots.FindObject(name.c_str()))}) {
1168 coutI(ObjectHandling) << "RooWorkspace::saveSnapshot(" << GetName() << ") replacing previous snapshot with name " << name << std::endl ;
1169 _snapshots.Remove(oldSnap.get()) ;
1170 }
1171
1172 _snapshots.Add(snapshot) ;
1173
1174 return true ;
1175}
1176
1177
1178
1179
1180////////////////////////////////////////////////////////////////////////////////
1181/// Load the values and attributes of the parameters in the snapshot saved with
1182/// the given name
1183/// \return true on success, false on failure
1185{
1186 RooArgSet* snap = static_cast<RooArgSet*>(_snapshots.find(name)) ;
1187 if (!snap) {
1188 coutE(ObjectHandling) << "RooWorkspace::loadSnapshot(" << GetName() << ") no snapshot with name " << name << " is available" << std::endl ;
1189 return false ;
1190 }
1191
1194 actualParams.assign(*snap) ;
1195
1196 return true ;
1197}
1198
1199
1200////////////////////////////////////////////////////////////////////////////////
1201/// Return the RooArgSet containing a snapshot of variables contained in the workspace
1202///
1203/// Note that the variables of the objects in the snapshots are **copies** of the
1204/// variables in the workspace. To load the values of a snapshot in the workspace
1205/// variables, use loadSnapshot() instead.
1206
1208{
1209 return static_cast<RooArgSet*>(_snapshots.find(name));
1210}
1211
1212
1213////////////////////////////////////////////////////////////////////////////////
1214/// Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found
1215
1217{
1218 return dynamic_cast<RooAbsPdf*>(_allOwnedNodes.find(name.c_str())) ;
1219}
1220
1221
1222////////////////////////////////////////////////////////////////////////////////
1223/// Retrieve function (RooAbsReal) with given name. Note that all RooAbsPdfs are also RooAbsReals. A null pointer is returned if not found.
1224
1226{
1227 return dynamic_cast<RooAbsReal*>(_allOwnedNodes.find(name.c_str())) ;
1228}
1229
1230
1231////////////////////////////////////////////////////////////////////////////////
1232/// Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found
1233
1235{
1236 return dynamic_cast<RooRealVar*>(_allOwnedNodes.find(name.c_str())) ;
1237}
1238
1239
1240////////////////////////////////////////////////////////////////////////////////
1241/// Retrieve discrete variable (RooCategory) with given name. A null pointer is returned if not found
1242
1244{
1245 return dynamic_cast<RooCategory*>(_allOwnedNodes.find(name.c_str())) ;
1246}
1247
1248
1249////////////////////////////////////////////////////////////////////////////////
1250/// Retrieve discrete function (RooAbsCategory) with given name. A null pointer is returned if not found
1251
1253{
1254 return dynamic_cast<RooAbsCategory*>(_allOwnedNodes.find(name.c_str())) ;
1255}
1256
1257
1258
1259////////////////////////////////////////////////////////////////////////////////
1260/// Return RooAbsArg with given name. A null pointer is returned if none is found.
1261
1263{
1264 return _allOwnedNodes.find(name.c_str()) ;
1265}
1266
1267
1268
1269////////////////////////////////////////////////////////////////////////////////
1270/// Return set of RooAbsArgs matching to given list of names
1271
1273{
1274 RooArgSet ret ;
1275
1276 for (const std::string& token : ROOT::Split(nameList, ",")) {
1277 RooAbsArg* oneArg = arg(token.c_str()) ;
1278 if (oneArg) {
1279 ret.add(*oneArg) ;
1280 } else {
1281 std::stringstream ss;
1282 ss << " RooWorkspace::argSet(" << GetName() << ") no RooAbsArg named \"" << token << "\" in workspace" ;
1283 const std::string errorMsg = ss.str();
1284 coutE(InputArguments) << errorMsg << std::endl;
1285 throw std::runtime_error(errorMsg);
1286 }
1287 }
1288 return ret ;
1289}
1290
1291
1292
1293////////////////////////////////////////////////////////////////////////////////
1294/// Return fundamental (i.e. non-derived) RooAbsArg with given name. Fundamental types
1295/// are e.g. RooRealVar, RooCategory. A null pointer is returned if none is found.
1296
1298{
1299 RooAbsArg* tmp = arg(name) ;
1300 if (!tmp) {
1301 return nullptr;
1302 }
1303 return tmp->isFundamental() ? tmp : nullptr;
1304}
1305
1306
1307
1308////////////////////////////////////////////////////////////////////////////////
1309/// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1310
1312{
1313 return static_cast<RooAbsData*>(_dataList.FindObject(name.c_str())) ;
1314}
1315
1316
1317////////////////////////////////////////////////////////////////////////////////
1318/// Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found
1319
1321{
1322 return static_cast<RooAbsData*>(_embeddedDataList.FindObject(name.c_str())) ;
1323}
1324
1325
1326
1327
1328////////////////////////////////////////////////////////////////////////////////
1329/// Return set with all variable objects
1330
1332{
1333 RooArgSet ret ;
1334
1335 // Split list of components in pdfs, functions and variables
1336 for(RooAbsArg* parg : _allOwnedNodes) {
1337 if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1338 ret.add(*parg) ;
1339 }
1340 }
1341
1342 return ret ;
1343}
1344
1345
1346////////////////////////////////////////////////////////////////////////////////
1347/// Return set with all category objects
1348
1350{
1351 RooArgSet ret ;
1352
1353 // Split list of components in pdfs, functions and variables
1354 for(RooAbsArg* parg : _allOwnedNodes) {
1355 if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
1356 ret.add(*parg) ;
1357 }
1358 }
1359
1360 return ret ;
1361}
1362
1363
1364
1365////////////////////////////////////////////////////////////////////////////////
1366/// Return set with all function objects
1367
1369{
1370 RooArgSet ret ;
1371
1372 // Split list of components in pdfs, functions and variables
1373 for(RooAbsArg* parg : _allOwnedNodes) {
1374 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
1375 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1376 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
1377 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
1378 ret.add(*parg) ;
1379 }
1380 }
1381
1382 return ret ;
1383}
1384
1385
1386////////////////////////////////////////////////////////////////////////////////
1387/// Return set with all category function objects
1388
1390{
1391 RooArgSet ret ;
1392
1393 // Split list of components in pdfs, functions and variables
1394 for(RooAbsArg* parg : _allOwnedNodes) {
1395 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
1396 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
1397 ret.add(*parg) ;
1398 }
1399 }
1400 return ret ;
1401}
1402
1403
1404
1405////////////////////////////////////////////////////////////////////////////////
1406/// Return set with all resolution model objects
1407
1409{
1410 RooArgSet ret ;
1411
1412 // Split list of components in pdfs, functions and variables
1413 for(RooAbsArg* parg : _allOwnedNodes) {
1414 if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1415 if (!(static_cast<RooResolutionModel*>(parg))->isConvolved()) {
1416 ret.add(*parg) ;
1417 }
1418 }
1419 }
1420 return ret ;
1421}
1422
1423
1424////////////////////////////////////////////////////////////////////////////////
1425/// Return set with all probability density function objects
1426
1428{
1429 RooArgSet ret ;
1430
1431 // Split list of components in pdfs, functions and variables
1432 for(RooAbsArg* parg : _allOwnedNodes) {
1433 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
1434 !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
1435 ret.add(*parg) ;
1436 }
1437 }
1438 return ret ;
1439}
1440
1441
1442
1443////////////////////////////////////////////////////////////////////////////////
1444/// Return list of all dataset in the workspace
1445
1446std::list<RooAbsData*> RooWorkspace::allData() const
1447{
1448 std::list<RooAbsData*> ret ;
1450 ret.push_back(dat) ;
1451 }
1452 return ret ;
1453}
1454
1455
1456////////////////////////////////////////////////////////////////////////////////
1457/// Return list of all dataset in the workspace
1458
1459std::list<RooAbsData*> RooWorkspace::allEmbeddedData() const
1460{
1461 std::list<RooAbsData*> ret ;
1463 ret.push_back(dat) ;
1464 }
1465 return ret ;
1466}
1467
1468
1469
1470////////////////////////////////////////////////////////////////////////////////
1471/// Return list of all generic objects in the workspace
1472
1473std::list<TObject*> RooWorkspace::allGenericObjects() const
1474{
1475 std::list<TObject*> ret ;
1476 for(TObject * gobj : _genObjects) {
1477
1478 // If found object is wrapper, return payload
1479 if (gobj->IsA()==RooTObjWrap::Class()) {
1480 ret.push_back((static_cast<RooTObjWrap*>(gobj))->obj()) ;
1481 } else {
1482 ret.push_back(gobj) ;
1483 }
1484 }
1485 return ret ;
1486}
1487
1488
1489namespace {
1490
1491std::string findFileInPath(std::string const &file, std::list<std::string> const &dirList)
1492{
1493 // Check list of additional paths
1494 for (std::string const &diter : dirList) {
1495 TString temp = file.c_str();
1496 const char *cpath = gSystem->PrependPathName(diter.c_str(), temp);
1497 std::string path = cpath;
1498 if (!gSystem->AccessPathName(path.c_str())) {
1499 // found file
1500 return path;
1501 }
1502 }
1503 return "";
1504}
1505
1506} // namespace
1507
1508
1509////////////////////////////////////////////////////////////////////////////////
1510/// Import code of class 'tc' into the repository. If code is already in repository it is only imported
1511/// again if doReplace is false. The names and location of the source files is determined from the information
1512/// in TClass. If no location is found in the TClass information, the files are searched in the workspace
1513/// search path, defined by addClassDeclImportDir() and addClassImplImportDir() for declaration and implementation
1514/// files respectively. If files cannot be found, abort with error status, otherwise update the internal
1515/// class-to-file map and import the contents of the files, if they are not imported yet.
1516/// \return true on success, false on failure
1518{
1519
1520 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") request to import code of class " << tc->GetName() << std::endl ;
1521
1522 // *** PHASE 1 *** Check if file needs to be imported, or is in ROOT distribution, and check if it can be persisted
1523
1524 // Check if we already have the class (i.e. it is in the classToFile map)
1525 if (!doReplace && _c2fmap.find(tc->GetName())!=_c2fmap.end()) {
1526 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " already imported, skipping" << std::endl ;
1527 return true ;
1528 }
1529
1530 // Check if class is listed in a ROOTMAP file - if so we can skip it because it is in the root distribution
1531 const char* mapEntry = gInterpreter->GetClassSharedLibs(tc->GetName()) ;
1532 if (mapEntry && strlen(mapEntry)>0) {
1533 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << std::endl ;
1534 return true ;
1535 }
1536
1537 // Retrieve file names through ROOT TClass interface
1538 string implfile = tc->GetImplFileName() ;
1539 string declfile = tc->GetDeclFileName() ;
1540
1541 // Check that file names are not empty
1542 if (implfile.empty() || declfile.empty()) {
1543 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") ERROR: cannot retrieve code file names for class "
1544 << tc->GetName() << " through ROOT TClass interface, unable to import code" << std::endl ;
1545 return false ;
1546 }
1547
1548 // Check if header filename is found in ROOT distribution, if so, do not import class
1549 TString rootsys = gSystem->Getenv("ROOTSYS") ;
1550 if (TString(implfile.c_str()).Index(rootsys)>=0) {
1551 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo(" << _wspace->GetName() << ") code of class " << tc->GetName() << " is in ROOT distribution, skipping " << std::endl ;
1552 return true ;
1553 }
1554
1555 // Require that class meets technical criteria to be persistable (i.e it has a default constructor)
1556 // (We also need a default constructor of abstract classes, but cannot check that through is interface
1557 // as TClass::HasDefaultCtor only returns true for callable default constructors)
1558 if (!(tc->Property() & kIsAbstract) && !tc->HasDefaultConstructor()) {
1559 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING cannot import class "
1560 << tc->GetName() << " : it cannot be persisted because it doesn't have a default constructor. Please fix " << std::endl ;
1561 return false ;
1562 }
1563
1564
1565 // *** PHASE 2 *** Check if declaration and implementation files can be located
1566
1567 std::string declpath;
1568 std::string implpath;
1569
1570 // Check if header file can be found in specified location
1571 // If not, scan through list of 'class declaration' paths in RooWorkspace
1572 if (gSystem->AccessPathName(declfile.c_str())) {
1573
1575
1576 // Header file cannot be found anywhere, warn user and abort operation
1577 if (declpath.empty()) {
1578 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1579 << tc->GetName() << " because header file " << declfile << " is not found in current directory nor in $ROOTSYS" ;
1580 if (!_classDeclDirList.empty()) {
1581 ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1583
1584 while(diter!= RooWorkspace::_classDeclDirList.end()) {
1585
1587 ooccoutW(_wspace,ObjectHandling) << "," ;
1588 }
1589 ooccoutW(_wspace,ObjectHandling) << diter->c_str() ;
1590 ++diter ;
1591 }
1592 }
1593 ooccoutW(_wspace,ObjectHandling) << ". To fix this problem, add the required directory to the search "
1594 << "path using RooWorkspace::addClassDeclImportDir(const char* dir)" << std::endl ;
1595
1596 return false ;
1597 }
1598 }
1599
1600
1601 // Check if implementation file can be found in specified location
1602 // If not, scan through list of 'class implementation' paths in RooWorkspace
1603 if (gSystem->AccessPathName(implfile.c_str())) {
1604
1606
1607 // Implementation file cannot be found anywhere, warn user and abort operation
1608 if (implpath.empty()) {
1609 oocoutW(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName() << ") WARNING Cannot access code of class "
1610 << tc->GetName() << " because implementation file " << implfile << " is not found in current directory nor in $ROOTSYS" ;
1611 if (!_classDeclDirList.empty()) {
1612 ooccoutW(_wspace,ObjectHandling) << ", nor in the search path " ;
1614
1615 while(iiter!= RooWorkspace::_classImplDirList.end()) {
1616
1618 ooccoutW(_wspace,ObjectHandling) << "," ;
1619 }
1620 ooccoutW(_wspace,ObjectHandling) << iiter->c_str() ;
1621 ++iiter ;
1622 }
1623 }
1624 ooccoutW(_wspace,ObjectHandling) << ". To fix this problem add the required directory to the search "
1625 << "path using RooWorkspace::addClassImplImportDir(const char* dir)" << std::endl;
1626 return false;
1627 }
1628 }
1629
1630 char buf[64000];
1631
1632 // *** Phase 3 *** Prepare to import code from files into STL string buffer
1633 //
1634 // Code storage is organized in two linked maps
1635 //
1636 // _fmap contains stl strings with code, indexed on declaration file name
1637 //
1638 // _c2fmap contains list of declaration file names and list of base classes
1639 // and is indexed on class name
1640 //
1641 // Phase 3 is skipped if fmap already contains an entry with given filebasename
1642
1643 const std::string declfilename = !declpath.empty() ? gSystem->BaseName(declpath.c_str())
1644 : gSystem->BaseName(declfile.c_str());
1645
1646 // Split in base and extension
1647 int dotpos2 = strrchr(declfilename.c_str(),'.') - declfilename.c_str() ;
1648 string declfilebase = declfilename.substr(0,dotpos2) ;
1649 string declfileext = declfilename.substr(dotpos2+1) ;
1650
1652
1653 // If file has not been stored yet, enter stl strings with implementation and declaration in file map
1654 if (_fmap.find(declfilebase) == _fmap.end()) {
1655
1656 // Open declaration file
1657 std::fstream fdecl(!declpath.empty() ? declpath.c_str() : declfile.c_str());
1658
1659 // Abort import if declaration file cannot be opened
1660 if (!fdecl) {
1661 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1662 << ") ERROR opening declaration file " << declfile << std::endl ;
1663 return false ;
1664 }
1665
1666 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1667 << ") importing code of class " << tc->GetName()
1668 << " from " << (!implpath.empty() ? implpath.c_str() : implfile.c_str())
1669 << " and " << (!declpath.empty() ? declpath.c_str() : declfile.c_str()) << std::endl ;
1670
1671
1672 // Read entire file into an stl string
1673 string decl ;
1674 while(fdecl.getline(buf,1023)) {
1675
1676 // Look for include state of self
1677 bool processedInclude = false ;
1678 char* extincfile = nullptr ;
1679
1680 // Look for include of declaration file corresponding to this implementation file
1681 if (strstr(buf,"#include")) {
1682 // Process #include statements here
1683 char tmp[64000];
1684 strlcpy(tmp, buf, 64000);
1685 bool stdinclude = strchr(buf, '<');
1686 strtok(tmp, " <\"");
1687 char *incfile = strtok(nullptr, " <>\"");
1688
1689 if (!stdinclude) {
1690 // check if it lives in $ROOTSYS/include
1691 TString hpath = gSystem->Getenv("ROOTSYS");
1692 hpath += "/include/";
1693 hpath += incfile;
1694 if (gSystem->AccessPathName(hpath.Data())) {
1695 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1696 << ") scheduling include file " << incfile << " for import" << std::endl;
1697 extraHeaders.push_back(incfile);
1699 processedInclude = true;
1700 }
1701 }
1702 }
1703
1704 if (processedInclude) {
1705 decl += "// external include file below retrieved from workspace code storage\n" ;
1706 decl += Form("#include \"%s\"\n",extincfile) ;
1707 } else {
1708 decl += buf ;
1709 decl += '\n' ;
1710 }
1711 }
1712
1713 // Open implementation file
1714 fstream fimpl(!implpath.empty() ? implpath.c_str() : implfile.c_str()) ;
1715
1716 // Abort import if implementation file cannot be opened
1717 if (!fimpl) {
1718 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1719 << ") ERROR opening implementation file " << implfile << std::endl ;
1720 return false ;
1721 }
1722
1723
1724 // Import entire implementation file into stl string
1725 string impl ;
1726 while(fimpl.getline(buf,1023)) {
1727 // Process #include statements here
1728
1729 // Look for include state of self
1730 bool foundSelfInclude=false ;
1731 bool processedInclude = false ;
1732 char* extincfile = nullptr ;
1733
1734 // Look for include of declaration file corresponding to this implementation file
1735 if (strstr(buf,"#include")) {
1736 // Process #include statements here
1737 char tmp[64000];
1738 strlcpy(tmp, buf, 64000);
1739 bool stdinclude = strchr(buf, '<');
1740 strtok(tmp, " <\"");
1741 char *incfile = strtok(nullptr, " <>\"");
1742
1743 if (strstr(incfile, declfilename.c_str())) {
1744 foundSelfInclude = true;
1745 }
1746
1747 if (!stdinclude && !foundSelfInclude) {
1748 // check if it lives in $ROOTSYS/include
1749 TString hpath = gSystem->Getenv("ROOTSYS");
1750 hpath += "/include/";
1751 hpath += incfile;
1752
1753 if (gSystem->AccessPathName(hpath.Data())) {
1754 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1755 << ") scheduling include file " << incfile << " for import" << std::endl;
1756 extraHeaders.push_back(incfile);
1758 processedInclude = true;
1759 }
1760 }
1761 }
1762
1763 // Explicitly rewrite include of own declaration file to string
1764 // any directory prefixes, copy all other lines verbatim in stl string
1765 if (foundSelfInclude) {
1766 // If include of self is found, substitute original include
1767 // which may have directory structure with a plain include
1768 impl += "// class declaration include file below retrieved from workspace code storage\n" ;
1769 impl += Form("#include \"%s.%s\"\n",declfilebase.c_str(),declfileext.c_str()) ;
1770 } else if (processedInclude) {
1771 impl += "// external include file below retrieved from workspace code storage\n" ;
1772 impl += Form("#include \"%s\"\n",extincfile) ;
1773 } else {
1774 impl += buf ;
1775 impl += '\n' ;
1776 }
1777 }
1778
1779 // Create entry in file map
1780 _fmap[declfilebase]._hfile = decl ;
1781 _fmap[declfilebase]._cxxfile = impl ;
1782 _fmap[declfilebase]._hext = declfileext ;
1783
1784 // Process extra includes now
1785 for (list<string>::iterator ehiter = extraHeaders.begin() ; ehiter != extraHeaders.end() ; ++ehiter ) {
1786 if (_ehmap.find(*ehiter) == _ehmap.end()) {
1787
1788 ExtraHeader eh ;
1789 eh._hname = ehiter->c_str() ;
1790 fstream fehdr(ehiter->c_str()) ;
1791 string ehimpl ;
1792 char buf2[1024] ;
1793 while(fehdr.getline(buf2,1023)) {
1794
1795 // Look for include of declaration file corresponding to this implementation file
1796 if (strstr(buf2,"#include")) {
1797 // Process #include statements here
1798 char tmp[64000];
1799 strlcpy(tmp, buf2, 64000);
1800 bool stdinclude = strchr(buf, '<');
1801 strtok(tmp, " <\"");
1802 char *incfile = strtok(nullptr, " <>\"");
1803
1804 if (!stdinclude) {
1805 // check if it lives in $ROOTSYS/include
1806 TString hpath = gSystem->Getenv("ROOTSYS");
1807 hpath += "/include/";
1808 hpath += incfile;
1809 if (gSystem->AccessPathName(hpath.Data())) {
1810 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1811 << ") scheduling recursive include file " << incfile << " for import"
1812 << std::endl;
1813 extraHeaders.push_back(incfile);
1814 }
1815 }
1816 }
1817
1818 ehimpl += buf2;
1819 ehimpl += '\n';
1820 }
1821 eh._hfile = ehimpl.c_str();
1822
1823 _ehmap[ehiter->c_str()] = eh;
1824 }
1825 }
1826
1827 } else {
1828
1829 // Inform that existing file entry is being recycled because it already contained class code
1830 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::autoImportClass(" << _wspace->GetName()
1831 << ") code of class " << tc->GetName()
1832 << " was already imported from " << (!implpath.empty() ? implpath : implfile)
1833 << " and " << (!declpath.empty() ? declpath.c_str() : declfile.c_str()) << std::endl;
1834
1835 }
1836
1837
1838 // *** PHASE 4 *** Import stl strings with code into workspace
1839 //
1840 // If multiple classes are declared in a single code unit, there will be
1841 // multiple _c2fmap entries all pointing to the same _fmap entry.
1842
1843 // Make list of all immediate base classes of this class
1845 TList* bl = tc->GetListOfBases() ;
1846 std::list<TClass*> bases ;
1847 for(auto * base : static_range_cast<TBaseClass*>(*bl)) {
1848 if (baseNameList.Length()>0) {
1849 baseNameList += "," ;
1850 }
1851 baseNameList += base->GetClassPointer()->GetName() ;
1852 bases.push_back(base->GetClassPointer()) ;
1853 }
1854
1855 // Map class name to above _fmap entries, along with list of base classes
1856 // in _c2fmap
1857 _c2fmap[tc->GetName()]._baseName = baseNameList ;
1858 _c2fmap[tc->GetName()]._fileBase = declfilebase ;
1859
1860 // Recursive store all base classes.
1861 for(TClass* bclass : bases) {
1863 }
1864
1865 return true ;
1866}
1867
1868
1869////////////////////////////////////////////////////////////////////////////////
1870/// Create transient TDirectory representation of this workspace. This directory
1871/// will appear as a subdirectory of the directory that contains the workspace
1872/// and will have the name of the workspace suffixed with "Dir". The TDirectory
1873/// interface is read-only. Any attempt to insert objects into the workspace
1874/// directory representation will result in an error message. Note that some
1875/// ROOT object like TH1 automatically insert themselves into the current directory
1876/// when constructed. This will give error messages when done in a workspace
1877/// directory.
1878/// \return true (success) always
1880{
1881 if (_dir) return true ;
1882
1883 std::string title= "TDirectory representation of RooWorkspace " + std::string(GetName());
1884 _dir = new WSDir(GetName(),title.c_str(),this) ;
1885
1886 for (RooAbsArg * darg : _allOwnedNodes) {
1887 if (darg->IsA() != RooConstVar::Class()) {
1889 }
1890 }
1891
1892 return true ;
1893}
1894
1895
1896
1897////////////////////////////////////////////////////////////////////////////////
1898/// Import a clone of a generic TObject into workspace generic object container. Imported
1899/// object can be retrieved by name through the obj() method. The object is cloned upon
1900/// importation and the input argument does not need to live beyond the import call
1901///
1902/// \return due to historical reasons: false (0) on success and true (1) on failure
1903
1905{
1906 // First check if object with given name already exists
1907 std::unique_ptr<TObject> oldObj{_genObjects.FindObject(object.GetName())};
1908 if (oldObj && !replaceExisting) {
1909 coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
1910 << object.GetName() << " is already in workspace and replaceExisting flag is set to false" << std::endl ;
1911 return true ;
1912 }
1913
1914 // Grab the current state of the directory Auto-Add
1915 ROOT::DirAutoAdd_t func = object.IsA()->GetDirectoryAutoAdd();
1916 object.IsA()->SetDirectoryAutoAdd(nullptr);
1917 bool tmp = RooPlot::setAddDirectoryStatus(false) ;
1918
1919 if (oldObj) {
1920 _genObjects.Replace(oldObj.get(),object.Clone()) ;
1921 } else {
1922 _genObjects.Add(object.Clone()) ;
1923 }
1924
1925 // Reset the state of the directory Auto-Add
1926 object.IsA()->SetDirectoryAutoAdd(func);
1928
1929 return false ;
1930}
1931
1932
1933
1934
1935////////////////////////////////////////////////////////////////////////////////
1936/// Import a clone of a generic TObject into workspace generic object container.
1937/// The imported object will be stored under the given alias name rather than its
1938/// own name. Imported object can be retrieved its alias name through the obj() method.
1939/// The object is cloned upon importation and the input argument does not need to live beyond the import call
1940/// This method is mostly useful for importing objects that do not have a settable name such as TMatrix
1941///
1942/// \return due to historical reasons: false (0) on success and true (1) on failure
1943
1944bool RooWorkspace::import(TObject const& object, const char* aliasName, bool replaceExisting)
1945{
1946 // First check if object with given name already exists
1947 std::unique_ptr<TObject> oldObj{_genObjects.FindObject(aliasName)};
1948 if (oldObj && !replaceExisting) {
1949 coutE(InputArguments) << "RooWorkspace::import(" << GetName() << ") generic object with name "
1950 << aliasName << " is already in workspace and replaceExisting flag is set to false" << std::endl ;
1951 return true ;
1952 }
1953
1954 TH1::AddDirectory(false) ;
1955 auto wrapper = new RooTObjWrap(object.Clone()) ;
1956 TH1::AddDirectory(true) ;
1957 wrapper->setOwning(true) ;
1958 wrapper->SetName(aliasName) ;
1959 wrapper->SetTitle(aliasName) ;
1960
1961 if (oldObj) {
1963 } else {
1965 }
1966 return false ;
1967}
1968
1969
1970
1971
1972////////////////////////////////////////////////////////////////////////////////
1973/// Insert RooStudyManager module
1974/// \return due to historical reasons: false (0) on success and true (1) on failure
1976{
1977 RooAbsStudy* clone = static_cast<RooAbsStudy*>(study.Clone()) ;
1978 _studyMods.Add(clone) ;
1979 return false ;
1980}
1981
1982
1983
1984
1985////////////////////////////////////////////////////////////////////////////////
1986/// Remove all RooStudyManager modules
1987
1992
1993
1994
1995
1996////////////////////////////////////////////////////////////////////////////////
1997/// Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
1998
2000{
2001 // Try RooAbsArg first
2002 TObject* ret = arg(name) ;
2003 if (ret) return ret ;
2004
2005 // Then try RooAbsData
2006 ret = data(name) ;
2007 if (ret) return ret ;
2008
2009 // Finally try generic object store
2010 return genobj(name) ;
2011}
2012
2013
2014
2015////////////////////////////////////////////////////////////////////////////////
2016/// Return generic object with given name
2017
2019{
2020 // Find object by name
2021 TObject* gobj = _genObjects.FindObject(name.c_str()) ;
2022
2023 // Exit here if not found
2024 if (!gobj) return nullptr;
2025
2026 // If found object is wrapper, return payload
2027 if (gobj->IsA()==RooTObjWrap::Class()) return (static_cast<RooTObjWrap*>(gobj))->obj() ;
2028
2029 return gobj ;
2030}
2031
2032
2033
2034////////////////////////////////////////////////////////////////////////////////
2035/// \return true on success, false on failure
2036/// \see TDirectoryFile::cd
2037bool RooWorkspace::cd(const char* path)
2038{
2039 makeDir() ;
2040 return _dir->cd(path) ;
2041}
2042
2043
2044
2045////////////////////////////////////////////////////////////////////////////////
2046/// Save this current workspace into given file
2047/// \return true if file correctly written, false in case of error
2048
2049bool RooWorkspace::writeToFile(const char* fileName, bool recreate)
2050{
2051 std::unique_ptr<TFile> f{ TFile::Open(fileName, recreate ? "RECREATE" : "UPDATE") };
2052 if (!f || f->IsZombie())
2053 return false;
2054 auto bytes = Write();
2055 return bytes > 0;
2056}
2057
2058
2059
2060////////////////////////////////////////////////////////////////////////////////
2061/// Return instance to factory tool
2062
2064{
2065 if (_factory) {
2066 return *_factory;
2067 }
2068 cxcoutD(ObjectHandling) << "INFO: Creating RooFactoryWSTool associated with this workspace" << std::endl ;
2070 return *_factory;
2071}
2072
2073
2074
2075
2076////////////////////////////////////////////////////////////////////////////////
2077/// Short-hand function for `factory()->process(expr);`
2078///
2079/// \copydoc RooFactoryWSTool::process(const char*)
2084
2085
2086
2087
2088////////////////////////////////////////////////////////////////////////////////
2089/// Print contents of the workspace
2090
2092{
2093 bool treeMode(false) ;
2094 bool verbose(false);
2095 if (TString(opts).Contains("t")) {
2096 treeMode=true ;
2097 }
2098 if (TString(opts).Contains("v")) {
2099 verbose = true;
2100 }
2101
2102 std::cout << std::endl << "RooWorkspace(" << GetName() << ") " << GetTitle() << " contents" << std::endl << std::endl ;
2103
2106 RooArgSet varSet ;
2110
2111
2112 // Split list of components in pdfs, functions and variables
2113 for(RooAbsArg* parg : _allOwnedNodes) {
2114
2115 //---------------
2116
2117 if (treeMode) {
2118
2119 // In tree mode, only add nodes with no clients to the print lists
2120
2121 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class())) {
2122 if (!parg->hasClients()) {
2123 pdfSet.add(*parg) ;
2124 }
2125 }
2126
2127 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2128 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2129 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2130 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2131 if (!parg->hasClients()) {
2132 funcSet.add(*parg) ;
2133 }
2134 }
2135
2136
2137 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2138 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2139 if (!parg->hasClients()) {
2140 catfuncSet.add(*parg) ;
2141 }
2142 }
2143
2144 } else {
2145
2146 if (parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2147 if ((static_cast<RooResolutionModel*>(parg))->isConvolved()) {
2148 convResoSet.add(*parg) ;
2149 } else {
2150 resoSet.add(*parg) ;
2151 }
2152 }
2153
2154 if (parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2155 !parg->IsA()->InheritsFrom(RooResolutionModel::Class())) {
2156 pdfSet.add(*parg) ;
2157 }
2158
2159 if (parg->IsA()->InheritsFrom(RooAbsReal::Class()) &&
2160 !parg->IsA()->InheritsFrom(RooAbsPdf::Class()) &&
2161 !parg->IsA()->InheritsFrom(RooConstVar::Class()) &&
2162 !parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2163 funcSet.add(*parg) ;
2164 }
2165
2166 if (parg->IsA()->InheritsFrom(RooAbsCategory::Class()) &&
2167 !parg->IsA()->InheritsFrom(RooCategory::Class())) {
2168 catfuncSet.add(*parg) ;
2169 }
2170 }
2171
2172 if (parg->IsA()->InheritsFrom(RooRealVar::Class())) {
2173 varSet.add(*parg) ;
2174 }
2175
2176 if (parg->IsA()->InheritsFrom(RooCategory::Class())) {
2177 varSet.add(*parg) ;
2178 }
2179
2180 }
2181
2182
2183 RooFit::MsgLevel oldLevel = RooMsgService::instance().globalKillBelow() ;
2184 RooMsgService::instance().setGlobalKillBelow(RooFit::WARNING) ;
2185
2186 if (!varSet.empty()) {
2187 varSet.sort() ;
2188 std::cout << "variables" << std::endl ;
2189 std::cout << "---------" << std::endl ;
2190 std::cout << varSet << std::endl ;
2191 std::cout << std::endl ;
2192 }
2193
2194 if (!pdfSet.empty()) {
2195 std::cout << "p.d.f.s" << std::endl ;
2196 std::cout << "-------" << std::endl ;
2197 pdfSet.sort() ;
2198 for(RooAbsArg* parg : pdfSet) {
2199 if (treeMode) {
2200 parg->printComponentTree() ;
2201 } else {
2202 parg->Print() ;
2203 }
2204 }
2205 std::cout << std::endl ;
2206 }
2207
2208 if (!treeMode) {
2209 if (!resoSet.empty()) {
2210 std::cout << "analytical resolution models" << std::endl ;
2211 std::cout << "----------------------------" << std::endl ;
2212 resoSet.sort() ;
2213 for(RooAbsArg* parg : resoSet) {
2214 parg->Print() ;
2215 }
2216 std::cout << std::endl ;
2217 }
2218 }
2219
2220 if (!funcSet.empty()) {
2221 std::cout << "functions" << std::endl ;
2222 std::cout << "--------" << std::endl ;
2223 funcSet.sort() ;
2224 for(RooAbsArg * parg : funcSet) {
2225 if (treeMode) {
2226 parg->printComponentTree() ;
2227 } else {
2228 parg->Print() ;
2229 }
2230 }
2231 std::cout << std::endl ;
2232 }
2233
2234 if (!catfuncSet.empty()) {
2235 std::cout << "category functions" << std::endl ;
2236 std::cout << "------------------" << std::endl ;
2237 catfuncSet.sort() ;
2238 for(RooAbsArg* parg : catfuncSet) {
2239 if (treeMode) {
2240 parg->printComponentTree() ;
2241 } else {
2242 parg->Print() ;
2243 }
2244 }
2245 std::cout << std::endl ;
2246 }
2247
2248 if (!_dataList.empty()) {
2249 std::cout << "datasets" << std::endl ;
2250 std::cout << "--------" << std::endl ;
2252 std::cout << data2->ClassName() << "::" << data2->GetName() << *data2->get() << std::endl;
2253 }
2254 std::cout << std::endl ;
2255 }
2256
2257 if (!_embeddedDataList.empty()) {
2258 std::cout << "embedded datasets (in pdfs and functions)" << std::endl ;
2259 std::cout << "-----------------------------------------" << std::endl ;
2261 std::cout << data2->ClassName() << "::" << data2->GetName() << *data2->get() << std::endl ;
2262 }
2263 std::cout << std::endl ;
2264 }
2265
2266 if (!_snapshots.empty()) {
2267 std::cout << "parameter snapshots" << std::endl ;
2268 std::cout << "-------------------" << std::endl ;
2270 std::cout << snap->GetName() << " = (" ;
2271 bool first(true) ;
2272 for(RooAbsArg* a : *snap) {
2273 if (first) { first=false ; } else { std::cout << "," ; }
2274 std::cout << a->GetName() << "=" ;
2275 a->printValue(std::cout) ;
2276 if (a->isConstant()) {
2277 std::cout << "[C]" ;
2278 }
2279 }
2280 std::cout << ")" << std::endl ;
2281 }
2282 std::cout << std::endl ;
2283 }
2284
2285
2286 if (!_namedSets.empty()) {
2287 std::cout << "named sets" << std::endl ;
2288 std::cout << "----------" << std::endl ;
2289 for (map<string,RooArgSet>::const_iterator it = _namedSets.begin() ; it != _namedSets.end() ; ++it) {
2290 if (verbose || !isCacheSet(it->first)) {
2291 std::cout << it->first << ":" << it->second << std::endl;
2292 }
2293 }
2294
2295 std::cout << std::endl ;
2296 }
2297
2298
2299 if (!_genObjects.empty()) {
2300 std::cout << "generic objects" << std::endl ;
2301 std::cout << "---------------" << std::endl ;
2302 for(TObject* gobj : _genObjects) {
2303 if (gobj->IsA()==RooTObjWrap::Class()) {
2304 std::cout << (static_cast<RooTObjWrap*>(gobj))->obj()->ClassName() << "::" << gobj->GetName() << std::endl ;
2305 } else {
2306 std::cout << gobj->ClassName() << "::" << gobj->GetName() << std::endl ;
2307 }
2308 }
2309 std::cout << std::endl ;
2310
2311 }
2312
2313 if (!_studyMods.empty()) {
2314 std::cout << "study modules" << std::endl ;
2315 std::cout << "-------------" << std::endl ;
2316 for(TObject* smobj : _studyMods) {
2317 std::cout << smobj->ClassName() << "::" << smobj->GetName() << std::endl ;
2318 }
2319 std::cout << std::endl ;
2320
2321 }
2322
2323 if (!_classes.listOfClassNames().empty()) {
2324 std::cout << "embedded class code" << std::endl ;
2325 std::cout << "-------------------" << std::endl ;
2326 std::cout << _classes.listOfClassNames() << std::endl ;
2327 std::cout << std::endl ;
2328 }
2329
2330 if (!_eocache.empty()) {
2331 std::cout << "embedded precalculated expensive components" << std::endl ;
2332 std::cout << "-------------------------------------------" << std::endl ;
2333 _eocache.print() ;
2334 }
2335
2336 RooMsgService::instance().setGlobalKillBelow(oldLevel) ;
2337
2338 return ;
2339}
2340
2341
2342////////////////////////////////////////////////////////////////////////////////
2343/// Custom streamer for the workspace. Stream contents of workspace
2344/// and code repository. When reading, read code repository first
2345/// and compile missing classes before proceeding with streaming
2346/// of workspace contents to avoid errors.
2347
2349{
2351
2352 // Stream an object of class RooWorkspace::CodeRepo.
2353 if (R__b.IsReading()) {
2354
2355 UInt_t R__s;
2356 UInt_t R__c;
2357 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2358
2359 // Stream contents of ClassFiles map
2360 Int_t count(0);
2361 R__b >> count;
2362 while (count--) {
2363 TString name;
2364 name.Streamer(R__b);
2365 _fmap[name]._hext.Streamer(R__b);
2366 _fmap[name]._hfile.Streamer(R__b);
2367 _fmap[name]._cxxfile.Streamer(R__b);
2368 }
2369
2370 // Stream contents of ClassRelInfo map
2371 count = 0;
2372 R__b >> count;
2373 while (count--) {
2374 TString name;
2375 name.Streamer(R__b);
2376 _c2fmap[name]._baseName.Streamer(R__b);
2377 _c2fmap[name]._fileBase.Streamer(R__b);
2378 }
2379
2380 if (R__v == 2) {
2381
2382 count = 0;
2383 R__b >> count;
2384 while (count--) {
2385 TString name;
2386 name.Streamer(R__b);
2387 _ehmap[name]._hname.Streamer(R__b);
2388 _ehmap[name]._hfile.Streamer(R__b);
2389 }
2390 }
2391
2392 R__b.CheckByteCount(R__s, R__c, thisClass::IsA());
2393
2394 // Instantiate any classes that are not defined in current session
2395 _compiledOK = !compileClasses();
2396
2397 } else {
2398
2399 UInt_t R__c;
2400 R__c = R__b.WriteVersion(thisClass::IsA(), true);
2401
2402 // Stream contents of ClassFiles map
2403 UInt_t count = _fmap.size();
2404 R__b << count;
2405 map<TString, ClassFiles>::iterator iter = _fmap.begin();
2406 while (iter != _fmap.end()) {
2407 TString key_copy(iter->first);
2408 key_copy.Streamer(R__b);
2409 iter->second._hext.Streamer(R__b);
2410 iter->second._hfile.Streamer(R__b);
2411 iter->second._cxxfile.Streamer(R__b);
2412
2413 ++iter;
2414 }
2415
2416 // Stream contents of ClassRelInfo map
2417 count = _c2fmap.size();
2418 R__b << count;
2419 map<TString, ClassRelInfo>::iterator iter2 = _c2fmap.begin();
2420 while (iter2 != _c2fmap.end()) {
2421 TString key_copy(iter2->first);
2422 key_copy.Streamer(R__b);
2423 iter2->second._baseName.Streamer(R__b);
2424 iter2->second._fileBase.Streamer(R__b);
2425 ++iter2;
2426 }
2427
2428 // Stream contents of ExtraHeader map
2429 count = _ehmap.size();
2430 R__b << count;
2431 map<TString, ExtraHeader>::iterator iter3 = _ehmap.begin();
2432 while (iter3 != _ehmap.end()) {
2433 TString key_copy(iter3->first);
2434 key_copy.Streamer(R__b);
2435 iter3->second._hname.Streamer(R__b);
2436 iter3->second._hfile.Streamer(R__b);
2437 ++iter3;
2438 }
2439
2440 R__b.SetByteCount(R__c, true);
2441 }
2442}
2443
2444
2445////////////////////////////////////////////////////////////////////////////////
2446/// Stream an object of class RooWorkspace. This is a standard ROOT streamer for the
2447/// I/O part. This custom function exists to detach all external client links
2448/// from the payload prior to writing the payload so that these client links
2449/// are not persisted. (Client links occur if external function objects use
2450/// objects contained in the workspace as input)
2451/// After the actual writing, these client links are restored.
2452
2454{
2455 if (R__b.IsReading()) {
2456
2457 R__b.ReadClassBuffer(RooWorkspace::Class(), this);
2458
2459 // Perform any pass-2 schema evolution here
2460 for (RooAbsArg *node : _allOwnedNodes) {
2461 node->ioStreamerPass2();
2462 }
2464
2465 // Make expensive object cache of all objects point to intermal copy.
2466 // Somehow this doesn't work OK automatically
2467 for (RooAbsArg *node : _allOwnedNodes) {
2468 node->setExpensiveObjectCache(_eocache);
2469 node->setWorkspace(*this);
2470#ifdef ROOFIT_LEGACY_EVAL_BACKEND
2471 if (dynamic_cast<RooAbsOptTestStatistic *>(node)) {
2472 RooAbsOptTestStatistic *tmp = static_cast<RooAbsOptTestStatistic *>(node);
2473 if (tmp->isSealed() && tmp->sealNotice() && strlen(tmp->sealNotice()) > 0) {
2474 std::cout << "RooWorkspace::Streamer(" << GetName() << ") " << node->ClassName() << "::" << node->GetName()
2475 << " : " << tmp->sealNotice() << std::endl;
2476 }
2477 }
2478#endif
2479 }
2480
2481 for(TObject * gobj : allGenericObjects()) {
2482 if (auto handle = dynamic_cast<RooWorkspaceHandle*>(gobj)) {
2483 handle->ReplaceWS(this);
2484 }
2485 }
2486
2487 } else {
2488
2489 // Make lists of external clients of WS objects, and remove those links temporarily
2490
2494
2496
2497 // Loop over client list of this arg
2498 std::vector<RooAbsArg *> clientsTmp{tmparg->_clientList.begin(), tmparg->_clientList.end()};
2499 for (auto client : clientsTmp) {
2500 if (!_allOwnedNodes.containsInstance(*client)) {
2501
2502 const auto refCount = tmparg->_clientList.refCount(client);
2503 auto &bufferVec = extClients[tmparg];
2504
2505 bufferVec.insert(bufferVec.end(), refCount, client);
2506 tmparg->_clientList.Remove(client, true);
2507 }
2508 }
2509
2510 // Loop over value client list of this arg
2511 clientsTmp.assign(tmparg->_clientListValue.begin(), tmparg->_clientListValue.end());
2512 for (auto vclient : clientsTmp) {
2514 cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2515 << " has external value client link to " << vclient << " (" << vclient->GetName()
2516 << ") with ref count " << tmparg->_clientListValue.refCount(vclient) << std::endl;
2517
2518 const auto refCount = tmparg->_clientListValue.refCount(vclient);
2520
2521 bufferVec.insert(bufferVec.end(), refCount, vclient);
2522 tmparg->_clientListValue.Remove(vclient, true);
2523 }
2524 }
2525
2526 // Loop over shape client list of this arg
2527 clientsTmp.assign(tmparg->_clientListShape.begin(), tmparg->_clientListShape.end());
2528 for (auto sclient : clientsTmp) {
2530 cxcoutD(ObjectHandling) << "RooWorkspace::Streamer(" << GetName() << ") element " << tmparg->GetName()
2531 << " has external shape client link to " << sclient << " (" << sclient->GetName()
2532 << ") with ref count " << tmparg->_clientListShape.refCount(sclient) << std::endl;
2533
2534 const auto refCount = tmparg->_clientListShape.refCount(sclient);
2536
2537 bufferVec.insert(bufferVec.end(), refCount, sclient);
2538 tmparg->_clientListShape.Remove(sclient, true);
2539 }
2540 }
2541 }
2542
2543 R__b.WriteClassBuffer(RooWorkspace::Class(), this);
2544
2545 // Reinstate clients here
2546
2547 for (auto &iterx : extClients) {
2548 for (auto client : iterx.second) {
2549 iterx.first->_clientList.Add(client);
2550 }
2551 }
2552
2553 for (auto &iterx : extValueClients) {
2554 for (auto client : iterx.second) {
2555 iterx.first->_clientListValue.Add(client);
2556 }
2557 }
2558
2559 for (auto &iterx : extShapeClients) {
2560 for (auto client : iterx.second) {
2561 iterx.first->_clientListShape.Add(client);
2562 }
2563 }
2564 }
2565}
2566
2567
2568
2569
2570////////////////////////////////////////////////////////////////////////////////
2571/// Return STL string with last of class names contained in the code repository
2572
2574{
2575 string ret ;
2576 map<TString,ClassRelInfo>::const_iterator iter = _c2fmap.begin() ;
2577 while(iter!=_c2fmap.end()) {
2578 if (!ret.empty()) {
2579 ret += ", " ;
2580 }
2581 ret += iter->first ;
2582 ++iter ;
2583 }
2584
2585 return ret ;
2586}
2587
2588namespace {
2589UInt_t crc32(const char* data, ULong_t sz, UInt_t crc)
2590{
2591 // update CRC32 with new data
2592
2593 // use precomputed table, rather than computing it on the fly
2594 static const UInt_t crctab[256] = { 0x00000000,
2595 0x04c11db7, 0x09823b6e, 0x0d4326d9, 0x130476dc, 0x17c56b6b,
2596 0x1a864db2, 0x1e475005, 0x2608edb8, 0x22c9f00f, 0x2f8ad6d6,
2597 0x2b4bcb61, 0x350c9b64, 0x31cd86d3, 0x3c8ea00a, 0x384fbdbd,
2598 0x4c11db70, 0x48d0c6c7, 0x4593e01e, 0x4152fda9, 0x5f15adac,
2599 0x5bd4b01b, 0x569796c2, 0x52568b75, 0x6a1936c8, 0x6ed82b7f,
2600 0x639b0da6, 0x675a1011, 0x791d4014, 0x7ddc5da3, 0x709f7b7a,
2601 0x745e66cd, 0x9823b6e0, 0x9ce2ab57, 0x91a18d8e, 0x95609039,
2602 0x8b27c03c, 0x8fe6dd8b, 0x82a5fb52, 0x8664e6e5, 0xbe2b5b58,
2603 0xbaea46ef, 0xb7a96036, 0xb3687d81, 0xad2f2d84, 0xa9ee3033,
2604 0xa4ad16ea, 0xa06c0b5d, 0xd4326d90, 0xd0f37027, 0xddb056fe,
2605 0xd9714b49, 0xc7361b4c, 0xc3f706fb, 0xceb42022, 0xca753d95,
2606 0xf23a8028, 0xf6fb9d9f, 0xfbb8bb46, 0xff79a6f1, 0xe13ef6f4,
2607 0xe5ffeb43, 0xe8bccd9a, 0xec7dd02d, 0x34867077, 0x30476dc0,
2608 0x3d044b19, 0x39c556ae, 0x278206ab, 0x23431b1c, 0x2e003dc5,
2609 0x2ac12072, 0x128e9dcf, 0x164f8078, 0x1b0ca6a1, 0x1fcdbb16,
2610 0x018aeb13, 0x054bf6a4, 0x0808d07d, 0x0cc9cdca, 0x7897ab07,
2611 0x7c56b6b0, 0x71159069, 0x75d48dde, 0x6b93dddb, 0x6f52c06c,
2612 0x6211e6b5, 0x66d0fb02, 0x5e9f46bf, 0x5a5e5b08, 0x571d7dd1,
2613 0x53dc6066, 0x4d9b3063, 0x495a2dd4, 0x44190b0d, 0x40d816ba,
2614 0xaca5c697, 0xa864db20, 0xa527fdf9, 0xa1e6e04e, 0xbfa1b04b,
2615 0xbb60adfc, 0xb6238b25, 0xb2e29692, 0x8aad2b2f, 0x8e6c3698,
2616 0x832f1041, 0x87ee0df6, 0x99a95df3, 0x9d684044, 0x902b669d,
2617 0x94ea7b2a, 0xe0b41de7, 0xe4750050, 0xe9362689, 0xedf73b3e,
2618 0xf3b06b3b, 0xf771768c, 0xfa325055, 0xfef34de2, 0xc6bcf05f,
2619 0xc27dede8, 0xcf3ecb31, 0xcbffd686, 0xd5b88683, 0xd1799b34,
2620 0xdc3abded, 0xd8fba05a, 0x690ce0ee, 0x6dcdfd59, 0x608edb80,
2621 0x644fc637, 0x7a089632, 0x7ec98b85, 0x738aad5c, 0x774bb0eb,
2622 0x4f040d56, 0x4bc510e1, 0x46863638, 0x42472b8f, 0x5c007b8a,
2623 0x58c1663d, 0x558240e4, 0x51435d53, 0x251d3b9e, 0x21dc2629,
2624 0x2c9f00f0, 0x285e1d47, 0x36194d42, 0x32d850f5, 0x3f9b762c,
2625 0x3b5a6b9b, 0x0315d626, 0x07d4cb91, 0x0a97ed48, 0x0e56f0ff,
2626 0x1011a0fa, 0x14d0bd4d, 0x19939b94, 0x1d528623, 0xf12f560e,
2627 0xf5ee4bb9, 0xf8ad6d60, 0xfc6c70d7, 0xe22b20d2, 0xe6ea3d65,
2628 0xeba91bbc, 0xef68060b, 0xd727bbb6, 0xd3e6a601, 0xdea580d8,
2629 0xda649d6f, 0xc423cd6a, 0xc0e2d0dd, 0xcda1f604, 0xc960ebb3,
2630 0xbd3e8d7e, 0xb9ff90c9, 0xb4bcb610, 0xb07daba7, 0xae3afba2,
2631 0xaafbe615, 0xa7b8c0cc, 0xa379dd7b, 0x9b3660c6, 0x9ff77d71,
2632 0x92b45ba8, 0x9675461f, 0x8832161a, 0x8cf30bad, 0x81b02d74,
2633 0x857130c3, 0x5d8a9099, 0x594b8d2e, 0x5408abf7, 0x50c9b640,
2634 0x4e8ee645, 0x4a4ffbf2, 0x470cdd2b, 0x43cdc09c, 0x7b827d21,
2635 0x7f436096, 0x7200464f, 0x76c15bf8, 0x68860bfd, 0x6c47164a,
2636 0x61043093, 0x65c52d24, 0x119b4be9, 0x155a565e, 0x18197087,
2637 0x1cd86d30, 0x029f3d35, 0x065e2082, 0x0b1d065b, 0x0fdc1bec,
2638 0x3793a651, 0x3352bbe6, 0x3e119d3f, 0x3ad08088, 0x2497d08d,
2639 0x2056cd3a, 0x2d15ebe3, 0x29d4f654, 0xc5a92679, 0xc1683bce,
2640 0xcc2b1d17, 0xc8ea00a0, 0xd6ad50a5, 0xd26c4d12, 0xdf2f6bcb,
2641 0xdbee767c, 0xe3a1cbc1, 0xe760d676, 0xea23f0af, 0xeee2ed18,
2642 0xf0a5bd1d, 0xf464a0aa, 0xf9278673, 0xfde69bc4, 0x89b8fd09,
2643 0x8d79e0be, 0x803ac667, 0x84fbdbd0, 0x9abc8bd5, 0x9e7d9662,
2644 0x933eb0bb, 0x97ffad0c, 0xafb010b1, 0xab710d06, 0xa6322bdf,
2645 0xa2f33668, 0xbcb4666d, 0xb8757bda, 0xb5365d03, 0xb1f740b4
2646 };
2647
2648 crc = ~crc;
2649 while (sz--) crc = (crc << 8) ^ UInt_t(*data++) ^ crctab[crc >> 24];
2650
2651 return ~crc;
2652}
2653
2654UInt_t crc32(const char* data)
2655{
2656 // Calculate crc32 checksum on given string
2657 unsigned long sz = strlen(data);
2658 switch (strlen(data)) {
2659 case 0:
2660 return 0;
2661 case 1:
2662 return data[0];
2663 case 2:
2664 return (data[0] << 8) | data[1];
2665 case 3:
2666 return (data[0] << 16) | (data[1] << 8) | data[2];
2667 case 4:
2668 return (data[0] << 24) | (data[1] << 16) | (data[2] << 8) | data[3];
2669 default:
2670 return crc32(data + 4, sz - 4, (data[0] << 24) | (data[1] << 16) |
2671 (data[2] << 8) | data[3]);
2672 }
2673}
2674
2675}
2676
2677////////////////////////////////////////////////////////////////////////////////
2678/// For all classes in the workspace for which no class definition is
2679/// found in the ROOT class table extract source code stored in code
2680/// repository into temporary directory set by
2681/// setClassFileExportDir(), compile classes and link them with
2682/// current ROOT session. If a compilation error occurs print
2683/// instructions for user how to fix errors and recover workspace and
2684/// abort import procedure.
2685/// \return true on success, false on failure
2686
2688{
2689 bool haveDir=false ;
2690
2691 // Retrieve name of directory in which to export code files
2692 string dirName = Form(_classFileExportDir.c_str(),_wspace->uuid().AsString(),_wspace->GetName()) ;
2693
2694 bool writeExtraHeaders(false) ;
2695
2696 // Process all class entries in repository
2697 map<TString,ClassRelInfo>::iterator iter = _c2fmap.begin() ;
2698 while(iter!=_c2fmap.end()) {
2699
2700 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing class " << iter->first.Data() << std::endl ;
2701
2702 // If class is already known, don't load
2703 if (gClassTable->GetDict(iter->first.Data())) {
2704 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Embedded class "
2705 << iter->first << " already in ROOT class table, skipping" << std::endl ;
2706 ++iter ;
2707 continue ;
2708 }
2709
2710 // Check that export directory exists
2711 if (!haveDir) {
2712
2713 // If not, make local directory to extract files
2714 if (!gSystem->AccessPathName(dirName.c_str())) {
2715 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() reusing code export directory " << dirName.c_str()
2716 << " to extract coded embedded in workspace" << std::endl ;
2717 } else {
2718 if (gSystem->MakeDirectory(dirName.c_str())==0) {
2719 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() creating code export directory " << dirName.c_str()
2720 << " to extract coded embedded in workspace" << std::endl ;
2721 } else {
2722 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR creating code export directory " << dirName.c_str()
2723 << " to extract coded embedded in workspace" << std::endl ;
2724 return false ;
2725 }
2726 }
2727 haveDir=true ;
2728
2729 }
2730
2731 // First write any extra header files
2732 if (!writeExtraHeaders) {
2734
2735 map<TString,ExtraHeader>::iterator extraIter = _ehmap.begin() ;
2736 while(extraIter!=_ehmap.end()) {
2737
2738 // Check if identical declaration file (header) is already written
2739 bool needEHWrite=true ;
2740 string fdname = Form("%s/%s",dirName.c_str(),extraIter->second._hname.Data()) ;
2741 ifstream ifdecl(fdname.c_str()) ;
2742 if (ifdecl) {
2743 TString contents ;
2744 char buf[64000];
2745 while (ifdecl.getline(buf, 64000)) {
2746 contents += buf;
2747 contents += "\n";
2748 }
2749 UInt_t crcFile = crc32(contents.Data());
2750 UInt_t crcWS = crc32(extraIter->second._hfile.Data());
2751 needEHWrite = (crcFile != crcWS);
2752 }
2753
2754 // Write declaration file if required
2755 if (needEHWrite) {
2756 oocoutI(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting extra header file "
2757 << fdname << std::endl;
2758
2759 // Extra headers may contain non-existing path - create first to be sure
2761
2762 ofstream fdecl(fdname.c_str());
2763 if (!fdecl) {
2764 oocoutE(_wspace, ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file " << fdname
2765 << " for writing" << std::endl;
2766 return false;
2767 }
2768 fdecl << extraIter->second._hfile.Data();
2769 fdecl.close();
2770 }
2771 ++extraIter;
2772 }
2773 }
2774
2775
2776 // Navigate from class to file
2777 ClassFiles& cfinfo = _fmap[iter->second._fileBase] ;
2778
2779 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() now processing file with base " << iter->second._fileBase << std::endl ;
2780
2781 // If file is already processed, skip to next class
2782 if (cfinfo._extracted) {
2783 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() file with base name " << iter->second._fileBase
2784 << " has already been extracted, skipping to next class" << std::endl ;
2785 continue ;
2786 }
2787
2788 // Check if identical declaration file (header) is already written
2789 bool needDeclWrite=true ;
2790 string fdname = Form("%s/%s.%s",dirName.c_str(),iter->second._fileBase.Data(),cfinfo._hext.Data()) ;
2791 ifstream ifdecl(fdname.c_str()) ;
2792 if (ifdecl) {
2793 TString contents ;
2794 char buf[64000];
2795 while (ifdecl.getline(buf, 64000)) {
2796 contents += buf;
2797 contents += "\n";
2798 }
2799 UInt_t crcFile = crc32(contents.Data()) ;
2800 UInt_t crcWS = crc32(cfinfo._hfile.Data()) ;
2802 }
2803
2804 // Write declaration file if required
2805 if (needDeclWrite) {
2806 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting declaration code of class " << iter->first << ", file " << fdname << std::endl ;
2807 ofstream fdecl(fdname.c_str()) ;
2808 if (!fdecl) {
2809 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file "
2810 << fdname << " for writing" << std::endl ;
2811 return false ;
2812 }
2813 fdecl << cfinfo._hfile ;
2814 fdecl.close() ;
2815 }
2816
2817 // Check if identical implementation file is already written
2818 bool needImplWrite=true ;
2819 string finame = Form("%s/%s.cxx",dirName.c_str(),iter->second._fileBase.Data()) ;
2820 ifstream ifimpl(finame.c_str()) ;
2821 if (ifimpl) {
2822 TString contents ;
2823 char buf[64000];
2824 while (ifimpl.getline(buf, 64000)) {
2825 contents += buf;
2826 contents += "\n";
2827 }
2828 UInt_t crcFile = crc32(contents.Data()) ;
2829 UInt_t crcWS = crc32(cfinfo._cxxfile.Data()) ;
2831 }
2832
2833 // Write implementation file if required
2834 if (needImplWrite) {
2835 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Extracting implementation code of class " << iter->first << ", file " << finame << std::endl ;
2836 ofstream fimpl(finame.c_str()) ;
2837 if (!fimpl) {
2838 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR opening file"
2839 << finame << " for writing" << std::endl ;
2840 return false ;
2841 }
2842 fimpl << cfinfo._cxxfile ;
2843 fimpl.close() ;
2844 }
2845
2846 // Mark this file as extracted
2847 cfinfo._extracted = true ;
2848 oocxcoutD(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() marking code unit " << iter->second._fileBase << " as extracted" << std::endl ;
2849
2850 // Compile class
2851 oocoutI(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() Compiling code unit " << iter->second._fileBase.Data() << " to define class " << iter->first << std::endl ;
2852 bool ok = gSystem->CompileMacro(finame.c_str(),"k") ;
2853
2854 if (!ok) {
2855 oocoutE(_wspace,ObjectHandling) << "RooWorkspace::CodeRepo::compileClasses() ERROR compiling class " << iter->first.Data() << ", to fix this you can do the following: " << std::endl
2856 << " 1) Fix extracted source code files in directory " << dirName.c_str() << "/" << std::endl
2857 << " 2) In clean ROOT session compiled fixed classes by hand using '.x " << dirName.c_str() << "/ClassName.cxx+'" << std::endl
2858 << " 3) Reopen file with RooWorkspace with broken source code in UPDATE mode. Access RooWorkspace to force loading of class" << std::endl
2859 << " Broken instances in workspace will _not_ be compiled, instead precompiled fixed instances will be used." << std::endl
2860 << " 4) Reimport fixed code in workspace using 'RooWorkspace::importClassCode(\"*\",true)' method, Write() updated workspace to file and close file" << std::endl
2861 << " 5) Reopen file in clean ROOT session to confirm that problems are fixed" << std::endl ;
2862 return false ;
2863 }
2864
2865 ++iter ;
2866 }
2867
2868 return true ;
2869}
2870
2871
2872
2873////////////////////////////////////////////////////////////////////////////////
2874/// Internal access to TDirectory append method
2875
2880
2881
2882////////////////////////////////////////////////////////////////////////////////
2883/// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
2884
2886{
2887 if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
2888 coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << std::endl ;
2889 } else {
2890 InternalAppend(obj) ;
2891 }
2892}
2893
2894
2895////////////////////////////////////////////////////////////////////////////////
2896/// Overload TDirectory interface method to prohibit insertion of objects in read-only directory workspace representation
2897
2899{
2900 if (dynamic_cast<RooAbsArg*>(obj) || dynamic_cast<RooAbsData*>(obj)) {
2901 coutE(ObjectHandling) << "RooWorkspace::WSDir::Add(" << GetName() << ") ERROR: Directory is read-only representation of a RooWorkspace, use RooWorkspace::import() to add objects" << std::endl ;
2902 } else {
2903 InternalAppend(obj) ;
2904 }
2905}
2906
2907
2908////////////////////////////////////////////////////////////////////////////////
2909/// If one of the TObject we have a referenced to is deleted, remove the
2910/// reference.
2911
2913{
2915 if (removedObj == _dir) _dir = nullptr;
2916
2918
2925
2926 std::vector<std::string> invalidSets;
2927
2928 for(auto &c : _namedSets) {
2929 auto const& setName = c.first;
2930 auto& set = c.second;
2931 std::size_t oldSize = set.size();
2932 set.RecursiveRemove(removedObj);
2933 // If the set is used internally by RooFit to cache parameters or
2934 // constraints, it is invalidated by object removal. We will keep track
2935 // of its name to remove the cache set later.
2936 if(set.size() < oldSize && isCacheSet(setName)) {
2937 invalidSets.emplace_back(setName);
2938 }
2939 }
2940
2941 // Remove the sets that got invalidated by the object removal
2942 for(std::string const& setName : invalidSets) {
2943 removeSet(setName.c_str());
2944 }
2945
2946 _eocache.RecursiveRemove(removedObj); // RooExpensiveObjectCache
2947}
#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:572
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:57
virtual const RooArgSet * get() const
Definition RooAbsData.h:101
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.
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:3764
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition TH1.cxx:1263
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:41
virtual void RecursiveRemove(TObject *obj)
Recursively remove this object from a list.
Definition TObject.cxx:678
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:981
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.