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