Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooJSONFactoryWSTool.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 * Authors:
4 * Carsten D. Burgard, DESY/ATLAS, Dec 2021
5 *
6 * Copyright (c) 2022, CERN
7 *
8 * Redistribution and use in source and binary forms,
9 * with or without modification, are permitted according to the terms
10 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
11 */
12
13#include <RooFitHS3/JSONIO.h>
15
16#include <RooConstVar.h>
17#include <RooRealVar.h>
18#include <RooBinning.h>
19#include <RooAbsCategory.h>
20#include <RooRealProxy.h>
21#include <RooListProxy.h>
22#include <RooAbsProxy.h>
23#include <RooCategory.h>
24#include <RooDataSet.h>
25#include <RooDataHist.h>
26#include <RooSimultaneous.h>
27#include <RooFormulaVar.h>
28#include <RooFit/ModelConfig.h>
29#include <RooFitImplHelpers.h>
30#include <RooAbsCollection.h>
31
32#include "JSONIOUtils.h"
33#include "Domains.h"
34
35#include "RooFitImplHelpers.h"
36
37#include <TROOT.h>
38
39#include <algorithm>
40#include <fstream>
41#include <iostream>
42#include <stack>
43#include <stdexcept>
44
45/** \class RooJSONFactoryWSTool
46\ingroup roofit_dev_docs_hs3
47
48When using \ref Roofitmain, statistical models can be conveniently handled and
49stored as a RooWorkspace. However, for the sake of interoperability
50with other statistical frameworks, and also ease of manipulation, it
51may be useful to store statistical models in text form.
52
53The RooJSONFactoryWSTool is a helper class to achieve exactly this,
54exporting to and importing from JSON.
55
56In order to import a workspace from a JSON file, you can do
57
58~~~ {.py}
59ws = ROOT.RooWorkspace("ws")
60tool = ROOT.RooJSONFactoryWSTool(ws)
61tool.importJSON("myjson.json")
62~~~
63
64Similarly, in order to export a workspace to a JSON file, you can do
65
66~~~ {.py}
67tool = ROOT.RooJSONFactoryWSTool(ws)
68tool.exportJSON("myjson.json")
69~~~
70
71Analogously, in C++, you can do
72
73~~~ {.cxx}
74#include "RooFitHS3/RooJSONFactoryWSTool.h"
75// ...
76RooWorkspace ws("ws");
77RooJSONFactoryWSTool tool(ws);
78tool.importJSON("myjson.json");
79~~~
80
81and
82
83~~~ {.cxx}
84#include "RooFitHS3/RooJSONFactoryWSTool.h"
85// ...
86RooJSONFactoryWSTool tool(ws);
87tool.exportJSON("myjson.json");
88~~~
89
90For more details, consult the tutorial <a href="rf515__hfJSON_8py.html">rf515_hfJSON</a>.
91
92The RooJSONFactoryWSTool only knows about a limited set of classes for
93import and export. If import or export of a class you're interested in
94fails, you might need to add your own importer or exporter. Please
95consult the relevant section in the \ref roofit_dev_docs to learn how to do that (\ref roofit_dev_docs_hs3).
96
97You can always get a list of all the available importers and exporters by calling the following functions:
98~~~ {.py}
99ROOT.RooFit.JSONIO.printImporters()
100ROOT.RooFit.JSONIO.printExporters()
101ROOT.RooFit.JSONIO.printFactoryExpressions()
102ROOT.RooFit.JSONIO.printExportKeys()
103~~~
104
105Alternatively, you can generate a LaTeX version of the available importers and exporters by calling
106~~~ {.py}
107tool = ROOT.RooJSONFactoryWSTool(ws)
108tool.writedoc("hs3.tex")
109~~~
110*/
111
112constexpr auto hs3VersionTag = "0.2";
113
116
117namespace {
118
119std::vector<std::string> valsToStringVec(JSONNode const &node)
120{
121 std::vector<std::string> out;
122 out.reserve(node.num_children());
123 for (JSONNode const &elem : node.children()) {
124 out.push_back(elem.val());
125 }
126 return out;
127}
128
129/**
130 * @brief Check if the number of components in CombinedData matches the number of categories in the RooSimultaneous PDF.
131 *
132 * This function checks whether the number of components in the provided CombinedData 'data' matches the number of
133 * categories in the provided RooSimultaneous PDF 'pdf'.
134 *
135 * @param data The reference to the CombinedData to be checked.
136 * @param pdf The pointer to the RooSimultaneous PDF for comparison.
137 * @return bool Returns true if the number of components in 'data' matches the number of categories in 'pdf'; otherwise,
138 * returns false.
139 */
140bool matches(const RooJSONFactoryWSTool::CombinedData &data, const RooSimultaneous *pdf)
141{
142 return data.components.size() == pdf->indexCat().size();
143}
144
145/**
146 * @struct Var
147 * @brief Structure to store variable information.
148 *
149 * This structure represents variable information such as the number of bins, minimum and maximum values,
150 * and a vector of binning edges for a variable.
151 */
152struct Var {
153 int nbins; // Number of bins
154 double min; // Minimum value
155 double max; // Maximum value
156 std::vector<double> edges; // Vector of edges
157
158 /**
159 * @brief Constructor for Var.
160 * @param n Number of bins.
161 */
162 Var(int n) : nbins(n), min(0), max(n) {}
163
164 /**
165 * @brief Constructor for Var from JSONNode.
166 * @param val JSONNode containing variable information.
167 */
168 Var(const JSONNode &val);
169};
170
171/**
172 * @brief Check if a string represents a valid number.
173 *
174 * This function checks whether the provided string 'str' represents a valid number.
175 * The function returns true if the entire string can be parsed as a number (integer or floating-point); otherwise, it
176 * returns false.
177 *
178 * @param str The string to be checked.
179 * @return bool Returns true if the string 'str' represents a valid number; otherwise, returns false.
180 */
181bool isNumber(const std::string &str)
182{
183 bool seen_digit = false;
184 bool seen_dot = false;
185 bool seen_e = false;
186 bool after_e = false;
187 bool sign_allowed = true;
188
189 for (size_t i = 0; i < str.size(); ++i) {
190 char c = str[i];
191
192 if (std::isdigit(c)) {
193 seen_digit = true;
194 sign_allowed = false;
195 } else if ((c == '+' || c == '-') && sign_allowed) {
196 // Sign allowed at the beginning or right after 'e'/'E'
197 sign_allowed = false;
198 } else if (c == '.' && !seen_dot && !after_e) {
199 seen_dot = true;
200 sign_allowed = false;
201 } else if ((c == 'e' || c == 'E') && seen_digit && !seen_e) {
202 seen_e = true;
203 after_e = true;
204 sign_allowed = true; // allow sign immediately after 'e'
205 seen_digit = false; // reset: we now expect digits after e
206 } else {
207 return false;
208 }
209 }
210
211 return seen_digit;
212}
213
214/**
215 * @brief Configure a RooRealVar based on information from a JSONNode.
216 *
217 * This function configures the provided RooRealVar 'v' based on the information provided in the JSONNode 'p'.
218 * The JSONNode 'p' contains information about various properties of the RooRealVar, such as its value, error, number of
219 * bins, etc. The function reads these properties from the JSONNode and sets the corresponding properties of the
220 * RooRealVar accordingly.
221 *
222 * @param domains The reference to the RooFit::JSONIO::Detail::Domains containing domain information for variables (not
223 * used in this function).
224 * @param p The JSONNode containing information about the properties of the RooRealVar 'v'.
225 * @param v The reference to the RooRealVar to be configured.
226 * @return void
227 */
229{
230 if (!p.has_child("name")) {
231 RooJSONFactoryWSTool::error("cannot instantiate variable without \"name\"!");
232 }
233 if (auto n = p.find("value"))
234 v.setVal(n->val_double());
235 domains.writeVariable(v);
236 if (auto n = p.find("nbins"))
237 v.setBins(n->val_int());
238 if (auto n = p.find("relErr"))
239 v.setError(v.getVal() * n->val_double());
240 if (auto n = p.find("err"))
241 v.setError(n->val_double());
242 if (auto n = p.find("const")) {
243 v.setConstant(n->val_bool());
244 } else {
245 v.setConstant(false);
246 }
247}
248
250{
251 auto paramPointsNode = rootNode.find("parameter_points");
252 if (!paramPointsNode)
253 return nullptr;
254 auto out = RooJSONFactoryWSTool::findNamedChild(*paramPointsNode, "default_values");
255 if (out == nullptr)
256 return nullptr;
257 return &((*out)["parameters"]);
258}
259
260Var::Var(const JSONNode &val)
261{
262 if (val.find("edges")) {
263 for (auto const &child : val.children()) {
264 this->edges.push_back(child.val_double());
265 }
266 this->nbins = this->edges.size();
267 this->min = this->edges[0];
268 this->max = this->edges[this->nbins - 1];
269 } else {
270 if (!val.find("nbins")) {
271 this->nbins = 1;
272 } else {
273 this->nbins = val["nbins"].val_int();
274 }
275 if (!val.find("min")) {
276 this->min = 0;
277 } else {
278 this->min = val["min"].val_double();
279 }
280 if (!val.find("max")) {
281 this->max = 1;
282 } else {
283 this->max = val["max"].val_double();
284 }
285 }
286}
287
288std::string genPrefix(const JSONNode &p, bool trailing_underscore)
289{
290 std::string prefix;
291 if (!p.is_map())
292 return prefix;
293 if (auto node = p.find("namespaces")) {
294 for (const auto &ns : node->children()) {
295 if (!prefix.empty())
296 prefix += "_";
297 prefix += ns.val();
298 }
299 }
300 if (trailing_underscore && !prefix.empty())
301 prefix += "_";
302 return prefix;
303}
304
305// helpers for serializing / deserializing binned datasets
306void genIndicesHelper(std::vector<std::vector<int>> &combinations, std::vector<int> &curr_comb,
307 const std::vector<int> &vars_numbins, size_t curridx)
308{
309 if (curridx == vars_numbins.size()) {
310 // we have filled a combination. Copy it.
311 combinations.emplace_back(curr_comb);
312 } else {
313 for (int i = 0; i < vars_numbins[curridx]; ++i) {
314 curr_comb[curridx] = i;
316 }
317 }
318}
319
320/**
321 * @brief Import attributes from a JSONNode into a RooAbsArg.
322 *
323 * This function imports attributes, represented by the provided JSONNode 'node', into the provided RooAbsArg 'arg'.
324 * The attributes are read from the JSONNode and applied to the RooAbsArg.
325 *
326 * @param arg The pointer to the RooAbsArg to which the attributes will be imported.
327 * @param node The JSONNode containing information about the attributes to be imported.
328 * @return void
329 */
330void importAttributes(RooAbsArg *arg, JSONNode const &node)
331{
332 if (auto seq = node.find("dict")) {
333 for (const auto &attr : seq->children()) {
334 arg->setStringAttribute(attr.key().c_str(), attr.val().c_str());
335 }
336 }
337 if (auto seq = node.find("tags")) {
338 for (const auto &attr : seq->children()) {
339 arg->setAttribute(attr.val().c_str());
340 }
341 }
342}
343
344// RooWSFactoryTool expression handling
345std::string generate(const RooFit::JSONIO::ImportExpression &ex, const JSONNode &p, RooJSONFactoryWSTool *tool)
346{
347 std::stringstream expression;
348 std::string classname(ex.tclass->GetName());
349 size_t colon = classname.find_last_of(':');
350 expression << (colon < classname.size() ? classname.substr(colon + 1) : classname);
351 bool first = true;
352 const auto &name = RooJSONFactoryWSTool::name(p);
353 for (auto k : ex.arguments) {
354 expression << (first ? "::" + name + "(" : ",");
355 first = false;
356 if (k == "true" || k == "false") {
357 expression << (k == "true" ? "1" : "0");
358 } else if (!p.has_child(k)) {
359 std::stringstream errMsg;
360 errMsg << "node '" << name << "' is missing key '" << k << "'";
362 } else if (p[k].is_seq()) {
363 bool firstInner = true;
364 expression << "{";
365 for (RooAbsArg *arg : tool->requestArgList<RooAbsReal>(p, k)) {
366 expression << (firstInner ? "" : ",") << arg->GetName();
367 firstInner = false;
368 }
369 expression << "}";
370 } else {
371 tool->requestArg<RooAbsReal>(p, p[k].key());
372 expression << p[k].val();
373 }
374 }
375 expression << ")";
376 return expression.str();
377}
378
379/**
380 * @brief Generate bin indices for a set of RooRealVars.
381 *
382 * This function generates all possible combinations of bin indices for the provided RooArgSet 'vars' containing
383 * RooRealVars. Each bin index represents a possible bin selection for the corresponding RooRealVar. The bin indices are
384 * stored in a vector of vectors, where each inner vector represents a combination of bin indices for all RooRealVars.
385 *
386 * @param vars The RooArgSet containing the RooRealVars for which bin indices will be generated.
387 * @return std::vector<std::vector<int>> A vector of vectors containing all possible combinations of bin indices.
388 */
389std::vector<std::vector<int>> generateBinIndices(const RooArgSet &vars)
390{
391 std::vector<std::vector<int>> combinations;
392 std::vector<int> vars_numbins;
393 vars_numbins.reserve(vars.size());
394 for (const auto *absv : static_range_cast<RooRealVar *>(vars)) {
395 vars_numbins.push_back(absv->getBins());
396 }
397 std::vector<int> curr_comb(vars.size());
399 return combinations;
400}
401
402template <typename... Keys_t>
403JSONNode const *findRooFitInternal(JSONNode const &node, Keys_t const &...keys)
404{
405 return node.find("misc", "ROOT_internal", keys...);
406}
407
408/**
409 * @brief Check if a RooAbsArg is a literal constant variable.
410 *
411 * This function checks whether the provided RooAbsArg 'arg' is a literal constant variable.
412 * A literal constant variable is a RooConstVar with a numeric value as a name.
413 *
414 * @param arg The reference to the RooAbsArg to be checked.
415 * @return bool Returns true if 'arg' is a literal constant variable; otherwise, returns false.
416 */
417bool isLiteralConstVar(RooAbsArg const &arg)
418{
419 bool isRooConstVar = dynamic_cast<RooConstVar const *>(&arg);
420 return isRooConstVar && isNumber(arg.GetName());
421}
422
423/**
424 * @brief Export attributes of a RooAbsArg to a JSONNode.
425 *
426 * This function exports the attributes of the provided RooAbsArg 'arg' to the JSONNode 'rootnode'.
427 *
428 * @param arg The pointer to the RooAbsArg from which attributes will be exported.
429 * @param rootnode The JSONNode to which the attributes will be exported.
430 * @return void
431 */
432void exportAttributes(const RooAbsArg *arg, JSONNode &rootnode)
433{
434 // If this RooConst is a literal number, we don't need to export the attributes.
435 if (isLiteralConstVar(*arg)) {
436 return;
437 }
438
439 JSONNode *node = nullptr;
440
441 auto initializeNode = [&]() {
442 if (node)
443 return;
444
445 node = &RooJSONFactoryWSTool::getRooFitInternal(rootnode, "attributes").set_map()[arg->GetName()].set_map();
446 };
447
448 // RooConstVars are not a thing in HS3, and also for RooFit they are not
449 // that important: they are just constants. So we don't need to remember
450 // any information about them.
451 if (dynamic_cast<RooConstVar const *>(arg)) {
452 return;
453 }
454
455 // export all string attributes of an object
456 if (!arg->stringAttributes().empty()) {
457 for (const auto &it : arg->stringAttributes()) {
458 // Skip some RooFit internals
459 if (it.first == "factory_tag" || it.first == "PROD_TERM_TYPE")
460 continue;
462 (*node)["dict"].set_map()[it.first] << it.second;
463 }
464 }
465 if (!arg->attributes().empty()) {
466 for (auto const &attr : arg->attributes()) {
467 // Skip some RooFit internals
468 if (attr == "SnapShot_ExtRefClone" || attr == "RooRealConstant_Factory_Object")
469 continue;
471 (*node)["tags"].set_seq().append_child() << attr;
472 }
473 }
474}
475
476/**
477 * @brief Create several observables in the workspace.
478 *
479 * This function obtains a list of observables from the provided
480 * RooWorkspace 'ws' based on their names given in the 'axes" field of
481 * the JSONNode 'node'. The observables are added to the RooArgSet
482 * 'out'.
483 *
484 * @param ws The RooWorkspace in which the observables will be created.
485 * @param node The JSONNode containing information about the observables to be created.
486 * @param out The RooArgSet to which the created observables will be added.
487 * @return void
488 */
489void getObservables(RooWorkspace const &ws, const JSONNode &node, RooArgSet &out)
490{
491 std::map<std::string, Var> vars;
492 for (const auto &p : node["axes"].children()) {
493 vars.emplace(RooJSONFactoryWSTool::name(p), Var(p));
494 }
495
496 for (auto v : vars) {
497 std::string name(v.first);
498 if (ws.var(name)) {
499 out.add(*ws.var(name));
500 } else {
501 std::stringstream errMsg;
502 errMsg << "The observable \"" << name << "\" could not be found in the workspace!";
504 }
505 }
506}
507
508/**
509 * @brief Import data from the JSONNode into the workspace.
510 *
511 * This function imports data, represented by the provided JSONNode 'p', into the workspace represented by the provided
512 * RooWorkspace. The data information is read from the JSONNode and added to the workspace.
513 *
514 * @param p The JSONNode representing the data to be imported.
515 * @param workspace The RooWorkspace to which the data will be imported.
516 * @return std::unique_ptr<RooAbsData> A unique pointer to the RooAbsData object representing the imported data.
517 * The caller is responsible for managing the memory of the returned object.
518 */
519std::unique_ptr<RooAbsData> loadData(const JSONNode &p, RooWorkspace &workspace)
520{
521 std::string name(RooJSONFactoryWSTool::name(p));
522
524
525 std::string const &type = p["type"].val();
526 if (type == "binned") {
527 // binned
529 } else if (type == "unbinned") {
530 // unbinned
531 RooArgSet vars;
532 getObservables(workspace, p, vars);
533 RooArgList varlist(vars);
534 auto data = std::make_unique<RooDataSet>(name, name, vars, RooFit::WeightVar());
535 auto &coords = p["entries"];
536 if (!coords.is_seq()) {
537 RooJSONFactoryWSTool::error("key 'entries' is not a list!");
538 }
539 std::vector<double> weightVals;
540 if (p.has_child("weights")) {
541 auto &weights = p["weights"];
542 if (coords.num_children() != weights.num_children()) {
543 RooJSONFactoryWSTool::error("inconsistent number of entries and weights!");
544 }
545 for (auto const &weight : weights.children()) {
546 weightVals.push_back(weight.val_double());
547 }
548 }
549 std::size_t i = 0;
550 for (auto const &point : coords.children()) {
551 if (!point.is_seq()) {
552 std::stringstream errMsg;
553 errMsg << "coordinate point '" << i << "' is not a list!";
555 }
556 if (point.num_children() != varlist.size()) {
557 RooJSONFactoryWSTool::error("inconsistent number of entries and observables!");
558 }
559 std::size_t j = 0;
560 for (auto const &pointj : point.children()) {
561 auto *v = static_cast<RooRealVar *>(varlist.at(j));
562 v->setVal(pointj.val_double());
563 ++j;
564 }
565 if (weightVals.size() > 0) {
566 data->add(vars, weightVals[i]);
567 } else {
568 data->add(vars, 1.);
569 }
570 ++i;
571 }
572 return data;
573 }
574
575 std::stringstream ss;
576 ss << "RooJSONFactoryWSTool() failed to create dataset " << name << std::endl;
578 return nullptr;
579}
580
581/**
582 * @brief Import an analysis from the JSONNode into the workspace.
583 *
584 * This function imports an analysis, represented by the provided JSONNodes 'analysisNode' and 'likelihoodsNode',
585 * into the workspace represented by the provided RooWorkspace. The analysis information is read from the JSONNodes
586 * and added to the workspace as one or more RooStats::ModelConfig objects.
587 *
588 * @param rootnode The root JSONNode representing the entire JSON file.
589 * @param analysisNode The JSONNode representing the analysis to be imported.
590 * @param likelihoodsNode The JSONNode containing information about likelihoods associated with the analysis.
591 * @param domainsNode The JSONNode containing information about domains associated with the analysis.
592 * @param workspace The RooWorkspace to which the analysis will be imported.
593 * @param datasets A vector of unique pointers to RooAbsData objects representing the data associated with the analysis.
594 * @return void
595 */
596void importAnalysis(const JSONNode &rootnode, const JSONNode &analysisNode, const JSONNode &likelihoodsNode,
597 const JSONNode &domainsNode, RooWorkspace &workspace,
598 const std::vector<std::unique_ptr<RooAbsData>> &datasets)
599{
600 // if this is a toplevel pdf, also create a modelConfig for it
602 JSONNode const *mcAuxNode = findRooFitInternal(rootnode, "ModelConfigs", analysisName);
603
604 JSONNode const *mcNameNode = mcAuxNode ? mcAuxNode->find("mcName") : nullptr;
605 std::string mcname = mcNameNode ? mcNameNode->val() : analysisName;
606 if (workspace.obj(mcname))
607 return;
608
609 workspace.import(RooStats::ModelConfig{mcname.c_str(), mcname.c_str()});
610 auto *mc = static_cast<RooStats::ModelConfig *>(workspace.obj(mcname));
611 mc->SetWS(workspace);
612
614 if (!nllNode) {
615 throw std::runtime_error("likelihood node not found!");
616 }
617 if (!nllNode->has_child("distributions")) {
618 throw std::runtime_error("likelihood node has no distributions attached!");
619 }
620 if (!nllNode->has_child("data")) {
621 throw std::runtime_error("likelihood node has no data attached!");
622 }
623 std::vector<std::string> nllDistNames = valsToStringVec((*nllNode)["distributions"]);
625 for (auto &nameNode : (*nllNode)["aux_distributions"].children()) {
626 if (RooAbsArg *extConstraint = workspace.arg(nameNode.val())) {
628 }
629 }
630 RooArgSet observables;
631 for (auto &nameNode : (*nllNode)["data"].children()) {
632 bool found = false;
633 for (const auto &d : datasets) {
634 if (d->GetName() == nameNode.val()) {
635 found = true;
636 observables.add(*d->get());
637 }
638 }
639 if (nameNode.val() != "0" && !found)
640 throw std::runtime_error("dataset '" + nameNode.val() + "' cannot be found!");
641 }
642
643 JSONNode const *pdfNameNode = mcAuxNode ? mcAuxNode->find("pdfName") : nullptr;
644 std::string const pdfName = pdfNameNode ? pdfNameNode->val() : "simPdf";
645
646 RooAbsPdf *pdf = static_cast<RooSimultaneous *>(workspace.pdf(pdfName));
647
648 if (!pdf) {
649 // if there is no simultaneous pdf, we can check whether there is only one pdf in the list
650 if (nllDistNames.size() == 1) {
651 // if so, we can use that one to populate the ModelConfig
652 pdf = workspace.pdf(nllDistNames[0]);
653 } else {
654 // otherwise, we have no choice but to build a simPdf by hand
655 std::string simPdfName = analysisName + "_simPdf";
656 std::string indexCatName = analysisName + "_categoryIndex";
657 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
658 std::map<std::string, RooAbsPdf *> pdfMap;
659 for (std::size_t i = 0; i < nllDistNames.size(); ++i) {
660 indexCat.defineType(nllDistNames[i], i);
661 pdfMap[nllDistNames[i]] = workspace.pdf(nllDistNames[i]);
662 }
663 RooSimultaneous simPdf{simPdfName.c_str(), simPdfName.c_str(), pdfMap, indexCat};
665 pdf = static_cast<RooSimultaneous *>(workspace.pdf(simPdfName));
666 }
667 }
668
669 mc->SetPdf(*pdf);
670
671 if (!extConstraints.empty())
672 mc->SetExternalConstraints(extConstraints);
673
674 auto readArgSet = [&](std::string const &name) {
675 RooArgSet out;
676 for (auto const &child : analysisNode[name].children()) {
677 out.add(*workspace.arg(child.val()));
678 }
679 return out;
680 };
681
682 mc->SetParametersOfInterest(readArgSet("parameters_of_interest"));
683 mc->SetObservables(observables);
684 RooArgSet pars;
685 pdf->getParameters(&observables, pars);
686
687 // Figure out the set parameters that appear in the main measurement:
688 // getAllConstraints() has the side effect to remove all parameters from
689 // "mainPars" that are not part of any pdf over observables.
690 RooArgSet mainPars{pars};
691 pdf->getAllConstraints(observables, mainPars, /*stripDisconnected*/ true);
692
694 for (auto &domain : analysisNode["domains"].children()) {
696 if (!thisDomain || !thisDomain->has_child("axes"))
697 continue;
698 for (auto &var : (*thisDomain)["axes"].children()) {
699 auto *wsvar = workspace.var(RooJSONFactoryWSTool::name(var));
700 if (wsvar)
701 domainPars.add(*wsvar);
702 }
703 }
704
706 RooArgSet globs;
707 for (const auto &p : pars) {
708 if (mc->GetParametersOfInterest()->find(*p))
709 continue;
710 if (p->isConstant() && !mainPars.find(*p) && domainPars.find(*p)) {
711 globs.add(*p);
712 } else if (domainPars.find(*p)) {
713 nps.add(*p);
714 }
715 }
716
717 mc->SetGlobalObservables(globs);
718 mc->SetNuisanceParameters(nps);
719
720 if (mcAuxNode) {
721 if (auto found = mcAuxNode->find("combined_data_name")) {
722 pdf->setStringAttribute("combined_data_name", found->val().c_str());
723 }
724 }
725
726 if (analysisNode.has_child("init") && workspace.getSnapshot(analysisNode["init"].val().c_str())) {
727 mc->SetSnapshot(*workspace.getSnapshot(analysisNode["init"].val().c_str()));
728 }
729}
730
731void combinePdfs(const JSONNode &rootnode, RooWorkspace &ws)
732{
733 auto *combinedPdfInfoNode = findRooFitInternal(rootnode, "combined_distributions");
734
735 // If there is no info on combining pdfs
736 if (combinedPdfInfoNode == nullptr) {
737 return;
738 }
739
740 for (auto &info : combinedPdfInfoNode->children()) {
741
742 // parse the information
743 std::string combinedName = info.key();
744 std::string indexCatName = info["index_cat"].val();
745 std::vector<std::string> labels = valsToStringVec(info["labels"]);
746 std::vector<int> indices;
747 std::vector<std::string> pdfNames = valsToStringVec(info["distributions"]);
748 for (auto &n : info["indices"].children()) {
749 indices.push_back(n.val_int());
750 }
751
752 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
753 std::map<std::string, RooAbsPdf *> pdfMap;
754
755 for (std::size_t iChannel = 0; iChannel < labels.size(); ++iChannel) {
756 indexCat.defineType(labels[iChannel], indices[iChannel]);
757 pdfMap[labels[iChannel]] = ws.pdf(pdfNames[iChannel]);
758 }
759
760 RooSimultaneous simPdf{combinedName.c_str(), combinedName.c_str(), pdfMap, indexCat};
762 }
763}
764
765void combineDatasets(const JSONNode &rootnode, std::vector<std::unique_ptr<RooAbsData>> &datasets)
766{
767 auto *combinedDataInfoNode = findRooFitInternal(rootnode, "combined_datasets");
768
769 // If there is no info on combining datasets
770 if (combinedDataInfoNode == nullptr) {
771 return;
772 }
773
774 for (auto &info : combinedDataInfoNode->children()) {
775
776 // parse the information
777 std::string combinedName = info.key();
778 std::string indexCatName = info["index_cat"].val();
779 std::vector<std::string> labels = valsToStringVec(info["labels"]);
780 std::vector<int> indices;
781 for (auto &n : info["indices"].children()) {
782 indices.push_back(n.val_int());
783 }
784 if (indices.size() != labels.size()) {
785 RooJSONFactoryWSTool::error("mismatch in number of indices and labels!");
786 }
787
788 // Create the combined dataset for RooFit
789 std::map<std::string, std::unique_ptr<RooAbsData>> dsMap;
790 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
791 RooArgSet allVars{indexCat};
792 for (std::size_t iChannel = 0; iChannel < labels.size(); ++iChannel) {
793 auto componentName = combinedName + "_" + labels[iChannel];
794 // We move the found channel data out of the "datasets" vector, such that
795 // the data components don't get imported anymore.
796 std::unique_ptr<RooAbsData> &component = *std::find_if(
797 datasets.begin(), datasets.end(), [&](auto &d) { return d && d->GetName() == componentName; });
798 if (!component)
799 RooJSONFactoryWSTool::error("unable to obtain component matching component name '" + componentName + "'");
800 allVars.add(*component->get());
801 dsMap.insert({labels[iChannel], std::move(component)});
802 indexCat.defineType(labels[iChannel], indices[iChannel]);
803 }
804
805 auto combined = std::make_unique<RooDataSet>(combinedName, combinedName, allVars, RooFit::Import(dsMap),
806 RooFit::Index(indexCat));
807 datasets.emplace_back(std::move(combined));
808 }
809}
810
811template <class T>
812void sortByName(T &coll)
813{
814 std::sort(coll.begin(), coll.end(), [](auto &l, auto &r) { return strcmp(l->GetName(), r->GetName()) < 0; });
815}
816
817} // namespace
818
820
822
824{
825 const size_t old_children = node.num_children();
826 node.set_seq();
827 size_t n = 0;
828 for (RooAbsArg const *arg : coll) {
829 if (n >= nMax)
830 break;
831 if (isLiteralConstVar(*arg)) {
832 node.append_child() << static_cast<RooConstVar const *>(arg)->getVal();
833 } else {
834 node.append_child() << arg->GetName();
835 }
836 ++n;
837 }
838 if (node.num_children() != old_children + coll.size()) {
839 error("unable to stream collection " + std::string(coll.GetName()) + " to " + node.key());
840 }
841}
842
844{
845 const size_t old_children = node.num_children();
846 node.set_seq();
847 size_t n = 0;
848 for (RooAbsArg const *arg : coll) {
849 if (n >= nMax)
850 break;
851 if (isLiteralConstVar(*arg)) {
852 node.append_child() << static_cast<RooConstVar const *>(arg)->getVal();
853 } else {
854 node.append_child() << sanitizeName(arg->GetName());
855 }
856 ++n;
857 }
858 if (node.num_children() != old_children + coll.size()) {
859 error("unable to stream collection " + std::string(coll.GetName()) + " to " + node.key());
860 }
861}
862
864{
866 return node.set_map()[name].set_map();
867 }
868 JSONNode &child = node.set_seq().append_child().set_map();
869 child["name"] << name;
870 return child;
871}
872
873JSONNode const *RooJSONFactoryWSTool::findNamedChild(JSONNode const &node, std::string const &name)
874{
876 if (!node.is_map())
877 return nullptr;
878 return node.find(name);
879 }
880 if (!node.is_seq())
881 return nullptr;
882 for (JSONNode const &child : node.children()) {
883 if (child["name"].val() == name)
884 return &child;
885 }
886
887 return nullptr;
888}
889
890/**
891 * @brief Check if a string is a valid name.
892 *
893 * A valid name should start with a letter or an underscore, followed by letters, digits, or underscores.
894 * Only characters from the ASCII character set are allowed.
895 *
896 * @param str The string to be checked.
897 * @return bool Returns true if the string is a valid name; otherwise, returns false.
898 */
899bool RooJSONFactoryWSTool::isValidName(const std::string &str)
900{
901 // Check if the string is empty or starts with a non-letter/non-underscore character
902 if (str.empty() || !(std::isalpha(str[0]) || str[0] == '_')) {
903 return false;
904 }
905
906 // Check the remaining characters in the string
907 for (char c : str) {
908 // Allow letters, digits, and underscore
909 if (!(std::isalnum(c) || c == '_')) {
910 return false;
911 }
912 }
913
914 // If all characters are valid, the string is a valid name
915 return true;
916}
917
921{
923 std::stringstream ss;
924 ss << "RooJSONFactoryWSTool() name '" << name << "' is not valid!" << std::endl
925 << "Sanitize names by setting RooJSONFactoryWSTool::allowSanitizeNames = True." << std::endl;
928 return false;
929 } else {
931 }
932 }
933 return true;
934}
935
937{
938 return useListsInsteadOfDicts ? n["name"].val() : n.key();
939}
940
942{
943 return appendNamedChild(rootNode["parameter_points"], "default_values")["parameters"];
944}
945
946template <>
947RooRealVar *RooJSONFactoryWSTool::requestImpl<RooRealVar>(const std::string &objname)
948{
950 return retval;
951 if (const auto *vars = getVariablesNode(*_rootnodeInput)) {
952 if (const auto &node = vars->find(objname)) {
953 this->importVariable(*node);
955 return retval;
956 }
957 }
958 return nullptr;
959}
960
961template <>
962RooAbsPdf *RooJSONFactoryWSTool::requestImpl<RooAbsPdf>(const std::string &objname)
963{
965 return retval;
966 if (const auto &distributionsNode = _rootnodeInput->find("distributions")) {
967 if (const auto &child = findNamedChild(*distributionsNode, objname)) {
968 this->importFunction(*child, true);
970 return retval;
971 }
972 }
973 return nullptr;
974}
975
976template <>
977RooAbsReal *RooJSONFactoryWSTool::requestImpl<RooAbsReal>(const std::string &objname)
978{
980 return retval;
981 if (isNumber(objname))
984 return pdf;
986 return var;
987 if (const auto &functionNode = _rootnodeInput->find("functions")) {
988 if (const auto &child = findNamedChild(*functionNode, objname)) {
989 this->importFunction(*child, true);
991 return retval;
992 }
993 }
994 return nullptr;
995}
996
997/**
998 * @brief Export a variable from the workspace to a JSONNode.
999 *
1000 * This function exports a variable, represented by the provided RooAbsArg pointer 'v', from the workspace to a
1001 * JSONNode. The variable's information is added to the JSONNode as key-value pairs.
1002 *
1003 * @param v The pointer to the RooAbsArg representing the variable to be exported.
1004 * @param node The JSONNode to which the variable will be exported.
1005 * @return void
1006 */
1008{
1009 auto *cv = dynamic_cast<const RooConstVar *>(v);
1010 auto *rrv = dynamic_cast<const RooRealVar *>(v);
1011 if (!cv && !rrv)
1012 return;
1013
1014 // for RooConstVar, if name and value are the same, we don't need to do anything
1015 if (cv && strcmp(cv->GetName(), TString::Format("%g", cv->getVal()).Data()) == 0) {
1016 return;
1017 }
1018
1019 // this variable was already exported
1020 if (findNamedChild(node, v->GetName())) {
1021 return;
1022 }
1023
1024 JSONNode &var = appendNamedChild(node, v->GetName());
1025
1026 if (cv) {
1027 var["value"] << cv->getVal();
1028 var["const"] << true;
1029 } else if (rrv) {
1030 var["value"] << rrv->getVal();
1031 if (rrv->isConstant()) {
1032 var["const"] << rrv->isConstant();
1033 }
1034 if (rrv->getBins() != 100) {
1035 var["nbins"] << rrv->getBins();
1036 }
1037 _domains->readVariable(*rrv);
1038 }
1039}
1040
1041/**
1042 * @brief Export variables from the workspace to a JSONNode.
1043 *
1044 * This function exports variables, represented by the provided RooArgSet, from the workspace to a JSONNode.
1045 * The variables' information is added to the JSONNode as key-value pairs.
1046 *
1047 * @param allElems The RooArgSet representing the variables to be exported.
1048 * @param n The JSONNode to which the variables will be exported.
1049 * @return void
1050 */
1052{
1053 // export a list of RooRealVar objects
1054 n.set_seq();
1055 for (RooAbsArg *arg : allElems) {
1056 exportVariable(arg, n);
1057 }
1058}
1059
1061 const std::string &formula)
1062{
1063 std::string newname = std::string(original->GetName()) + suffix;
1065 trafo_node["type"] << "generic_function";
1066 trafo_node["expression"] << TString::Format(formula.c_str(), original->GetName()).Data();
1067 this->setAttribute(newname, "roofit_skip"); // this function should not be imported back in
1068 return newname;
1069}
1070
1071/**
1072 * @brief Export an object from the workspace to a JSONNode.
1073 *
1074 * This function exports an object, represented by the provided RooAbsArg, from the workspace to a JSONNode.
1075 * The object's information is added to the JSONNode as key-value pairs.
1076 *
1077 * @param func The RooAbsArg representing the object to be exported.
1078 * @param exportedObjectNames A set of strings containing names of previously exported objects to avoid duplicates.
1079 * This set is updated with the name of the newly exported object.
1080 * @return void
1081 */
1083{
1084 // const std::string name = sanitizeName(func.GetName());
1085 std::string name = func.GetName();
1086
1087 // if this element was already exported, skip
1089 return;
1090
1091 exportedObjectNames.insert(name);
1092
1093 if (auto simPdf = dynamic_cast<RooSimultaneous const *>(&func)) {
1094 // RooSimultaneous is not used in the HS3 standard, we only export the
1095 // dependents and some ROOT internal information.
1097
1098 std::vector<std::string> channelNames;
1099 for (auto const &item : simPdf->indexCat()) {
1100 channelNames.push_back(item.first);
1101 }
1102
1103 auto &infoNode = getRooFitInternal(*_rootnodeOutput, "combined_distributions").set_map();
1104 auto &child = infoNode[simPdf->GetName()].set_map();
1105 child["index_cat"] << simPdf->indexCat().GetName();
1106 exportCategory(simPdf->indexCat(), child);
1107 child["distributions"].set_seq();
1108 for (auto const &item : simPdf->indexCat()) {
1109 child["distributions"].append_child() << simPdf->getPdf(item.first.c_str())->GetName();
1110 }
1111
1112 return;
1113 } else if (dynamic_cast<RooAbsCategory const *>(&func)) {
1114 // categories are created by the respective RooSimultaneous, so we're skipping the export here
1115 return;
1116 } else if (dynamic_cast<RooRealVar const *>(&func) || dynamic_cast<RooConstVar const *>(&func)) {
1117 exportVariable(&func, *_varsNode);
1118 return;
1119 }
1120
1121 auto &collectionNode = (*_rootnodeOutput)[dynamic_cast<RooAbsPdf const *>(&func) ? "distributions" : "functions"];
1122
1123 auto const &exporters = RooFit::JSONIO::exporters();
1124 auto const &exportKeys = RooFit::JSONIO::exportKeys();
1125
1126 TClass *cl = func.IsA();
1127
1129
1130 auto it = exporters.find(cl);
1131 if (it != exporters.end()) { // check if we have a specific exporter available
1132 for (auto &exp : it->second) {
1133 _serversToExport.clear();
1134 _serversToDelete.clear();
1135 if (!exp->exportObject(this, &func, elem)) {
1136 // The exporter might have messed with the content of the node
1137 // before failing. That's why we clear it and only reset the name.
1138 elem.clear();
1139 elem.set_map();
1141 elem["name"] << name;
1142 }
1143 continue;
1144 }
1145 if (exp->autoExportDependants()) {
1147 } else {
1149 }
1150 for (auto &s : _serversToDelete) {
1151 delete s;
1152 }
1153 return;
1154 }
1155 }
1156
1157 // generic export using the factory expressions
1158 const auto &dict = exportKeys.find(cl);
1159 if (dict == exportKeys.end()) {
1160 std::cerr << "unable to export class '" << cl->GetName() << "' - no export keys available!\n"
1161 << "there are several possible reasons for this:\n"
1162 << " 1. " << cl->GetName() << " is a custom class that you or some package you are using added.\n"
1163 << " 2. " << cl->GetName()
1164 << " is a ROOT class that nobody ever bothered to write a serialization definition for.\n"
1165 << " 3. something is wrong with your setup, e.g. you might have called "
1166 "RooFit::JSONIO::clearExportKeys() and/or never successfully read a file defining these "
1167 "keys with RooFit::JSONIO::loadExportKeys(filename)\n"
1168 << "either way, please make sure that:\n"
1169 << " 3: you are reading a file with export keys - call RooFit::JSONIO::printExportKeys() to "
1170 "see what is available\n"
1171 << " 2 & 1: you might need to write a serialization definition yourself. check "
1172 "https://root.cern/doc/master/group__roofit__dev__docs__hs3.html to "
1173 "see how to do this!\n";
1174 return;
1175 }
1176
1177 elem["type"] << dict->second.type;
1178
1179 size_t nprox = func.numProxies();
1180
1181 for (size_t i = 0; i < nprox; ++i) {
1182 RooAbsProxy *p = func.getProxy(i);
1183 if (!p)
1184 continue;
1185
1186 // some proxies start with a "!". This is a magic symbol that we don't want to stream
1187 std::string pname(p->name());
1188 if (pname[0] == '!')
1189 pname.erase(0, 1);
1190
1191 auto k = dict->second.proxies.find(pname);
1192 if (k == dict->second.proxies.end()) {
1193 std::cerr << "failed to find key matching proxy '" << pname << "' for type '" << dict->second.type
1194 << "', encountered in '" << func.GetName() << "', skipping" << std::endl;
1195 return;
1196 }
1197
1198 // empty string is interpreted as an instruction to ignore this value
1199 if (k->second.empty())
1200 continue;
1201
1202 if (auto l = dynamic_cast<RooAbsCollection *>(p)) {
1203 fillSeq(elem[k->second], *l);
1204 }
1205 if (auto r = dynamic_cast<RooArgProxy *>(p)) {
1206 if (isLiteralConstVar(*r->absArg())) {
1207 elem[k->second] << static_cast<RooConstVar *>(r->absArg())->getVal();
1208 } else {
1209 elem[k->second] << r->absArg()->GetName();
1210 }
1211 }
1212 }
1213
1214 // export all the servers of a given RooAbsArg
1215 for (RooAbsArg *s : func.servers()) {
1216 if (!s) {
1217 std::cerr << "unable to locate server of " << func.GetName() << std::endl;
1218 continue;
1219 }
1221 }
1222}
1223
1224/**
1225 * @brief Import a function from the JSONNode into the workspace.
1226 *
1227 * This function imports a function from the given JSONNode into the workspace.
1228 * The function's information is read from the JSONNode and added to the workspace.
1229 *
1230 * @param p The JSONNode representing the function to be imported.
1231 * @param importAllDependants A boolean flag indicating whether to import all dependants (servers) of the function.
1232 * @return void
1233 */
1235{
1236 std::string name(RooJSONFactoryWSTool::name(p));
1237
1238 // If this node if marked to be skipped by RooFit, exit
1239 if (hasAttribute(name, "roofit_skip")) {
1240 return;
1241 }
1242
1243 auto const &importers = RooFit::JSONIO::importers();
1245
1246 // some preparations: what type of function are we dealing with here?
1248
1249 // if the RooAbsArg already exists, we don't need to do anything
1250 if (_workspace.arg(name)) {
1251 return;
1252 }
1253 // if the key we found is not a map, it's an error
1254 if (!p.is_map()) {
1255 std::stringstream ss;
1256 ss << "RooJSONFactoryWSTool() function node " + name + " is not a map!";
1258 return;
1259 }
1260 std::string prefix = genPrefix(p, true);
1261 if (!prefix.empty())
1262 name = prefix + name;
1263 if (!p.has_child("type")) {
1264 std::stringstream ss;
1265 ss << "RooJSONFactoryWSTool() no type given for function '" << name << "', skipping." << std::endl;
1267 return;
1268 }
1269
1270 std::string functype(p["type"].val());
1271
1272 // import all dependents if importing a workspace, not for creating new objects
1273 if (!importAllDependants) {
1274 this->importDependants(p);
1275 }
1276
1277 // check for specific implementations
1278 auto it = importers.find(functype);
1279 bool ok = false;
1280 if (it != importers.end()) {
1281 for (auto &imp : it->second) {
1282 ok = imp->importArg(this, p);
1283 if (ok)
1284 break;
1285 }
1286 }
1287 if (!ok) { // generic import using the factory expressions
1288 auto expr = factoryExpressions.find(functype);
1289 if (expr != factoryExpressions.end()) {
1290 std::string expression = ::generate(expr->second, p, this);
1291 if (!_workspace.factory(expression)) {
1292 std::stringstream ss;
1293 ss << "RooJSONFactoryWSTool() failed to create " << expr->second.tclass->GetName() << " '" << name
1294 << "', skipping. expression was\n"
1295 << expression << std::endl;
1297 }
1298 } else {
1299 std::stringstream ss;
1300 ss << "RooJSONFactoryWSTool() no handling for type '" << functype << "' implemented, skipping."
1301 << "\n"
1302 << "there are several possible reasons for this:\n"
1303 << " 1. " << functype << " is a custom type that is not available in RooFit.\n"
1304 << " 2. " << functype
1305 << " is a ROOT class that nobody ever bothered to write a deserialization definition for.\n"
1306 << " 3. something is wrong with your setup, e.g. you might have called "
1307 "RooFit::JSONIO::clearFactoryExpressions() and/or never successfully read a file defining "
1308 "these expressions with RooFit::JSONIO::loadFactoryExpressions(filename)\n"
1309 << "either way, please make sure that:\n"
1310 << " 3: you are reading a file with factory expressions - call "
1311 "RooFit::JSONIO::printFactoryExpressions() "
1312 "to see what is available\n"
1313 << " 2 & 1: you might need to write a deserialization definition yourself. check "
1314 "https://root.cern/doc/master/group__roofit__dev__docs__hs3.html to see "
1315 "how to do this!"
1316 << std::endl;
1318 return;
1319 }
1320 }
1322 if (!func) {
1323 std::stringstream err;
1324 err << "something went wrong importing function '" << name << "'.";
1325 RooJSONFactoryWSTool::error(err.str());
1326 }
1327}
1328
1329/**
1330 * @brief Import a function from a JSON string into the workspace.
1331 *
1332 * This function imports a function from the provided JSON string into the workspace.
1333 * The function's information is read from the JSON string and added to the workspace.
1334 *
1335 * @param jsonString The JSON string containing the function information.
1336 * @param importAllDependants A boolean flag indicating whether to import all dependants (servers) of the function.
1337 * @return void
1338 */
1340{
1341 this->importFunction((JSONTree::create(jsonString))->rootnode(), importAllDependants);
1342}
1343
1344/**
1345 * @brief Export histogram data to a JSONNode.
1346 *
1347 * This function exports histogram data, represented by the provided variables and contents, to a JSONNode.
1348 * The histogram's axes information and bin contents are added as key-value pairs to the JSONNode.
1349 *
1350 * @param vars The RooArgSet representing the variables associated with the histogram.
1351 * @param n The number of bins in the histogram.
1352 * @param contents A pointer to the array containing the bin contents of the histogram.
1353 * @param output The JSONNode to which the histogram data will be exported.
1354 * @return void
1355 */
1356void RooJSONFactoryWSTool::exportHisto(RooArgSet const &vars, std::size_t n, double const *contents, JSONNode &output)
1357{
1358 auto &observablesNode = output["axes"].set_seq();
1359 // axes have to be ordered to get consistent bin indices
1360 for (auto *var : static_range_cast<RooRealVar *>(vars)) {
1361 std::string name = var->GetName();
1363 JSONNode &obsNode = observablesNode.append_child().set_map();
1364 obsNode["name"] << name;
1365 if (var->getBinning().isUniform()) {
1366 obsNode["min"] << var->getMin();
1367 obsNode["max"] << var->getMax();
1368 obsNode["nbins"] << var->getBins();
1369 } else {
1370 auto &edges = obsNode["edges"];
1371 edges.set_seq();
1372 double val = var->getBinning().binLow(0);
1373 edges.append_child() << val;
1374 for (int i = 0; i < var->getBinning().numBins(); ++i) {
1375 val = var->getBinning().binHigh(i);
1376 edges.append_child() << val;
1377 }
1378 }
1379 }
1380
1381 return exportArray(n, contents, output["contents"]);
1382}
1383
1384/**
1385 * @brief Export an array of doubles to a JSONNode.
1386 *
1387 * This function exports an array of doubles, represented by the provided size and contents,
1388 * to a JSONNode. The array elements are added to the JSONNode as a sequence of values.
1389 *
1390 * @param n The size of the array.
1391 * @param contents A pointer to the array containing the double values.
1392 * @param output The JSONNode to which the array will be exported.
1393 * @return void
1394 */
1395void RooJSONFactoryWSTool::exportArray(std::size_t n, double const *contents, JSONNode &output)
1396{
1397 output.set_seq();
1398 for (std::size_t i = 0; i < n; ++i) {
1399 double w = contents[i];
1400 // To make sure there are no unnecessary floating points in the JSON
1401 if (int(w) == w) {
1402 output.append_child() << int(w);
1403 } else {
1404 output.append_child() << w;
1405 }
1406 }
1407}
1408
1409/**
1410 * @brief Export a RooAbsCategory object to a JSONNode.
1411 *
1412 * This function exports a RooAbsCategory object, represented by the provided categories and indices,
1413 * to a JSONNode. The category labels and corresponding indices are added to the JSONNode as key-value pairs.
1414 *
1415 * @param cat The RooAbsCategory object to be exported.
1416 * @param node The JSONNode to which the category data will be exported.
1417 * @return void
1418 */
1420{
1421 auto &labels = node["labels"].set_seq();
1422 auto &indices = node["indices"].set_seq();
1423
1424 for (auto const &item : cat) {
1425 std::string label;
1426 if (std::isalpha(item.first[0])) {
1428 if (label != item.first) {
1429 oocoutW(nullptr, IO) << "RooFitHS3: changed '" << item.first << "' to '" << label
1430 << "' to become a valid name";
1431 }
1432 } else {
1433 RooJSONFactoryWSTool::error("refusing to change first character of string '" + item.first +
1434 "' to make a valid name!");
1435 label = item.first;
1436 }
1437 labels.append_child() << label;
1438 indices.append_child() << item.second;
1439 }
1440}
1441
1442/**
1443 * @brief Export combined data from the workspace to a custom struct.
1444 *
1445 * This function exports combined data from the workspace, represented by the provided RooAbsData object,
1446 * to a CombinedData struct. The struct contains information such as variables, categories,
1447 * and bin contents of the combined data.
1448 *
1449 * @param data The RooAbsData object representing the combined data to be exported.
1450 * @return CombinedData A custom struct containing the exported combined data.
1451 */
1453{
1454 // find category observables
1455 RooAbsCategory *cat = nullptr;
1456 for (RooAbsArg *obs : *data.get()) {
1457 if (dynamic_cast<RooAbsCategory *>(obs)) {
1458 if (cat) {
1459 RooJSONFactoryWSTool::error("dataset '" + std::string(data.GetName()) +
1460 " has several category observables!");
1461 }
1462 cat = static_cast<RooAbsCategory *>(obs);
1463 }
1464 }
1465
1466 // prepare return value
1468
1469 if (!cat)
1470 return datamap;
1471 // this is a combined dataset
1472
1473 datamap.name = data.GetName();
1474
1475 // Write information necessary to reconstruct the combined dataset upon import
1476 auto &child = getRooFitInternal(*_rootnodeOutput, "combined_datasets").set_map()[data.GetName()].set_map();
1477 child["index_cat"] << cat->GetName();
1478 exportCategory(*cat, child);
1479
1480 // Find a RooSimultaneous model that would fit to this dataset
1481 RooSimultaneous const *simPdf = nullptr;
1482 auto *combinedPdfInfoNode = findRooFitInternal(*_rootnodeOutput, "combined_distributions");
1483 if (combinedPdfInfoNode) {
1484 for (auto &info : combinedPdfInfoNode->children()) {
1485 if (info["index_cat"].val() == cat->GetName()) {
1486 simPdf = static_cast<RooSimultaneous const *>(_workspace.pdf(info.key()));
1487 }
1488 }
1489 }
1490
1491 // If there is an associated simultaneous pdf for the index category, we
1492 // use the RooAbsData::split() overload that takes the RooSimultaneous.
1493 // Like this, the observables that are not relevant for a given channel
1494 // are automatically split from the component datasets.
1495 std::vector<std::unique_ptr<RooAbsData>> dataList{simPdf ? data.split(*simPdf, true) : data.split(*cat, true)};
1496
1497 for (std::unique_ptr<RooAbsData> const &absData : dataList) {
1498 std::string catName(absData->GetName());
1499 std::string dataName;
1500 if (std::isalpha(catName[0])) {
1502 if (dataName != catName) {
1503 oocoutW(nullptr, IO) << "RooFitHS3: changed '" << catName << "' to '" << dataName
1504 << "' to become a valid name";
1505 }
1506 } else {
1507 RooJSONFactoryWSTool::error("refusing to change first character of string '" + catName +
1508 "' to make a valid name!");
1509 dataName = catName;
1510 }
1511 absData->SetName((std::string(data.GetName()) + "_" + dataName).c_str());
1512 datamap.components[catName] = absData->GetName();
1513 this->exportData(*absData);
1514 }
1515 return datamap;
1516}
1517
1518/**
1519 * @brief Export data from the workspace to a JSONNode.
1520 *
1521 * This function exports data represented by the provided RooAbsData object,
1522 * to a JSONNode. The data's information is added as key-value pairs to the JSONNode.
1523 *
1524 * @param data The RooAbsData object representing the data to be exported.
1525 * @return void
1526 */
1528{
1529 // find category observables
1530
1531 RooAbsCategory *cat = nullptr;
1532 for (RooAbsArg *obs : *data.get()) {
1533 if (dynamic_cast<RooAbsCategory *>(obs)) {
1534 if (cat) {
1535 RooJSONFactoryWSTool::error("dataset '" + std::string(data.GetName()) +
1536 " has several category observables!");
1537 }
1538 cat = static_cast<RooAbsCategory *>(obs);
1539 }
1540 }
1541
1542 if (cat)
1543 return;
1544
1545 JSONNode &output = appendNamedChild((*_rootnodeOutput)["data"], data.GetName());
1546 /*std::ofstream file("/home/scello/Data/ZvvH126_5.txt", std::ios::app);
1547 if (!file.is_open()) {
1548 std::cerr << "Error: Could not open file for writing.\n";
1549 }*/
1550
1551 // This works around a problem in RooStats/HistFactory that was only fixed
1552 // in ROOT 6.30: until then, the weight variable of the observed dataset,
1553 // called "weightVar", was added to the observables. Therefore, it also got
1554 // added to the Asimov dataset. But the Asimov has its own weight variable,
1555 // called "binWeightAsimov", making "weightVar" an actual observable in the
1556 // Asimov data. But this is only by accident and should be removed.
1557 RooArgSet variables = *data.get();
1558 if (auto weightVar = variables.find("weightVar")) {
1559 variables.remove(*weightVar);
1560 }
1561
1562 // this is a regular binned dataset
1563 if (auto dh = dynamic_cast<RooDataHist const *>(&data)) {
1564 output["type"] << "binned";
1565 for (auto *var : static_range_cast<RooRealVar *>(variables)) {
1566 _domains->readVariable(*var);
1567 }
1568 return exportHisto(variables, dh->numEntries(), dh->weightArray(), output);
1569 }
1570
1571 // Check if this actually represents a binned dataset, and then import it
1572 // like a RooDataHist. This happens frequently when people create combined
1573 // RooDataSets from binned data to fit HistFactory models. In this case, it
1574 // doesn't make sense to export them like an unbinned dataset, because the
1575 // coordinates are redundant information with the binning. We only do this
1576 // for 1D data for now.
1577 if (data.isWeighted() && variables.size() == 1) {
1578 bool isBinnedData = false;
1579 auto &x = static_cast<RooRealVar const &>(*variables[0]);
1580 std::vector<double> contents;
1581 int i = 0;
1582 for (; i < data.numEntries(); ++i) {
1583 data.get(i);
1584 if (x.getBin() != i)
1585 break;
1586 contents.push_back(data.weight());
1587 }
1588 if (i == x.getBins())
1589 isBinnedData = true;
1590 if (isBinnedData) {
1591 output["type"] << "binned";
1592 for (auto *var : static_range_cast<RooRealVar *>(variables)) {
1593 _domains->readVariable(*var);
1594 }
1595 return exportHisto(variables, data.numEntries(), contents.data(), output);
1596 }
1597 }
1598
1599 // this really is an unbinned dataset
1600 output["type"] << "unbinned";
1601 exportVariables(variables, output["axes"]);
1602 auto &coords = output["entries"].set_seq();
1603 std::vector<double> weightVals;
1604 bool hasNonUnityWeights = false;
1605 for (int i = 0; i < data.numEntries(); ++i) {
1606 data.get(i);
1607 coords.append_child().fill_seq(variables, [](auto x) { return static_cast<RooRealVar *>(x)->getVal(); });
1608 std::string datasetName = data.GetName();
1609 /*if (datasetName.find("combData_ZvvH126.5") != std::string::npos) {
1610 file << dynamic_cast<RooAbsReal *>(data.get(i)->find("atlas_invMass_PttEtaConvVBFCat1"))->getVal() <<
1611 std::endl;
1612 }*/
1613 if (data.isWeighted()) {
1614 weightVals.push_back(data.weight());
1615 if (data.weight() != 1.)
1616 hasNonUnityWeights = true;
1617 }
1618 }
1619 if (data.isWeighted() && hasNonUnityWeights) {
1620 output["weights"].fill_seq(weightVals);
1621 }
1622 // file.close();
1623}
1624
1625/**
1626 * @brief Read axes from the JSONNode and create a RooArgSet representing them.
1627 *
1628 * This function reads axes information from the given JSONNode and
1629 * creates a RooArgSet with variables representing these axes.
1630 *
1631 * @param topNode The JSONNode containing the axes information to be read.
1632 * @return RooArgSet A RooArgSet containing the variables created from the JSONNode.
1633 */
1635{
1636 RooArgSet vars;
1637
1638 for (JSONNode const &node : topNode["axes"].children()) {
1639 if (node.has_child("edges")) {
1640 std::vector<double> edges;
1641 for (auto const &bound : node["edges"].children()) {
1642 edges.push_back(bound.val_double());
1643 }
1644 auto obs = std::make_unique<RooRealVar>(node["name"].val().c_str(), node["name"].val().c_str(), edges[0],
1645 edges[edges.size() - 1]);
1646 RooBinning bins(obs->getMin(), obs->getMax());
1647 for (auto b : edges) {
1648 bins.addBoundary(b);
1649 }
1650 obs->setBinning(bins);
1651 vars.addOwned(std::move(obs));
1652 } else {
1653 auto obs = std::make_unique<RooRealVar>(node["name"].val().c_str(), node["name"].val().c_str(),
1654 node["min"].val_double(), node["max"].val_double());
1655 obs->setBins(node["nbins"].val_int());
1656 vars.addOwned(std::move(obs));
1657 }
1658 }
1659
1660 return vars;
1661}
1662
1663/**
1664 * @brief Read binned data from the JSONNode and create a RooDataHist object.
1665 *
1666 * This function reads binned data from the given JSONNode and creates a RooDataHist object.
1667 * The binned data is associated with the specified name and variables (RooArgSet) in the workspace.
1668 *
1669 * @param n The JSONNode representing the binned data to be read.
1670 * @param name The name to be associated with the created RooDataHist object.
1671 * @param vars The RooArgSet representing the variables associated with the binned data.
1672 * @return std::unique_ptr<RooDataHist> A unique pointer to the created RooDataHist object.
1673 */
1674std::unique_ptr<RooDataHist>
1675RooJSONFactoryWSTool::readBinnedData(const JSONNode &n, const std::string &name, RooArgSet const &vars)
1676{
1677 if (!n.has_child("contents"))
1678 RooJSONFactoryWSTool::error("no contents given");
1679
1680 JSONNode const &contents = n["contents"];
1681
1682 if (!contents.is_seq())
1683 RooJSONFactoryWSTool::error("contents are not in list form");
1684
1685 JSONNode const *errors = nullptr;
1686 if (n.has_child("errors")) {
1687 errors = &n["errors"];
1688 if (!errors->is_seq())
1689 RooJSONFactoryWSTool::error("errors are not in list form");
1690 }
1691
1692 auto bins = generateBinIndices(vars);
1693 if (contents.num_children() != bins.size()) {
1694 std::stringstream errMsg;
1695 errMsg << "inconsistent bin numbers: contents=" << contents.num_children() << ", bins=" << bins.size();
1697 }
1698 auto dh = std::make_unique<RooDataHist>(name, name, vars);
1699 std::vector<double> contentVals;
1700 contentVals.reserve(contents.num_children());
1701 for (auto const &cont : contents.children()) {
1702 contentVals.push_back(cont.val_double());
1703 }
1704 std::vector<double> errorVals;
1705 if (errors) {
1706 errorVals.reserve(errors->num_children());
1707 for (auto const &err : errors->children()) {
1708 errorVals.push_back(err.val_double());
1709 }
1710 }
1711 for (size_t ibin = 0; ibin < bins.size(); ++ibin) {
1712 const double err = errors ? errorVals[ibin] : -1;
1713 dh->set(ibin, contentVals[ibin], err);
1714 }
1715 return dh;
1716}
1717
1718/**
1719 * @brief Import a variable from the JSONNode into the workspace.
1720 *
1721 * This function imports a variable from the given JSONNode into the workspace.
1722 * The variable's information is read from the JSONNode and added to the workspace.
1723 *
1724 * @param p The JSONNode representing the variable to be imported.
1725 * @return void
1726 */
1728{
1729 // import a RooRealVar object
1730 std::string name(RooJSONFactoryWSTool::name(p));
1732
1733 if (_workspace.var(name))
1734 return;
1735 if (!p.is_map()) {
1736 std::stringstream ss;
1737 ss << "RooJSONFactoryWSTool() node '" << name << "' is not a map, skipping.";
1738 oocoutE(nullptr, InputArguments) << ss.str() << std::endl;
1739 return;
1740 }
1741 if (_attributesNode) {
1742 if (auto *attrNode = _attributesNode->find(name)) {
1743 // We should not create RooRealVar objects for RooConstVars!
1744 if (attrNode->has_child("is_const_var") && (*attrNode)["is_const_var"].val_int() == 1) {
1745 wsEmplace<RooConstVar>(name, p["value"].val_double());
1746 return;
1747 }
1748 }
1749 }
1751}
1752
1753/**
1754 * @brief Import all dependants (servers) of a node into the workspace.
1755 *
1756 * This function imports all the dependants (servers) of the given JSONNode into the workspace.
1757 * The dependants' information is read from the JSONNode and added to the workspace.
1758 *
1759 * @param n The JSONNode representing the node whose dependants are to be imported.
1760 * @return void
1761 */
1763{
1764 // import all the dependants of an object
1765 if (JSONNode const *varsNode = getVariablesNode(n)) {
1766 for (const auto &p : varsNode->children()) {
1768 }
1769 }
1770 if (auto seq = n.find("functions")) {
1771 for (const auto &p : seq->children()) {
1772 this->importFunction(p, true);
1773 }
1774 }
1775 if (auto seq = n.find("distributions")) {
1776 for (const auto &p : seq->children()) {
1777 this->importFunction(p, true);
1778 }
1779 }
1780}
1781
1783 const std::vector<CombinedData> &combDataSets,
1784 const std::vector<RooAbsData *> &singleDataSets)
1785{
1786 auto pdf = mc.GetPdf();
1787 auto simpdf = dynamic_cast<RooSimultaneous const *>(pdf);
1788 if (simpdf) {
1789 for (std::size_t i = 0; i < std::max(combDataSets.size(), std::size_t(1)); ++i) {
1790 const bool hasdata = i < combDataSets.size();
1791 if (hasdata && !matches(combDataSets.at(i), simpdf))
1792 continue;
1793
1794 std::string analysisName(simpdf->GetName());
1795 if (hasdata)
1796 analysisName += "_" + combDataSets[i].name;
1797
1798 exportSingleModelConfig(rootnode, mc, analysisName, hasdata ? &combDataSets[i].components : nullptr);
1799 }
1800 } else {
1801 RooArgSet observables(*mc.GetObservables());
1802 int founddata = 0;
1803 for (auto *data : singleDataSets) {
1804 if (observables.equals(*(data->get()))) {
1805 std::map<std::string, std::string> mapping;
1806 mapping[pdf->GetName()] = data->GetName();
1807 exportSingleModelConfig(rootnode, mc, std::string(pdf->GetName()) + "_" + data->GetName(), &mapping);
1808 ++founddata;
1809 }
1810 }
1811 if (founddata == 0) {
1812 exportSingleModelConfig(rootnode, mc, pdf->GetName(), nullptr);
1813 }
1814 }
1815}
1816
1818 std::string const &analysisName,
1819 std::map<std::string, std::string> const *dataComponents)
1820{
1821 auto pdf = mc.GetPdf();
1822
1823 JSONNode &analysisNode = appendNamedChild(rootnode["analyses"], analysisName);
1824
1825 auto &domains = analysisNode["domains"].set_seq();
1826
1827 analysisNode["likelihood"] << analysisName;
1828
1829 auto &nllNode = appendNamedChild(rootnode["likelihoods"], analysisName);
1830 nllNode["distributions"].set_seq();
1831 nllNode["data"].set_seq();
1832
1833 if (dataComponents) {
1834 auto simPdf = static_cast<RooSimultaneous const *>(pdf);
1835 if (simPdf) {
1836 for (auto const &item : simPdf->indexCat()) {
1837 const auto &dataComp = dataComponents->find(item.first);
1838 nllNode["distributions"].append_child() << simPdf->getPdf(item.first)->GetName();
1839 nllNode["data"].append_child() << dataComp->second;
1840 }
1841 } else {
1842 for (auto it : *dataComponents) {
1843 nllNode["distributions"].append_child() << it.first;
1844 nllNode["data"].append_child() << it.second;
1845 }
1846 }
1847 } else {
1848 nllNode["distributions"].append_child() << pdf->GetName();
1849 nllNode["data"].append_child() << 0;
1850 }
1851
1852 if (mc.GetExternalConstraints()) {
1853 auto &extConstrNode = nllNode["aux_distributions"];
1854 extConstrNode.set_seq();
1855 for (const auto &constr : *mc.GetExternalConstraints()) {
1856 extConstrNode.append_child() << constr->GetName();
1857 }
1858 }
1859
1860 auto writeList = [&](const char *name, RooArgSet const *args) {
1861 if (!args || !args->size())
1862 return;
1863
1864 std::vector<std::string> names;
1865 names.reserve(args->size());
1866 for (RooAbsArg const *arg : *args)
1867 names.push_back(arg->GetName());
1868 std::sort(names.begin(), names.end());
1869 analysisNode[name].fill_seq(names);
1870 };
1871
1872 writeList("parameters_of_interest", mc.GetParametersOfInterest());
1873
1874 auto &domainsNode = rootnode["domains"];
1875
1876 if (mc.GetNuisanceParameters() && mc.GetNuisanceParameters()->size() > 0) {
1877 std::string npDomainName = analysisName + "_nuisance_parameters";
1878 domains.append_child() << npDomainName;
1880 for (auto *np : static_range_cast<const RooRealVar *>(*mc.GetNuisanceParameters())) {
1881 npDomain.readVariable(*np);
1882 }
1884 }
1885
1886 if (mc.GetGlobalObservables() && mc.GetGlobalObservables()->size() > 0) {
1887 std::string globDomainName = analysisName + "_global_observables";
1888 domains.append_child() << globDomainName;
1890 for (auto *glob : static_range_cast<const RooRealVar *>(*mc.GetGlobalObservables())) {
1891 globDomain.readVariable(*glob);
1892 }
1894 }
1895
1896 if (mc.GetParametersOfInterest() && mc.GetParametersOfInterest()->size() > 0) {
1897 std::string poiDomainName = analysisName + "_parameters_of_interest";
1898 domains.append_child() << poiDomainName;
1900 for (auto *poi : static_range_cast<const RooRealVar *>(*mc.GetParametersOfInterest())) {
1901 poiDomain.readVariable(*poi);
1902 }
1904 }
1905
1906 auto &modelConfigAux = getRooFitInternal(rootnode, "ModelConfigs", analysisName);
1907 modelConfigAux.set_map();
1908 modelConfigAux["pdfName"] << pdf->GetName();
1909 modelConfigAux["mcName"] << mc.GetName();
1910}
1911
1912/**
1913 * @brief Export all objects in the workspace to a JSONNode.
1914 *
1915 * This function exports all the objects in the workspace to the provided JSONNode.
1916 * The objects' information is added as key-value pairs to the JSONNode.
1917 *
1918 * @param n The JSONNode to which the objects will be exported.
1919 * @return void
1920 */
1922{
1923 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
1925 _rootnodeOutput = &n;
1926
1927 // export all toplevel pdfs
1928 std::vector<RooAbsPdf *> allpdfs;
1929 for (auto &arg : _workspace.allPdfs()) {
1930 if (!arg->hasClients()) {
1931 if (auto *pdf = dynamic_cast<RooAbsPdf *>(arg)) {
1932 allpdfs.push_back(pdf);
1933 }
1934 }
1935 }
1937 std::set<std::string> exportedObjectNames;
1939
1940 // export all toplevel functions
1941 std::vector<RooAbsReal *> allfuncs;
1942 for (auto &arg : _workspace.allFunctions()) {
1943 if (!arg->hasClients()) {
1944 if (auto *func = dynamic_cast<RooAbsReal *>(arg)) {
1945 allfuncs.push_back(func);
1946 }
1947 }
1948 }
1951
1952 // export attributes of all objects
1953 for (RooAbsArg *arg : _workspace.components()) {
1954 exportAttributes(arg, n);
1955 }
1956
1957 // collect all datasets
1958 std::vector<RooAbsData *> alldata;
1959 for (auto &d : _workspace.allData()) {
1960 alldata.push_back(d);
1961 }
1963 // first, take care of combined datasets
1964 std::vector<RooAbsData *> singleData;
1965 std::vector<RooJSONFactoryWSTool::CombinedData> combData;
1966 for (auto &d : alldata) {
1967 auto data = this->exportCombinedData(*d);
1968 if (!data.components.empty())
1969 combData.push_back(data);
1970 else
1971 singleData.push_back(d);
1972 }
1973 // next, take care datasets
1974 for (auto &d : alldata) {
1975 this->exportData(*d);
1976 }
1977
1978 // export all ModelConfig objects and attached Pdfs
1979 for (TObject *obj : _workspace.allGenericObjects()) {
1980 if (auto mc = dynamic_cast<RooStats::ModelConfig *>(obj)) {
1982 }
1983 }
1984
1987 // We only want to add the variables that actually got exported and skip
1988 // the ones that the pdfs encoded implicitly (like in the case of
1989 // HistFactory).
1990 for (RooAbsArg *arg : *snsh) {
1991 if (exportedObjectNames.find(arg->GetName()) != exportedObjectNames.end()) {
1992 bool do_export = false;
1993 for (const auto &pdf : allpdfs) {
1994 if (pdf->dependsOn(*arg)) {
1995 do_export = true;
1996 }
1997 }
1998 if (do_export) {
1999 RooJSONFactoryWSTool::testValidName(arg->GetName(), true);
2000 snapshotSorted.add(*arg);
2001 }
2002 }
2003 }
2004 snapshotSorted.sort();
2005 std::string name(snsh->GetName());
2006 if (name != "default_values") {
2007 this->exportVariables(snapshotSorted, appendNamedChild(n["parameter_points"], name)["parameters"]);
2008 }
2009 }
2010 _varsNode = nullptr;
2011 _domains->writeJSON(n["domains"]);
2012 _domains.reset();
2013 _rootnodeOutput = nullptr;
2014}
2015
2016/**
2017 * @brief Import the workspace from a JSON string.
2018 *
2019 * @param s The JSON string containing the workspace data.
2020 * @return bool Returns true on successful import, false otherwise.
2021 */
2023{
2024 std::stringstream ss(s);
2025 return importJSON(ss);
2026}
2027
2028/**
2029 * @brief Import the workspace from a YML string.
2030 *
2031 * @param s The YML string containing the workspace data.
2032 * @return bool Returns true on successful import, false otherwise.
2033 */
2035{
2036 std::stringstream ss(s);
2037 return importYML(ss);
2038}
2039
2040/**
2041 * @brief Export the workspace to a JSON string.
2042 *
2043 * @return std::string The JSON string representing the exported workspace.
2044 */
2046{
2047 std::stringstream ss;
2048 exportJSON(ss);
2049 return ss.str();
2050}
2051
2052/**
2053 * @brief Export the workspace to a YML string.
2054 *
2055 * @return std::string The YML string representing the exported workspace.
2056 */
2058{
2059 std::stringstream ss;
2060 exportYML(ss);
2061 return ss.str();
2062}
2063
2064/**
2065 * @brief Create a new JSON tree with version information.
2066 *
2067 * @return std::unique_ptr<JSONTree> A unique pointer to the created JSON tree.
2068 */
2070{
2071 std::unique_ptr<JSONTree> tree = JSONTree::create();
2072 JSONNode &n = tree->rootnode();
2073 n.set_map();
2074 auto &metadata = n["metadata"].set_map();
2075
2076 // add the mandatory hs3 version number
2077 metadata["hs3_version"] << hs3VersionTag;
2078
2079 // Add information about the ROOT version that was used to generate this file
2080 auto &rootInfo = appendNamedChild(metadata["packages"], "ROOT");
2081 std::string versionName = gROOT->GetVersion();
2082 // We want to consistently use dots such that the version name can be easily
2083 // digested automatically.
2084 std::replace(versionName.begin(), versionName.end(), '/', '.');
2085 rootInfo["version"] << versionName;
2086
2087 return tree;
2088}
2089
2090/**
2091 * @brief Export the workspace to JSON format and write to the output stream.
2092 *
2093 * @param os The output stream to write the JSON data to.
2094 * @return bool Returns true on successful export, false otherwise.
2095 */
2097{
2098 std::unique_ptr<JSONTree> tree = createNewJSONTree();
2099 JSONNode &n = tree->rootnode();
2100 this->exportAllObjects(n);
2101 n.writeJSON(os);
2102 return true;
2103}
2104
2105/**
2106 * @brief Export the workspace to JSON format and write to the specified file.
2107 *
2108 * @param filename The name of the JSON file to create and write the data to.
2109 * @return bool Returns true on successful export, false otherwise.
2110 */
2112{
2113 std::ofstream out(filename.c_str());
2114 if (!out.is_open()) {
2115 std::stringstream ss;
2116 ss << "RooJSONFactoryWSTool() invalid output file '" << filename << "'." << std::endl;
2118 return false;
2119 }
2120 return this->exportJSON(out);
2121}
2122
2123/**
2124 * @brief Export the workspace to YML format and write to the output stream.
2125 *
2126 * @param os The output stream to write the YML data to.
2127 * @return bool Returns true on successful export, false otherwise.
2128 */
2130{
2131 std::unique_ptr<JSONTree> tree = createNewJSONTree();
2132 JSONNode &n = tree->rootnode();
2133 this->exportAllObjects(n);
2134 n.writeYML(os);
2135 return true;
2136}
2137
2138/**
2139 * @brief Export the workspace to YML format and write to the specified file.
2140 *
2141 * @param filename The name of the YML file to create and write the data to.
2142 * @return bool Returns true on successful export, false otherwise.
2143 */
2145{
2146 std::ofstream out(filename.c_str());
2147 if (!out.is_open()) {
2148 std::stringstream ss;
2149 ss << "RooJSONFactoryWSTool() invalid output file '" << filename << "'." << std::endl;
2151 return false;
2152 }
2153 return this->exportYML(out);
2154}
2155
2156bool RooJSONFactoryWSTool::hasAttribute(const std::string &obj, const std::string &attrib)
2157{
2158 if (!_attributesNode)
2159 return false;
2160 if (auto attrNode = _attributesNode->find(obj)) {
2161 if (auto seq = attrNode->find("tags")) {
2162 for (auto &a : seq->children()) {
2163 if (a.val() == attrib)
2164 return true;
2165 }
2166 }
2167 }
2168 return false;
2169}
2170void RooJSONFactoryWSTool::setAttribute(const std::string &obj, const std::string &attrib)
2171{
2172 auto node = &RooJSONFactoryWSTool::getRooFitInternal(*_rootnodeOutput, "attributes").set_map()[obj].set_map();
2173 auto &tags = (*node)["tags"];
2174 tags.set_seq();
2175 tags.append_child() << attrib;
2176}
2177
2178std::string RooJSONFactoryWSTool::getStringAttribute(const std::string &obj, const std::string &attrib)
2179{
2180 if (!_attributesNode)
2181 return "";
2182 if (auto attrNode = _attributesNode->find(obj)) {
2183 if (auto dict = attrNode->find("dict")) {
2184 if (auto *a = dict->find(attrib)) {
2185 return a->val();
2186 }
2187 }
2188 }
2189 return "";
2190}
2191void RooJSONFactoryWSTool::setStringAttribute(const std::string &obj, const std::string &attrib,
2192 const std::string &value)
2193{
2194 auto node = &RooJSONFactoryWSTool::getRooFitInternal(*_rootnodeOutput, "attributes").set_map()[obj].set_map();
2195 auto &dict = (*node)["dict"];
2196 dict.set_map();
2197 dict[attrib] << value;
2198}
2199
2200/**
2201 * @brief Imports all nodes of the JSON data and adds them to the workspace.
2202 *
2203 * @param n The JSONNode representing the root node of the JSON data.
2204 * @return void
2205 */
2207{
2208 // Per HS3 standard, the hs3_version in the metadata is required. So we
2209 // error out if it is missing. TODO: now we are only checking if the
2210 // hs3_version tag exists, but in the future when the HS3 specification
2211 // versions are actually frozen, we should also check if the hs3_version is
2212 // one that RooFit can actually read.
2213 auto metadata = n.find("metadata");
2214 if (!metadata || !metadata->find("hs3_version")) {
2215 std::stringstream ss;
2216 ss << "The HS3 version is missing in the JSON!\n"
2217 << "Please include the HS3 version in the metadata field, e.g.:\n"
2218 << " \"metadata\" :\n"
2219 << " {\n"
2220 << " \"hs3_version\" : \"" << hs3VersionTag << "\"\n"
2221 << " }";
2222 error(ss.str());
2223 }
2224
2225 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
2226 if (auto domains = n.find("domains")) {
2227 _domains->readJSON(*domains);
2228 }
2229 _domains->populate(_workspace);
2230
2231 _rootnodeInput = &n;
2232
2234
2235 this->importDependants(n);
2236
2237 if (auto paramPointsNode = n.find("parameter_points")) {
2238 for (const auto &snsh : paramPointsNode->children()) {
2239 std::string name = RooJSONFactoryWSTool::name(snsh);
2241
2242 RooArgSet vars;
2243 for (const auto &var : snsh["parameters"].children()) {
2246 vars.add(*rrv);
2247 }
2248 }
2250 }
2251 }
2252
2254
2255 // Import attributes
2256 if (_attributesNode) {
2257 for (const auto &elem : _attributesNode->children()) {
2258 if (RooAbsArg *arg = _workspace.arg(elem.key()))
2259 importAttributes(arg, elem);
2260 }
2261 }
2262
2263 _attributesNode = nullptr;
2264
2265 // We delay the import of the data to after combineDatasets(), because it
2266 // might be that some datasets are merged to combined datasets there. In
2267 // that case, we will remove the components from the "datasets" vector so they
2268 // don't get imported.
2269 std::vector<std::unique_ptr<RooAbsData>> datasets;
2270 if (auto dataNode = n.find("data")) {
2271 for (const auto &p : dataNode->children()) {
2272 datasets.push_back(loadData(p, _workspace));
2273 }
2274 }
2275
2276 // Now, read in analyses and likelihoods if there are any
2277
2278 if (auto analysesNode = n.find("analyses")) {
2279 for (JSONNode const &analysisNode : analysesNode->children()) {
2280 importAnalysis(*_rootnodeInput, analysisNode, n["likelihoods"], n["domains"], _workspace, datasets);
2281 }
2282 }
2283
2284 combineDatasets(*_rootnodeInput, datasets);
2285
2286 for (auto const &d : datasets) {
2287 if (d)
2289 }
2290
2291 _rootnodeInput = nullptr;
2292 _domains.reset();
2293}
2294
2295/**
2296 * @brief Imports a JSON file from the given input stream to the workspace.
2297 *
2298 * @param is The input stream containing the JSON data.
2299 * @return bool Returns true on successful import, false otherwise.
2300 */
2302{
2303 // import a JSON file to the workspace
2304 std::unique_ptr<JSONTree> tree = JSONTree::create(is);
2305 this->importAllNodes(tree->rootnode());
2306 if (this->workspace()->getSnapshot("default_values")) {
2307 this->workspace()->loadSnapshot("default_values");
2308 }
2309 return true;
2310}
2311
2312/**
2313 * @brief Imports a JSON file from the given filename to the workspace.
2314 *
2315 * @param filename The name of the JSON file to import.
2316 * @return bool Returns true on successful import, false otherwise.
2317 */
2319{
2320 // import a JSON file to the workspace
2321 std::ifstream infile(filename.c_str());
2322 if (!infile.is_open()) {
2323 std::stringstream ss;
2324 ss << "RooJSONFactoryWSTool() invalid input file '" << filename << "'." << std::endl;
2326 return false;
2327 }
2328 return this->importJSON(infile);
2329}
2330
2331/**
2332 * @brief Imports a YML file from the given input stream to the workspace.
2333 *
2334 * @param is The input stream containing the YML data.
2335 * @return bool Returns true on successful import, false otherwise.
2336 */
2338{
2339 // import a YML file to the workspace
2340 std::unique_ptr<JSONTree> tree = JSONTree::create(is);
2341 this->importAllNodes(tree->rootnode());
2342 return true;
2343}
2344
2345/**
2346 * @brief Imports a YML file from the given filename to the workspace.
2347 *
2348 * @param filename The name of the YML file to import.
2349 * @return bool Returns true on successful import, false otherwise.
2350 */
2352{
2353 // import a YML file to the workspace
2354 std::ifstream infile(filename.c_str());
2355 if (!infile.is_open()) {
2356 std::stringstream ss;
2357 ss << "RooJSONFactoryWSTool() invalid input file '" << filename << "'." << std::endl;
2359 return false;
2360 }
2361 return this->importYML(infile);
2362}
2363
2364void RooJSONFactoryWSTool::importJSONElement(const std::string &name, const std::string &jsonString)
2365{
2366 std::unique_ptr<RooFit::Detail::JSONTree> tree = RooFit::Detail::JSONTree::create(jsonString);
2367 JSONNode &n = tree->rootnode();
2368 n["name"] << name;
2369
2370 bool isVariable = true;
2371 if (n.find("type")) {
2372 isVariable = false;
2373 }
2374
2375 if (isVariable) {
2376 this->importVariableElement(n);
2377 } else {
2378 this->importFunction(n, false);
2379 }
2380}
2381
2383{
2384 std::unique_ptr<RooFit::Detail::JSONTree> tree = varJSONString(elementNode);
2385 JSONNode &n = tree->rootnode();
2386 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
2387 if (auto domains = n.find("domains"))
2388 _domains->readJSON(*domains);
2389
2390 _rootnodeInput = &n;
2392
2394 const auto &p = varsNode->child(0);
2396
2397 auto paramPointsNode = n.find("parameter_points");
2398 const auto &snsh = paramPointsNode->child(0);
2399 std::string name = RooJSONFactoryWSTool::name(snsh);
2400 RooArgSet vars;
2401 const auto &var = snsh["parameters"].child(0);
2404 vars.add(*rrv);
2405 }
2406
2407 // Import attributes
2408 if (_attributesNode) {
2409 for (const auto &elem : _attributesNode->children()) {
2410 if (RooAbsArg *arg = _workspace.arg(elem.key()))
2411 importAttributes(arg, elem);
2412 }
2413 }
2414
2415 _attributesNode = nullptr;
2416 _rootnodeInput = nullptr;
2417 _domains.reset();
2418}
2419
2420/**
2421 * @brief Writes a warning message to the RooFit message service.
2422 *
2423 * @param str The warning message to be logged.
2424 * @return std::ostream& A reference to the output stream.
2425 */
2426std::ostream &RooJSONFactoryWSTool::warning(std::string const &str)
2427{
2428 return RooMsgService::instance().log(nullptr, RooFit::MsgLevel::ERROR, RooFit::IO) << str << std::endl;
2429}
2430
2431/**
2432 * @brief Writes an error message to the RooFit message service and throws a runtime_error.
2433 *
2434 * @param s The error message to be logged and thrown.
2435 * @return void
2436 */
2438{
2439 RooMsgService::instance().log(nullptr, RooFit::MsgLevel::ERROR, RooFit::IO) << s << std::endl;
2440 throw std::runtime_error(s);
2441}
2442
2443/**
2444 * @brief Cleans up names to the HS3 standard
2445 *
2446 * @param str The string to be sanitized.
2447 * @return std::string
2448 */
2449std::string RooJSONFactoryWSTool::sanitizeName(const std::string str)
2450{
2451 std::string result;
2453 for (char c : str) {
2454 switch (c) {
2455 case '[':
2456 case '|':
2457 case ',':
2458 case '(': result += '_'; break;
2459 case ']':
2460 case ')':
2461 // skip these characters entirely
2462 break;
2463 case '.': result += "_dot_"; break;
2464 case '@': result += "at"; break;
2465 case '-': result += "minus"; break;
2466 case '/': result += "_div_"; break;
2467
2468 default: result += c; break;
2469 }
2470 }
2471 return result;
2472 }
2473 return str;
2474}
2475
2477{
2478 // Variables
2479
2481 if (onlyModelConfig) {
2482 for (auto *obj : ws.allGenericObjects()) {
2483 if (auto *mc = dynamic_cast<RooStats::ModelConfig *>(obj)) {
2484 tmpWS.import(*mc->GetPdf(), RooFit::RecycleConflictNodes(true));
2485 }
2486 }
2487
2488 } else {
2489
2490 for (auto *pdf : ws.allPdfs()) {
2491 if (!pdf->hasClients()) {
2492 tmpWS.import(*pdf, RooFit::RecycleConflictNodes(true));
2493 }
2494 }
2495
2496 for (auto *func : ws.allFunctions()) {
2497 if (!func->hasClients()) {
2498 tmpWS.import(*func, RooFit::RecycleConflictNodes(true));
2499 }
2500 }
2501 }
2502
2503 for (auto *data : ws.allData()) {
2504 tmpWS.import(*data);
2505 }
2506
2507 for (auto *obj : ws.allGenericObjects()) {
2508 tmpWS.import(*obj);
2509 }
2510
2511 /*
2512 if (auto* mc = dynamic_cast<RooStats::ModelConfig*>(obj)) {
2513 // Import the PDF
2514 tmpWS.import(*mc->GetPdf());
2515
2516 // Import all observables
2517 RooArgSet* obs = (RooArgSet*)mc->GetObservables()->snapshot();
2518 tmpWS.import(*obs);
2519
2520 // Import global observables
2521 RooArgSet* globObs = (RooArgSet*)mc->GetGlobalObservables()->snapshot();
2522 tmpWS.import(*globObs);
2523
2524 // Import POIs
2525 RooArgSet* pois = (RooArgSet*)mc->GetParametersOfInterest()->snapshot();
2526 tmpWS.import(*pois);
2527
2528 // Import nuisance parameters
2529 RooArgSet* nuis = (RooArgSet*)mc->GetNuisanceParameters()->snapshot();
2530 tmpWS.import(*nuis);
2531
2532
2533 RooStats::ModelConfig* mc_new = new RooStats::ModelConfig(mc->GetName(), mc->GetName());
2534
2535 mc_new->SetPdf(*tmpWS.pdf(mc->GetPdf()->GetName()));
2536 mc_new->SetObservables(*tmpWS.set(obs->GetName()));
2537 mc_new->SetGlobalObservables(*tmpWS.set(globObs->GetName()));
2538 mc_new->SetParametersOfInterest(*tmpWS.set(pois->GetName()));
2539 mc_new->SetNuisanceParameters(*tmpWS.set(nuis->GetName()));
2540
2541 // Import the ModelConfig into the new workspace
2542 tmpWS.import(*mc_new);
2543 }else {
2544
2545 tmpWS.import(*obj);
2546 }
2547 */
2548
2549 for (auto *snsh : ws.getSnapshots()) {
2550 auto *snshSet = dynamic_cast<RooArgSet *>(snsh);
2551 if (snshSet) {
2552 tmpWS.saveSnapshot(snshSet->GetName(), *snshSet, true);
2553 }
2554 }
2555
2556 return tmpWS;
2557}
2558
2559// Sanitize all names in the workspace to be HS3 compliant
2561{
2562 // Variables
2563
2564 RooWorkspace tmpWS = cleanWS(ws, false);
2565
2566 for (auto *obj : tmpWS.allVars()) {
2567 if (!isValidName(obj->GetName())) {
2568 obj->SetName(sanitizeName(obj->GetName()).c_str());
2569 }
2570 }
2571
2572 // Functions
2573 for (auto *obj : tmpWS.allFunctions()) {
2574 if (!isValidName(obj->GetName())) {
2575 obj->SetName(sanitizeName(obj->GetName()).c_str());
2576 }
2577 }
2578
2579 // PDFs
2580 for (auto *obj : tmpWS.allPdfs()) {
2581 if (!isValidName(obj->GetName())) {
2582 obj->SetName(sanitizeName(obj->GetName()).c_str());
2583 }
2584 }
2585
2586 // Datasets
2587 for (auto *data : tmpWS.allData()) {
2588 // Sanitize dataset name
2589 if (!isValidName(data->GetName())) {
2590 data->SetName(sanitizeName(data->GetName()).c_str());
2591 }
2592 for (auto *obj : *data->get()) {
2593 obj->SetName(sanitizeName(obj->GetName()).c_str());
2594 }
2595 }
2596 /* // Sanitize dataset observables
2597 const RooArgSet* obsSet = data->get();
2598 if (obsSet) {
2599 RooArgSet* mutableObs = const_cast<RooArgSet*>(obsSet);
2600 std::string oldSetName = mutableObs->GetName();
2601 std::string newSetName = sanitizeName(oldSetName);
2602 if (oldSetName != newSetName) {
2603 mutableObs->setName(newSetName.c_str());
2604 }
2605 }
2606
2607 for (auto* arg : *obsSet) {
2608 std::string oldObsName = arg->GetName();
2609 std::string newObsName = sanitizeName(oldObsName);
2610 if (oldObsName != newObsName) {
2611 arg->SetName(newObsName.c_str());
2612 data->changeObservableName(arg->GetName(), newObsName.c_str());
2613 }
2614 }
2615 */
2616 for (auto *data : tmpWS.allEmbeddedData()) {
2617 // Sanitize dataset name
2618 data->SetName(sanitizeName(data->GetName()).c_str());
2619 for (auto *obj : *data->get()) {
2620 obj->SetName(sanitizeName(obj->GetName()).c_str());
2621 }
2622 }
2623 for (auto *snshObj : tmpWS.getSnapshots()) {
2624 // Snapshots are stored as TObject*, but really they are RooArgSet*
2625 auto *snsh = dynamic_cast<RooArgSet *>(snshObj);
2626 if (!snsh) {
2627 std::cerr << "Warning: found snapshot that is not a RooArgSet, skipping\n";
2628 continue;
2629 }
2630
2631 // Sanitize snapshot name
2632 if (!isValidName(snsh->GetName())) {
2633 snsh->setName(sanitizeName(snsh->GetName()).c_str());
2634 }
2635
2636 // Sanitize the variables inside the snapshot
2637 for (auto *arg : *snsh) {
2638 if (!isValidName(arg->GetName())) {
2639 arg->SetName(sanitizeName(arg->GetName()).c_str());
2640 }
2641 }
2642 }
2643
2644 // Generic objects (ModelConfigs, attributes, etc.)
2645 for (auto *obj : tmpWS.allGenericObjects()) {
2646 if (!isValidName(obj->GetName())) {
2647 if (auto *named = dynamic_cast<TNamed *>(obj)) {
2648 named->SetName(sanitizeName(named->GetName()).c_str());
2649 } else {
2650 std::cerr << "Warning: object " << obj->GetName() << " is not TNamed, cannot rename.\n";
2651 }
2652 }
2653
2654 if (auto *mc = dynamic_cast<RooStats::ModelConfig *>(obj)) {
2655 // Sanitize ModelConfig name
2656 if (!isValidName(mc->GetName())) {
2657 mc->SetName(sanitizeName(mc->GetName()).c_str());
2658 }
2659
2660 // Sanitize the sets inside ModelConfig
2661 for (auto *obs : mc->GetObservables()->get()) {
2662 if (obs) {
2663 obs->SetName(sanitizeName(obs->GetName()).c_str());
2664 }
2665 }
2666 for (auto *poi : mc->GetParametersOfInterest()->get()) {
2667 if (poi) {
2668 poi->SetName(sanitizeName(poi->GetName()).c_str());
2669 }
2670 }
2671 for (auto *nuis : mc->GetNuisanceParameters()->get()) {
2672 if (nuis) {
2673 nuis->SetName(sanitizeName(nuis->GetName()).c_str());
2674 }
2675 }
2676 for (auto *glob : mc->GetGlobalObservables()->get()) {
2677 if (glob) {
2678 glob->SetName(sanitizeName(glob->GetName()).c_str());
2679 }
2680 }
2681 }
2682 }
2683 std::string wsName = std::string{ws.GetName()} + "_sanitized";
2684 RooWorkspace newWS = cleanWS(tmpWS, false);
2685 newWS.SetName(wsName.c_str());
2686
2687 return newWS;
2688}
std::unique_ptr< RooFit::Detail::JSONTree > varJSONString(const JSONNode &treeRoot)
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
double toDouble(const char *s)
constexpr auto hs3VersionTag
#define oocoutW(o, a)
#define oocoutE(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
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 r
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 result
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 child
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t attr
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
#define gROOT
Definition TROOT.h:411
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
TClass * IsA() const override
Definition RooAbsArg.h:678
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
const std::set< std::string > & attributes() const
Definition RooAbsArg.h:258
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:145
const std::map< std::string, std::string > & stringAttributes() const
Definition RooAbsArg.h:267
Int_t numProxies() const
Return the number of registered proxies.
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
RooAbsProxy * getProxy(Int_t index) const
Return the nth proxy from the proxy list.
A space to attach TBranches.
Abstract container object that can hold multiple RooAbsArg objects.
bool equals(const RooAbsCollection &otherColl) const
Check if this and other collection have identically-named contents.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
std::unique_ptr< RooArgSet > getAllConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, bool stripDisconnected=true) const
This helper function finds and collects all constraints terms of all component p.d....
Abstract interface for proxy classes.
Definition RooAbsProxy.h:37
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
Abstract interface for RooAbsArg proxy classes.
Definition RooArgProxy.h:24
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Implements a RooAbsBinning in terms of an array of boundary values, posing no constraints on the choi...
Definition RooBinning.h:27
bool addBoundary(double boundary)
Add bin boundary at given value.
Object to represent discrete states.
Definition RooCategory.h:28
Represents a constant real-valued object.
Definition RooConstVar.h:23
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
virtual double val_double() const
virtual JSONNode & set_map()=0
virtual JSONNode & append_child()=0
virtual children_view children()
virtual size_t num_children() const =0
virtual JSONNode & set_seq()=0
virtual bool is_seq() const =0
virtual bool is_map() const =0
virtual std::string key() const =0
JSONNode const * find(std::string const &key) const
virtual int val_int() const
static std::unique_ptr< JSONTree > create()
When using RooFit, statistical models can be conveniently handled and stored as a RooWorkspace.
void importFunction(const RooFit::Detail::JSONNode &n, bool importAllDependants)
Import a function from the JSONNode into the workspace.
static constexpr bool useListsInsteadOfDicts
std::string getStringAttribute(const std::string &obj, const std::string &attrib)
bool importYML(std::string const &filename)
Imports a YML file from the given filename to the workspace.
static void fillSeq(RooFit::Detail::JSONNode &node, RooAbsCollection const &coll, size_t nMax=-1)
void exportObjects(T const &args, std::set< std::string > &exportedObjectNames)
void exportCategory(RooAbsCategory const &cat, RooFit::Detail::JSONNode &node)
Export a RooAbsCategory object to a JSONNode.
RooJSONFactoryWSTool(RooWorkspace &ws)
void importVariable(const RooFit::Detail::JSONNode &n)
Import a variable from the JSONNode into the workspace.
void exportData(RooAbsData const &data)
Export data from the workspace to a JSONNode.
void exportModelConfig(RooFit::Detail::JSONNode &rootnode, RooStats::ModelConfig const &mc, const std::vector< RooJSONFactoryWSTool::CombinedData > &combined, const std::vector< RooAbsData * > &single)
bool hasAttribute(const std::string &obj, const std::string &attrib)
void exportVariables(const RooArgSet &allElems, RooFit::Detail::JSONNode &n)
Export variables from the workspace to a JSONNode.
bool importJSON(std::string const &filename)
Imports a JSON file from the given filename to the workspace.
static std::unique_ptr< RooDataHist > readBinnedData(const RooFit::Detail::JSONNode &n, const std::string &namecomp, RooArgSet const &vars)
Read binned data from the JSONNode and create a RooDataHist object.
static RooFit::Detail::JSONNode & appendNamedChild(RooFit::Detail::JSONNode &node, std::string const &name)
std::string exportYMLtoString()
Export the workspace to a YML string.
static RooFit::Detail::JSONNode & getRooFitInternal(RooFit::Detail::JSONNode &node, Keys_t const &...keys)
static void exportArray(std::size_t n, double const *contents, RooFit::Detail::JSONNode &output)
Export an array of doubles to a JSONNode.
bool importYMLfromString(const std::string &s)
Import the workspace from a YML string.
static bool testValidName(const std::string &str, bool forcError)
RooFit::Detail::JSONNode * _rootnodeOutput
static void exportHisto(RooArgSet const &vars, std::size_t n, double const *contents, RooFit::Detail::JSONNode &output)
Export histogram data to a JSONNode.
std::vector< RooAbsArg const * > _serversToDelete
void exportSingleModelConfig(RooFit::Detail::JSONNode &rootnode, RooStats::ModelConfig const &mc, std::string const &analysisName, std::map< std::string, std::string > const *dataComponents)
static void fillSeqSanitizedName(RooFit::Detail::JSONNode &node, RooAbsCollection const &coll, size_t nMax=-1)
static std::unique_ptr< RooFit::Detail::JSONTree > createNewJSONTree()
Create a new JSON tree with version information.
void exportVariable(const RooAbsArg *v, RooFit::Detail::JSONNode &n)
Export a variable from the workspace to a JSONNode.
const RooFit::Detail::JSONNode * _rootnodeInput
RooJSONFactoryWSTool::CombinedData exportCombinedData(RooAbsData const &data)
Export combined data from the workspace to a custom struct.
std::string exportJSONtoString()
Export the workspace to a JSON string.
static RooWorkspace cleanWS(const RooWorkspace &ws, bool onlyModelConfig=false)
std::string exportTransformed(const RooAbsReal *original, const std::string &suffix, const std::string &formula)
const RooFit::Detail::JSONNode * _attributesNode
static bool isValidName(const std::string &str)
Check if a string is a valid name.
void importDependants(const RooFit::Detail::JSONNode &n)
Import all dependants (servers) of a node into the workspace.
void importJSONElement(const std::string &name, const std::string &jsonString)
static RooWorkspace sanitizeWS(const RooWorkspace &ws)
static void error(const char *s)
Writes an error message to the RooFit message service and throws a runtime_error.
void setAttribute(const std::string &obj, const std::string &attrib)
bool exportYML(std::string const &fileName)
Export the workspace to YML format and write to the specified file.
bool importJSONfromString(const std::string &s)
Import the workspace from a JSON string.
RooFit::Detail::JSONNode * _varsNode
void exportObject(RooAbsArg const &func, std::set< std::string > &exportedObjectNames)
Export an object from the workspace to a JSONNode.
static RooFit::Detail::JSONNode & makeVariablesNode(RooFit::Detail::JSONNode &rootNode)
static std::string sanitizeName(const std::string str)
Cleans up names to the HS3 standard.
void importAllNodes(const RooFit::Detail::JSONNode &n)
Imports all nodes of the JSON data and adds them to the workspace.
static std::string name(const RooFit::Detail::JSONNode &n)
void exportAllObjects(RooFit::Detail::JSONNode &n)
Export all objects in the workspace to a JSONNode.
bool exportJSON(std::string const &fileName)
Export the workspace to JSON format and write to the specified file.
static RooFit::Detail::JSONNode const * findNamedChild(RooFit::Detail::JSONNode const &node, std::string const &name)
void setStringAttribute(const std::string &obj, const std::string &attrib, const std::string &value)
std::vector< RooAbsArg const * > _serversToExport
std::unique_ptr< RooFit::JSONIO::Detail::Domains > _domains
static std::ostream & warning(const std::string &s)
Writes a warning message to the RooFit message service.
static RooArgSet readAxes(const RooFit::Detail::JSONNode &node)
Read axes from the JSONNode and create a RooArgSet representing them.
void importVariableElement(const RooFit::Detail::JSONNode &n)
static RooMsgService & instance()
Return reference to singleton instance.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
< A class that holds configuration information for a model using a workspace as a store
Definition ModelConfig.h:34
Persistable container for RooFit projects.
TObject * obj(RooStringView name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
const RooArgSet * getSnapshot(const char *name) const
Return the RooArgSet containing a snapshot of variables contained in the workspace.
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
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.
std::list< RooAbsData * > allData() const
Return list of all dataset in the workspace.
RooLinkedList const & getSnapshots() const
std::list< TObject * > allGenericObjects() const
Return list of all generic objects in the workspace.
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.
const RooArgSet & components() const
RooArgSet allFunctions() const
Return set with all function objects.
RooFactoryWSTool & factory()
Return instance to factory tool.
RooRealVar * var(RooStringView name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
bool loadSnapshot(const char *name)
Load the values and attributes of the parameters in the snapshot saved with the given name.
bool import(const RooAbsArg &arg, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={})
Import a RooAbsArg object, e.g.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Mother of all ROOT objects.
Definition TObject.h:41
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
RooCmdArg RecycleConflictNodes(bool flag=true)
RooConstVar & RooConst(double val)
RooCmdArg Silence(bool flag=true)
RooCmdArg Index(RooCategory &icat)
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Import(const char *state, TH1 &histo)
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Double_t ex[n]
Definition legend1.C:17
std::string makeValidVarName(std::string const &in)
ImportMap & importers()
Definition JSONIO.cxx:53
ExportMap & exporters()
Definition JSONIO.cxx:59
ImportExpressionMap & importExpressions()
Definition JSONIO.cxx:65
ExportKeysMap & exportKeys()
Definition JSONIO.cxx:72
TLine l
Definition textangle.C:4
static void output()