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