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